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

Enter data of a 2D rectilinear grid.

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

Definition at line 450 of file vtkfigDataSet.h.

451  {
452  vtkSmartPointer<vtkDoubleArray> xcoord;
453  vtkSmartPointer<vtkDoubleArray> ycoord;
454  assert(spacedim!=3);
455  spacedim=2;
456  int Nx = x.size();
457  int Ny = y.size();
458 
459  if (this->data==NULL)
460  {
461  this->data=vtkSmartPointer<vtkRectilinearGrid>::New();
462  auto rdata=vtkRectilinearGrid::SafeDownCast(this->data);
463  assert(rdata);
464  xcoord = vtkSmartPointer<vtkDoubleArray>::New();
465  ycoord = vtkSmartPointer<vtkDoubleArray>::New();
466  rdata->SetXCoordinates(xcoord);
467  rdata->SetYCoordinates(ycoord);
468  rdata->SetDimensions(Nx, Ny, 1);
469  }
470  else
471  {
472  auto rdata=vtkRectilinearGrid::SafeDownCast(this->data);
473  assert(rdata);
474  xcoord=vtkDoubleArray::SafeDownCast(rdata->GetXCoordinates());
475  ycoord=vtkDoubleArray::SafeDownCast(rdata->GetYCoordinates());
476  xcoord->Initialize();
477  ycoord->Initialize();
478  rdata->SetDimensions(Nx, Ny, 1);
479  }
480 
481  xcoord->SetNumberOfComponents(1);
482  xcoord->SetNumberOfTuples(Nx);
483 
484  ycoord->SetNumberOfComponents(1);
485  ycoord->SetNumberOfTuples(Ny);
486 
487  for (int i=0; i<Nx; i++)
488  xcoord->InsertComponent(i, 0, x[i]*coordinate_scale_factor*coordinate_scale_factor_xyz[0]);
489  for (int i=0; i<Ny; i++)
490  ycoord->InsertComponent(i, 0, y[i]*coordinate_scale_factor*coordinate_scale_factor_xyz[1]);
491 
492  this->data->Modified();
493  }
vtkSmartPointer< vtkDataSet > data
double coordinate_scale_factor_xyz[3]
double coordinate_scale_factor