VTKFIG  0.20.20181114
Easy VTK based in situ visualization
template<class V , class IV >
void vtkfig::DataSet::SetSimplexGrid ( int  dim,
const V &  points,
const IV &  cells 
)
inline

Enter data of a simplex grid.

Template Parameters
VVector class counting from zero with member functions size() and operator[]. std::vector will work.
IVVector class counting from zero with member functions size() and operator[]. std::vector will work.
Parameters
dimSpace dimension (2 or 3)
pointsPoint coordinates, stored consecutively.
In 2D, coordinates of point i are (points[2*i],points[2*i+1])
In 3D, coordinates of point i are (points[3*i],points[3*i+1],,points[3*i+2])
cellsSimplex point indices (counting from 0!) In 2D, point indices of triangle i are (cells[3*i],cells[3*i+1],cells[3*i+2])
In 3D, point indices of tetrahedron i are (cells[3*i],cells[3*i+1],cells[3*i+2],cells[3*i+3])

Definition at line 343 of file vtkfigDataSet.h.

346  {
347  if (this->data==NULL)
348  this->data=vtkSmartPointer<vtkUnstructuredGrid>::New();
349 
350  auto udata=vtkUnstructuredGrid::SafeDownCast(this->data);
351  udata->Reset();
352 
353  this->spacedim=dim;
354  auto gridpoints = vtkSmartPointer<vtkPoints>::New();
355  udata->SetPoints(gridpoints);
356 
357  if (this->spacedim==2)
358  {
359  for (int icell=0;icell<cells.size(); icell+=3)
360  {
361  vtkIdType c[3]={cells[icell+0],cells[icell+1],cells[icell+2]};
362  udata->InsertNextCell(VTK_TRIANGLE,3,c);
363  }
364 
365  for (int ipoint=0;ipoint<points.size(); ipoint+=2)
366  {
367  gridpoints->InsertNextPoint(
369  points[ipoint+1]*coordinate_scale_factor*coordinate_scale_factor_xyz[1],
370  0);
371  }
372  }
373  else
374  {
375  for (int icell=0;icell<cells.size(); icell+=4)
376  {
377  vtkIdType c[4]={cells[icell+0],cells[icell+1],cells[icell+2],cells[icell+3]};
378  udata->InsertNextCell(VTK_TETRA,4,c);
379  }
380 
381  for (int ipoint=0;ipoint<points.size(); ipoint+=3)
382  {
383  gridpoints->InsertNextPoint(points[ipoint+0]*coordinate_scale_factor*coordinate_scale_factor_xyz[0],
384  points[ipoint+1]*coordinate_scale_factor*coordinate_scale_factor_xyz[1],
385  points[ipoint+2]*coordinate_scale_factor*coordinate_scale_factor_xyz[2]
386  );
387  }
388  }
389  }
vtkSmartPointer< vtkDataSet > data
double coordinate_scale_factor_xyz[3]
double coordinate_scale_factor