VTKFIG  0.20.20181114
Easy VTK based in situ visualization
template<typename V >
void vtkfig::DataSet::SetRectilinearGrid ( const V &  x,
const V &  y,
const V &  z 
)
inline

Enter data of a 3D rectilinear grid.

Template Parameters
VVector class counting from zero with member functions size() and operator[]. std::vector will work.
Parameters
xx-coordinates
yy-coordinates
zz-coordinates

Definition at line 498 of file vtkfigDataSet.h.

499  {
500  vtkSmartPointer<vtkDoubleArray> xcoord;
501  vtkSmartPointer<vtkDoubleArray> ycoord;
502  vtkSmartPointer<vtkDoubleArray> zcoord;
503  assert(spacedim!=2);
504  spacedim=3;
505  int Nx = x.size();
506  int Ny = y.size();
507  int Nz = z.size();
508 
509  if (this->data==NULL)
510  {
511  this->data=vtkSmartPointer<vtkRectilinearGrid>::New();
512  auto rdata=vtkRectilinearGrid::SafeDownCast(this->data);
513  assert(rdata);
514  xcoord = vtkSmartPointer<vtkDoubleArray>::New();
515  ycoord = vtkSmartPointer<vtkDoubleArray>::New();
516  zcoord = vtkSmartPointer<vtkDoubleArray>::New();
517  rdata->SetXCoordinates(xcoord);
518  rdata->SetYCoordinates(ycoord);
519  rdata->SetZCoordinates(zcoord);
520  rdata->SetDimensions(Nx, Ny, Nz );
521  }
522  else
523  {
524  auto rdata=vtkRectilinearGrid::SafeDownCast(this->data);
525  assert(rdata);
526  xcoord=vtkDoubleArray::SafeDownCast(rdata->GetXCoordinates());
527  ycoord=vtkDoubleArray::SafeDownCast(rdata->GetYCoordinates());
528  zcoord=vtkDoubleArray::SafeDownCast(rdata->GetZCoordinates());
529  xcoord->Initialize();
530  ycoord->Initialize();
531  zcoord->Initialize();
532  rdata->SetDimensions(Nx, Ny, Nz );
533  }
534 
535  xcoord->SetNumberOfComponents(1);
536  xcoord->SetNumberOfTuples(Nx);
537 
538  ycoord->SetNumberOfComponents(1);
539  ycoord->SetNumberOfTuples(Ny);
540 
541  zcoord->SetNumberOfComponents(1);
542  zcoord->SetNumberOfTuples(Nz);
543 
544  for (int i=0; i<Nx; i++)
545  xcoord->InsertComponent(i, 0, x[i]*coordinate_scale_factor*coordinate_scale_factor_xyz[0]);
546  for (int i=0; i<Ny; i++)
547  ycoord->InsertComponent(i, 0, y[i]*coordinate_scale_factor*coordinate_scale_factor_xyz[1]);
548  for (int i=0; i<Nz; i++)
549  zcoord->InsertComponent(i, 0, z[i]*coordinate_scale_factor*coordinate_scale_factor_xyz[2]);
550 
551  this->data->Modified();
552  }
vtkSmartPointer< vtkDataSet > data
double coordinate_scale_factor_xyz[3]
double coordinate_scale_factor