VTKFIG  0.25.0
Easy VTK based in situ visualization
vtkfigDataSet.h
Go to the documentation of this file.
1 
7 #ifndef VTKFIG_DATASET_H
8 #define VTKFIG_DATASET_H
9 
10 #include <memory>
11 #include <map>
12 #include <cmath>
13 
14 #include <vtkSmartPointer.h>
15 #include <vtkUnstructuredGrid.h>
16 #include <vtkRectilinearGrid.h>
17 #include <vtkDoubleArray.h>
18 #include <vtkPointData.h>
19 #include <vtkCellData.h>
20 #include <vtkCellType.h>
21 #include <vtkIdList.h>
22 
23 
24 namespace vtkfig
25 {
26 
36 
37  class DataSet
38  {
39 
40  public:
41 
42 
44  static std::shared_ptr<DataSet> New();
45 
46 
65 
66  template <class V, class IV>
67  void SetSimplexGrid(int dim, const V& points, const IV& cells);
68 
77  template<typename V>
78  void SetRectilinearGrid(const V &x, const V &y);
79 
89  template<typename V>
90  void SetRectilinearGrid(const V &x, const V &y, const V &z);
91 
92 
93 
103  template <class V>
104  void SetCellRegions(const V& cr);
105 
115  template <class V>
116  void SetBoundaryCellRegions(const V& cr);
117 
118 
128  template <class IV>
129  void SetSimplexGridBoundaryCells( const IV& cells);
130 
131 
132 
141  template <class IV>
142  void SetCellMaskByMask(const IV& mask_in, const std::string name);
143 
152  template <class IV>
153  void SetCellMaskByList(const IV& cells_in, const std::string name);
154 
155 
156 
167  template <class IV>
168  void SetCellMaskByRegionsOmitted(const IV& regions_omitted, const std::string name);
169 
170 
180  template <class V>
181  void SetPointScalar(const V&f , const std::string name);
182 
183 
195  template <class V>
196  void SetPointVector(const V&u, const V& v, const std::string name);
197 
211  template <class V>
212  void SetPointVector(const V&u, const V& v, const V& w, const std::string name);
213 
223  template <class V>
224  void SetPointVector(const V&uvw, const int dim, const std::string name);
225 
226 
236  template <class V>
237  void SetCellScalar(const V&f, const std::string name);
238 
239 
249  template <class V>
250  void SetBoundaryCellScalar(const V&f, const std::string name);
251 
259  void SetCoordinateScaleFactor(double factor) { coordinate_scale_factor=factor;}
260 
266  void SetCoordinateScaleFactorX(double factor) { coordinate_scale_factor_xyz[0]=factor;}
267 
273  void SetCoordinateScaleFactorY(double factor) { coordinate_scale_factor_xyz[1]=factor;}
274 
280  void SetCoordinateScaleFactorZ(double factor) { coordinate_scale_factor_xyz[2]=factor;}
281 
282 
289  void WriteVTK(const std::string fname, const std::string filetype);
290 
295  int GetSpaceDimension() { return spacedim;}
296 
300  enum class DataType
301  {
302  NoType=0,
303 
304  RectilinearGrid=1,
305 
306  UnstructuredGrid=2
307 
308  };
309 
315 
320  vtkSmartPointer<vtkDataSet> GetVTKDataSet() { return data;}
321 
326  vtkSmartPointer<vtkDataSet> GetVTKBoundaryDataSet() { return boundary_data;}
327 
328 
333  vtkSmartPointer<vtkIdList> GetCellList(std::string name) {return masks[name];}
334 
335 
336  protected:
337  int spacedim=0;
338 
339  DataSet();
340  ~DataSet(){};
341 
342  private:
343 
344  friend class Figure;
345  vtkSmartPointer<vtkDataSet> data=NULL;
346  vtkSmartPointer<vtkDataSet> boundary_data=NULL;
347  double coordinate_scale_factor=1.0;
348  double coordinate_scale_factor_xyz[3]={1.0,1.0,1.0};
349 
350 
351  template<class DATA, class WRITER>
352  void WriteVTK(vtkSmartPointer<DATA> data, const std::string fname, const std::string filetype);
353 
354  std::map<std::string,vtkSmartPointer<vtkIdList>> masks;
355 
356 
357  };
358 
359 
360  template <class V, class IV>
361  inline
363  const V& points,
364  const IV& cells)
365  {
366  if (this->data==NULL)
367  this->data=vtkSmartPointer<vtkUnstructuredGrid>::New();
368 
369  auto udata=vtkUnstructuredGrid::SafeDownCast(this->data);
370  udata->Reset();
371 
372  this->spacedim=dim;
373  auto gridpoints = vtkSmartPointer<vtkPoints>::New();
374  udata->SetPoints(gridpoints);
375 
376  if (this->spacedim==2)
377  {
378  for (size_t icell=0;icell<cells.size(); icell+=3)
379  {
380  vtkIdType c[3]={cells[icell+0],cells[icell+1],cells[icell+2]};
381  udata->InsertNextCell(VTK_TRIANGLE,3,c);
382  }
383 
384  for (size_t ipoint=0;ipoint<points.size(); ipoint+=2)
385  {
386  gridpoints->InsertNextPoint(
387  points[ipoint+0]*coordinate_scale_factor*coordinate_scale_factor_xyz[0],
388  points[ipoint+1]*coordinate_scale_factor*coordinate_scale_factor_xyz[1],
389  0);
390  }
391  }
392  else
393  {
394  for (size_t icell=0;icell<cells.size(); icell+=4)
395  {
396  vtkIdType c[4]={cells[icell+0],cells[icell+1],cells[icell+2],cells[icell+3]};
397  udata->InsertNextCell(VTK_TETRA,4,c);
398  }
399 
400  for (size_t ipoint=0;ipoint<points.size(); ipoint+=3)
401  {
402  gridpoints->InsertNextPoint(points[ipoint+0]*coordinate_scale_factor*coordinate_scale_factor_xyz[0],
403  points[ipoint+1]*coordinate_scale_factor*coordinate_scale_factor_xyz[1],
404  points[ipoint+2]*coordinate_scale_factor*coordinate_scale_factor_xyz[2]
405  );
406  }
407  }
408  }
409 
410  template <class IV>
411  inline
413  {
414  if (this->boundary_data==NULL)
415  this->boundary_data=vtkSmartPointer<vtkUnstructuredGrid>::New();
416 
417  auto budata=vtkUnstructuredGrid::SafeDownCast(this->boundary_data);
418  budata->Reset();
419 
420  auto bgridpoints = vtkSmartPointer<vtkPoints>::New();
421  budata->SetPoints(bgridpoints);
422 
423  auto udata=vtkUnstructuredGrid::SafeDownCast(this->data);
424  auto pdata=udata->GetPoints();
425 
426  int np=pdata->GetNumberOfPoints();
427  std::vector<int>pmask(np);
428  for (size_t i=0;i<pmask.size(); i++)
429  {
430  pmask[i]=-1;
431  }
432 
433 
434 
435  int ip=0;
436  for (int icell=0;icell<cells.size();icell+=this->spacedim)
437  {
438  for (int id=0;id<this->spacedim;id++)
439  {
440  if (pmask[cells[icell+id]]==-1)
441  {
442  pmask[cells[icell+id]]=ip++;
443  double point[3];
444  pdata->GetPoint(cells[icell+id],point);
445  bgridpoints->InsertNextPoint(point[0],point[1],point[2]);
446  }
447 
448  }
449  if (this->spacedim==2)
450  {
451  vtkIdType c[2]={pmask[cells[icell+0]],pmask[cells[icell+1]]};
452  budata->InsertNextCell(VTK_LINE,2,c);
453  }
454  else
455  {
456  vtkIdType c[3]={pmask[cells[icell+0]],pmask[cells[icell+1]],pmask[cells[icell+2]]};
457  budata->InsertNextCell(VTK_TRIANGLE,3,c);
458  }
459 
460  }
461 
462  }
463 
464 
465 
466 
467  template<typename V>
468  inline
469  void DataSet::SetRectilinearGrid(const V &x, const V &y)
470  {
471  vtkSmartPointer<vtkDoubleArray> xcoord;
472  vtkSmartPointer<vtkDoubleArray> ycoord;
473  assert(spacedim!=3);
474  spacedim=2;
475  int Nx = x.size();
476  int Ny = y.size();
477 
478  if (this->data==NULL)
479  {
480  this->data=vtkSmartPointer<vtkRectilinearGrid>::New();
481  auto rdata=vtkRectilinearGrid::SafeDownCast(this->data);
482  assert(rdata);
483  xcoord = vtkSmartPointer<vtkDoubleArray>::New();
484  ycoord = vtkSmartPointer<vtkDoubleArray>::New();
485  rdata->SetXCoordinates(xcoord);
486  rdata->SetYCoordinates(ycoord);
487  rdata->SetDimensions(Nx, Ny, 1);
488  }
489  else
490  {
491  auto rdata=vtkRectilinearGrid::SafeDownCast(this->data);
492  assert(rdata);
493  xcoord=vtkDoubleArray::SafeDownCast(rdata->GetXCoordinates());
494  ycoord=vtkDoubleArray::SafeDownCast(rdata->GetYCoordinates());
495  xcoord->Initialize();
496  ycoord->Initialize();
497  rdata->SetDimensions(Nx, Ny, 1);
498  }
499 
500  xcoord->SetNumberOfComponents(1);
501  xcoord->SetNumberOfTuples(Nx);
502 
503  ycoord->SetNumberOfComponents(1);
504  ycoord->SetNumberOfTuples(Ny);
505 
506  for (int i=0; i<Nx; i++)
507  xcoord->InsertComponent(i, 0, x[i]*coordinate_scale_factor*coordinate_scale_factor_xyz[0]);
508  for (int i=0; i<Ny; i++)
509  ycoord->InsertComponent(i, 0, y[i]*coordinate_scale_factor*coordinate_scale_factor_xyz[1]);
510 
511  this->data->Modified();
512  }
513 
514 
515  template<typename V>
516  inline
517  void DataSet::SetRectilinearGrid(const V &x, const V &y, const V &z)
518  {
519  vtkSmartPointer<vtkDoubleArray> xcoord;
520  vtkSmartPointer<vtkDoubleArray> ycoord;
521  vtkSmartPointer<vtkDoubleArray> zcoord;
522  assert(spacedim!=2);
523  spacedim=3;
524  int Nx = x.size();
525  int Ny = y.size();
526  int Nz = z.size();
527 
528  if (this->data==NULL)
529  {
530  this->data=vtkSmartPointer<vtkRectilinearGrid>::New();
531  auto rdata=vtkRectilinearGrid::SafeDownCast(this->data);
532  assert(rdata);
533  xcoord = vtkSmartPointer<vtkDoubleArray>::New();
534  ycoord = vtkSmartPointer<vtkDoubleArray>::New();
535  zcoord = vtkSmartPointer<vtkDoubleArray>::New();
536  rdata->SetXCoordinates(xcoord);
537  rdata->SetYCoordinates(ycoord);
538  rdata->SetZCoordinates(zcoord);
539  rdata->SetDimensions(Nx, Ny, Nz );
540  }
541  else
542  {
543  auto rdata=vtkRectilinearGrid::SafeDownCast(this->data);
544  assert(rdata);
545  xcoord=vtkDoubleArray::SafeDownCast(rdata->GetXCoordinates());
546  ycoord=vtkDoubleArray::SafeDownCast(rdata->GetYCoordinates());
547  zcoord=vtkDoubleArray::SafeDownCast(rdata->GetZCoordinates());
548  xcoord->Initialize();
549  ycoord->Initialize();
550  zcoord->Initialize();
551  rdata->SetDimensions(Nx, Ny, Nz );
552  }
553 
554  xcoord->SetNumberOfComponents(1);
555  xcoord->SetNumberOfTuples(Nx);
556 
557  ycoord->SetNumberOfComponents(1);
558  ycoord->SetNumberOfTuples(Ny);
559 
560  zcoord->SetNumberOfComponents(1);
561  zcoord->SetNumberOfTuples(Nz);
562 
563  for (int i=0; i<Nx; i++)
564  xcoord->InsertComponent(i, 0, x[i]*coordinate_scale_factor*coordinate_scale_factor_xyz[0]);
565  for (int i=0; i<Ny; i++)
566  ycoord->InsertComponent(i, 0, y[i]*coordinate_scale_factor*coordinate_scale_factor_xyz[1]);
567  for (int i=0; i<Nz; i++)
568  zcoord->InsertComponent(i, 0, z[i]*coordinate_scale_factor*coordinate_scale_factor_xyz[2]);
569 
570  this->data->Modified();
571  }
572 
573 
574 
575  template <class V>
576  inline
577  void DataSet::SetCellRegions(const V&values)
578  {
579  SetCellScalar(values, "cellregions");
580  }
581 
582  template <class V>
583  inline
584  void DataSet::SetBoundaryCellRegions(const V&values)
585  {
586  SetBoundaryCellScalar(values, "boundarycellregions");
587  }
588 
589 
590  template <class V>
591  inline
592  void DataSet::SetCellScalar(const V&values, const std::string name)
593  {
594  assert(this->data!=NULL);
595  auto ncells=this->data->GetNumberOfCells();
596  assert(ncells==values.size());
597  vtkSmartPointer<vtkDoubleArray>gridvalues;
598 
599  if (this->data->GetPointData()->HasArray(name.c_str()))
600  gridvalues=vtkDoubleArray::SafeDownCast(this->data->GetPointData()->GetAbstractArray(name.c_str()));
601  else
602  {
603  gridvalues=vtkSmartPointer<vtkDoubleArray>::New();
604  gridvalues->SetNumberOfComponents(1);
605  gridvalues->SetNumberOfTuples(ncells);
606  gridvalues->SetName(name.c_str());
607  this->data->GetCellData()->AddArray(gridvalues);
608  }
609 
610 
611  for (int i=0;i<ncells; i++)
612  gridvalues->InsertComponent(i,0,values[i]);
613 
614  gridvalues->Modified();
615  }
616 
617 
618  template <class V>
619  inline
620  void DataSet::SetBoundaryCellScalar(const V&values, const std::string name)
621  {
622  assert(this->boundary_data!=NULL);
623  auto ncells=this->boundary_data->GetNumberOfCells();
624  assert(ncells==values.size());
625  vtkSmartPointer<vtkDoubleArray>gridvalues;
626 
627  if (this->boundary_data->GetPointData()->HasArray(name.c_str()))
628  gridvalues=vtkDoubleArray::SafeDownCast(this->boundary_data->GetPointData()->GetAbstractArray(name.c_str()));
629  else
630  {
631  gridvalues=vtkSmartPointer<vtkDoubleArray>::New();
632  gridvalues->SetNumberOfComponents(1);
633  gridvalues->SetNumberOfTuples(ncells);
634  gridvalues->SetName(name.c_str());
635  this->boundary_data->GetCellData()->AddArray(gridvalues);
636  }
637 
638 
639  for (int i=0;i<ncells; i++)
640  gridvalues->InsertComponent(i,0,values[i]);
641 
642  gridvalues->Modified();
643  }
644 
645 
646 
647  template <class IV>
648  inline
649  void DataSet::SetCellMaskByRegionsOmitted(const IV& regions_omitted, const std::string name)
650  {
651  auto celllist=vtkSmartPointer<vtkIdList>::New();
652  auto cellregions=vtkDoubleArray::SafeDownCast(this->data->GetCellData()->GetAbstractArray("cellregions"));
653  assert(cellregions);
654  int ncells=this->data->GetNumberOfCells();
655  int icelllist=0;
656  for (int icell=0;icell<ncells;icell++)
657  {
658  int ireg=cellregions->GetComponent(icell,0);
659  bool use_cell=true;
660 
661  for (int io=0;io<regions_omitted.size();io++)
662  if (ireg==regions_omitted[io])
663  {
664  use_cell=false;
665  continue;
666  }
667  if (use_cell)
668  celllist->InsertId(icelllist++,icell);
669  }
670  masks[name]=celllist;
671  }
672 
673 
674 
675  template <class IV>
676  inline
677  void DataSet::SetCellMaskByMask(const IV& cells_in, const std::string name)
678  {
679  assert(cells_in.size()==this->data->GetNumberOfCells());
680  auto celllist=vtkSmartPointer<vtkIdList>::New();
681  int j=0;
682  for (int i=0;i<cells_in.size();i++)
683  if (cells_in[i]) celllist->InsertId(j++,i);
684  masks[name]=celllist;
685  }
686 
687  template <class IV>
688  inline
689  void DataSet::SetCellMaskByList(const IV& cells_in, const std::string name)
690  {
691  auto celllist=vtkSmartPointer<vtkIdList>::New();
692  int j=0;
693  for (int i=0;i<cells_in.size();i++)
694  celllist->InsertId(j++,cells_in[i]);
695  masks[name]=celllist;
696  }
697 
698 
699 
700 
701  template <class V>
702  inline
703  void DataSet::SetPointScalar(const V&values, const std::string name)
704  {
705  assert(this->data!=NULL);
706  size_t npoints=this->data->GetNumberOfPoints();
707  assert(npoints==values.size());
708  vtkSmartPointer<vtkDoubleArray>gridvalues;
709 
710  if (this->data->GetPointData()->HasArray(name.c_str()))
711  gridvalues=vtkDoubleArray::SafeDownCast(this->data->GetPointData()->GetAbstractArray(name.c_str()));
712  else
713  {
714  gridvalues=vtkSmartPointer<vtkDoubleArray>::New();
715  gridvalues->SetNumberOfComponents(1);
716  gridvalues->SetNumberOfTuples(npoints);
717  gridvalues->SetName(name.c_str());
718  this->data->GetPointData()->AddArray(gridvalues);
719  }
720 
721  for (size_t i=0;i<npoints; i++)
722  gridvalues->InsertComponent(i,0,values[i]);
723  gridvalues->Modified();
724  }
725 
726  template <class V>
727  inline
728  void DataSet::SetPointVector(const V&u, const V& v, const std::string name)
729  {
730  assert(this->spacedim==2);
731  assert(this->data!=NULL);
732  size_t npoints=this->data->GetNumberOfPoints();
733  assert(npoints==u.size());
734  assert(npoints==v.size());
735  vtkSmartPointer<vtkDoubleArray>gridvalues;
736 
737  if (this->data->GetPointData()->HasArray(name.c_str()))
738  gridvalues=vtkDoubleArray::SafeDownCast(this->data->GetPointData()->GetAbstractArray(name.c_str()));
739  else
740  {
741  gridvalues=vtkSmartPointer<vtkDoubleArray>::New();
742  gridvalues->SetNumberOfComponents(3);
743  gridvalues->SetNumberOfTuples(npoints);
744  gridvalues->SetName(name.c_str());
745  this->data->GetPointData()->AddArray(gridvalues);
746  }
747 
748  for (size_t i=0;i<npoints; i++)
749  gridvalues->InsertTuple3(i,u[i],v[i],0);
750 
751 
752  gridvalues->Modified();
753  }
754 
755 
756 
757  template <class V>
758  inline
759  void DataSet::SetPointVector(const V&u, const V& v, const V& w, const std::string name)
760  {
761  assert(this->spacedim==3);
762  assert(this->data!=NULL);
763  size_t npoints=this->data->GetNumberOfPoints();
764  assert(npoints==u.size());
765  assert(npoints==v.size());
766  assert(npoints==w.size());
767  vtkSmartPointer<vtkDoubleArray>gridvalues;
768 
769  if (this->data->GetPointData()->HasArray(name.c_str()))
770  gridvalues=vtkDoubleArray::SafeDownCast(this->data->GetPointData()->GetAbstractArray(name.c_str()));
771  else
772  {
773  gridvalues=vtkSmartPointer<vtkDoubleArray>::New();
774  gridvalues->SetNumberOfComponents(3);
775  gridvalues->SetNumberOfTuples(npoints);
776  gridvalues->SetName(name.c_str());
777  this->data->GetPointData()->AddArray(gridvalues);
778  }
779 
780 
781  for (size_t i=0;i<npoints; i++)
782  gridvalues->InsertTuple3(i,u[i],v[i],w[i]);
783  gridvalues->Modified();
784  }
785 
786 
787  template <class V>
788  inline
789  void DataSet::SetPointVector(const V&uvw, int dim, const std::string name)
790  {
791  assert(this->spacedim==dim);
792  assert(this->data!=NULL);
793  size_t npoints=this->data->GetNumberOfPoints();
794  assert(npoints==uvw.size()/dim);
795 
796  vtkSmartPointer<vtkDoubleArray>gridvalues;
797 
798  if (this->data->GetPointData()->HasArray(name.c_str()))
799  {
800  gridvalues=vtkDoubleArray::SafeDownCast(this->data->GetPointData()->GetAbstractArray(name.c_str()));
801  }
802  else
803  {
804  gridvalues=vtkSmartPointer<vtkDoubleArray>::New();
805  gridvalues->SetNumberOfComponents(3);
806  gridvalues->SetNumberOfTuples(npoints);
807  gridvalues->SetName(name.c_str());
808  this->data->GetPointData()->AddArray(gridvalues);
809  }
810 
811  switch(dim)
812  {
813  case 2:
814  for (size_t i=0,j=0;i<npoints; i++,j+=2)
815  gridvalues->InsertTuple3(i,uvw[j],uvw[j+1],0);
816  break;
817  case 3:
818  for (size_t i=0,j=0;i<npoints; i++,j+=3)
819  gridvalues->InsertTuple3(i,uvw[j],uvw[j+1],uvw[j+2]);
820  break;
821  default: break;
822  }
823 
824  gridvalues->Modified();
825  }
826 
827 }
828 
829 
830 #endif
vtkfig::DataSet::GetVTKDataSet
vtkSmartPointer< vtkDataSet > GetVTKDataSet()
Request vtkDataset which contains all the data.
Definition: vtkfigDataSet.h:320
vtkfig::DataSet::WriteVTK
void WriteVTK(const std::string fname, const std::string filetype)
Write dataset to disk in VTK format.
vtkfig::DataSet::SetCoordinateScaleFactorZ
void SetCoordinateScaleFactorZ(double factor)
Set scale factor for z coordinates.
Definition: vtkfigDataSet.h:280
vtkfig::DataSet::SetCoordinateScaleFactorX
void SetCoordinateScaleFactorX(double factor)
Set scale factor for x coordinates.
Definition: vtkfigDataSet.h:266
vtkfig::DataSet::SetRectilinearGrid
void SetRectilinearGrid(const V &x, const V &y)
Enter data of a 2D rectilinear grid.
Definition: vtkfigDataSet.h:469
vtkfig::DataSet::DataType
DataType
Enum describing different possible data types.
Definition: vtkfigDataSet.h:300
vtkfig::DataSet::SetCoordinateScaleFactorY
void SetCoordinateScaleFactorY(double factor)
Set scale factor for y coordinates.
Definition: vtkfigDataSet.h:273
vtkfig::DataSet::SetCellRegions
void SetCellRegions(const V &cr)
Set cell region data (for multiregion grids)
Definition: vtkfigDataSet.h:577
vtkfig::DataSet::SetPointScalar
void SetPointScalar(const V &f, const std::string name)
Set data of a scalar function defined on the points of the grid.
Definition: vtkfigDataSet.h:703
vtkfig::DataSet::SetBoundaryCellScalar
void SetBoundaryCellScalar(const V &f, const std::string name)
Set data of a scalar function defined on the cells of the grid.
Definition: vtkfigDataSet.h:620
vtkfig::DataSet::SetCellMaskByRegionsOmitted
void SetCellMaskByRegionsOmitted(const IV &regions_omitted, const std::string name)
Mask cells which shall be shown.
Definition: vtkfigDataSet.h:649
vtkfig::DataSet::SetSimplexGridBoundaryCells
void SetSimplexGridBoundaryCells(const IV &cells)
Set boundary cells for simplex grids.
Definition: vtkfigDataSet.h:412
vtkfig::DataSet::SetCellMaskByList
void SetCellMaskByList(const IV &cells_in, const std::string name)
Mask cells which shall be shown.
Definition: vtkfigDataSet.h:689
vtkfig::DataSet::SetSimplexGrid
void SetSimplexGrid(int dim, const V &points, const IV &cells)
Enter data of a simplex grid.
Definition: vtkfigDataSet.h:362
vtkfig::DataSet::SetPointVector
void SetPointVector(const V &u, const V &v, const std::string name)
Set data of a vector function defined on the points on a 2D grid.
Definition: vtkfigDataSet.h:728
vtkfig::DataSet
Class to collect all data given on one grid.
Definition: vtkfigDataSet.h:37
vtkfig::DataSet::SetCellMaskByMask
void SetCellMaskByMask(const IV &mask_in, const std::string name)
Mask cells which shall be shown.
Definition: vtkfigDataSet.h:677
vtkfig::DataSet::SetCoordinateScaleFactor
void SetCoordinateScaleFactor(double factor)
Set scale factor for coordinates.
Definition: vtkfigDataSet.h:259
vtkfig::DataSet::GetDataType
DataType GetDataType()
Request the data type of the dataset.
vtkfig::DataSet::GetVTKBoundaryDataSet
vtkSmartPointer< vtkDataSet > GetVTKBoundaryDataSet()
Request boundary dataset.
Definition: vtkfigDataSet.h:326
vtkfig::DataSet::SetBoundaryCellRegions
void SetBoundaryCellRegions(const V &cr)
Set boundary cell region data (for multiregion grids)
Definition: vtkfigDataSet.h:584
vtkfig::DataSet::SetCellScalar
void SetCellScalar(const V &f, const std::string name)
Set data of a scalar function defined on the cells of the grid.
Definition: vtkfigDataSet.h:592
vtkfig::DataSet::New
static std::shared_ptr< DataSet > New()
Static constructor of smart pointer to an empty instance.
vtkfig::DataSet::GetCellList
vtkSmartPointer< vtkIdList > GetCellList(std::string name)
Request celllist.
Definition: vtkfigDataSet.h:333
vtkfig::DataSet::GetSpaceDimension
int GetSpaceDimension()
Request the space dimension of the dataset.
Definition: vtkfigDataSet.h:295