VTKFIG  0.25.0
Easy VTK based in situ visualization
example-simplexquiver3d.cxx

Vector function on 3D simplex grid.

#include "vtkMath.h"
#include "vtkfigFrame.h"
#include "vtkfigTools.h"
inline double G(double x,double y, double z, double t)
{
return exp(-(x*x+y*y+z*z))*sin(t+x)*cos(y-t)*cos(z-2*t);
}
inline double dGdx(double x,double y, double z, double t)
{
return cos(y-t)*cos(z-2*t)*exp(-(x*x+y*y+z*z))*(cos(t+x)-2*x*sin(t+x));
}
inline double dGdy(double x,double y, double z, double t)
{
return sin(t+x)*cos(z-2*t)*exp(-(x*x+y*y+z*z))*(-2*y*cos(y-t) -sin(y-t));
}
inline double dGdz(double x,double y, double z, double t)
{
return sin(t+x)*cos(y-t)*exp(-(x*x+y*y+z*z))*(-2*z*cos(z-2*t) -sin(z-2*t));
}
int main(void)
{
size_t nspin=vtkfig::NSpin();
std::vector<double> inpoints;
for(double x = -2; x < 2; x+=0.2)
for(double y = -2; y < 2; y+=0.2)
for(double z = -2; z < 2; z+=0.2)
{
inpoints.push_back(x + vtkMath::Random(-.1, .1));
inpoints.push_back(y + vtkMath::Random(-.1, .1));
inpoints.push_back(z + vtkMath::Random(-.1, .1));
}
std::vector<double>points;
std::vector<int>cells;
vtkfig::Delaunay3D(inpoints,points,cells);
int npoints=points.size()/3;
std::vector<double> u(npoints);
std::vector<double> v(npoints);
std::vector<double> w(npoints);
std::vector<double>val(npoints);
auto colors=vtkfig::RGBTable
{
{0.0, 0.3, 0.3, 1.0},
{0.5, 0.3, 1.0, 0.3},
{1.0, 1.0, 0.3, 0.3}
};
auto qcolors=vtkfig::RGBTable
{
{0.0, 0.0, 0.0, 0.0},
{1.0, 0.0, 0.0, 0.0}
};
auto frame=vtkfig::Frame::New();
auto dataset=vtkfig::DataSet::New();
dataset->SetSimplexGrid(3,points,cells);
dataset->SetPointScalar(val ,"val");
dataset->SetPointVector(u,v,w ,"grad");
auto contour=vtkfig::ScalarView::New();
contour->SetData(dataset,"val");
contour->SetSurfaceRGBTable(colors,255);
contour->ShowIsolines(false);
contour->SetValueRange(-1,1);
auto quiver=vtkfig::VectorView::New();
quiver->SetData(dataset,"grad");
quiver->SetQuiverGrid(10,10,10);
frame->AddFigure(contour);
frame->AddFigure(quiver);
double t=0;
double dt=0.1;
size_t ii=0;
double t0=(double)clock()/(double)CLOCKS_PER_SEC;
double i0=ii;
while (ii<nspin)
{
for (size_t ipoint=0, ival=0;ipoint<points.size(); ipoint+=3,ival++)
{
val[ival]=G(points[ipoint+0],points[ipoint+1],points[ipoint+2],t);
u[ival]=dGdx(points[ipoint+0],points[ipoint+1],points[ipoint+2],t);
v[ival]=dGdy(points[ipoint+0],points[ipoint+1],points[ipoint+2],t);
w[ival]=dGdz(points[ipoint+0],points[ipoint+1],points[ipoint+2],t);
}
dataset->SetPointScalar(val ,"val");
dataset->SetPointVector(u,v ,w,"grad");
frame->Show();
if (ii==3)
frame->WritePNG("example-simplexquiver3d.png");
t+=dt;
double t1=(double)clock()/(double)CLOCKS_PER_SEC;
double i1=ii;
if (t1-t0>4.0)
{
printf("Frame rate: %.2f fps\n",(double)(i1-i0)/4.0);
t0=(double)clock()/(double)CLOCKS_PER_SEC;
i0=ii;
}
ii++;
}
}
vtkfig::VectorView::New
static std::shared_ptr< VectorView > New()
Smart pointer constructor.
vtkfigScalarView.h
vtkfig::Delaunay3D
void Delaunay3D(const V0 &points_in, V &points, IV &cells)
Create a 3D Delaunay triangulation from point set.
Definition: vtkfigTools.h:115
vtkfig::ScalarView::New
static std::shared_ptr< ScalarView > New()
Create smart pointer to ScalarView instance.
vtkfig::NSpin
size_t NSpin()
Return content of environment variable VTKFIG_NSPIN or large value.
vtkfig::RGBTable
std::vector< RGBPoint > RGBTable
Color table.
Definition: vtkfigTools.h:43
vtkfigVectorView.h
vtkfigTools.h
vtkfig::DataSet::New
static std::shared_ptr< DataSet > New()
Static constructor of smart pointer to an empty instance.
vtkfigFrame.h
vtkfig::Frame::New
static std::shared_ptr< Frame > New()
Create smart pointer to frame.