VTKFIG  0.25.0
Easy VTK based in situ visualization
Public Types | Public Member Functions | Static Public Member Functions | List of all members
vtkfig::DataSet Class Reference

Class to collect all data given on one grid. More...

Public Types

enum  DataType
 Enum describing different possible data types. More...
 

Public Member Functions

template<class V , class IV >
void SetSimplexGrid (int dim, const V &points, const IV &cells)
 Enter data of a simplex grid. More...
 
template<typename V >
void SetRectilinearGrid (const V &x, const V &y)
 Enter data of a 2D rectilinear grid. More...
 
template<typename V >
void SetRectilinearGrid (const V &x, const V &y, const V &z)
 Enter data of a 3D rectilinear grid. More...
 
template<class V >
void SetCellRegions (const V &cr)
 Set cell region data (for multiregion grids) More...
 
template<class V >
void SetBoundaryCellRegions (const V &cr)
 Set boundary cell region data (for multiregion grids) More...
 
template<class IV >
void SetSimplexGridBoundaryCells (const IV &cells)
 Set boundary cells for simplex grids. More...
 
template<class IV >
void SetCellMaskByMask (const IV &mask_in, const std::string name)
 Mask cells which shall be shown. More...
 
template<class IV >
void SetCellMaskByList (const IV &cells_in, const std::string name)
 Mask cells which shall be shown. More...
 
template<class IV >
void SetCellMaskByRegionsOmitted (const IV &regions_omitted, const std::string name)
 Mask cells which shall be shown. More...
 
template<class V >
void SetPointScalar (const V &f, const std::string name)
 Set data of a scalar function defined on the points of the grid. More...
 
template<class V >
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. More...
 
template<class V >
void SetPointVector (const V &u, const V &v, const V &w, const std::string name)
 Set data of a vector function defined on the points on a 2D grid. More...
 
template<class V >
void SetPointVector (const V &uvw, const int dim, const std::string name)
 Set data of a vector function defined on the points on a 2D grid. More...
 
template<class V >
void SetCellScalar (const V &f, const std::string name)
 Set data of a scalar function defined on the cells of the grid. More...
 
template<class V >
void SetBoundaryCellScalar (const V &f, const std::string name)
 Set data of a scalar function defined on the cells of the grid. More...
 
void SetCoordinateScaleFactor (double factor)
 Set scale factor for coordinates. More...
 
void SetCoordinateScaleFactorX (double factor)
 Set scale factor for x coordinates. More...
 
void SetCoordinateScaleFactorY (double factor)
 Set scale factor for y coordinates. More...
 
void SetCoordinateScaleFactorZ (double factor)
 Set scale factor for z coordinates. More...
 
void WriteVTK (const std::string fname, const std::string filetype)
 Write dataset to disk in VTK format. More...
 
int GetSpaceDimension ()
 Request the space dimension of the dataset. More...
 
DataType GetDataType ()
 Request the data type of the dataset. More...
 
vtkSmartPointer< vtkDataSet > GetVTKDataSet ()
 Request vtkDataset which contains all the data. More...
 
vtkSmartPointer< vtkDataSet > GetVTKBoundaryDataSet ()
 Request boundary dataset. More...
 
vtkSmartPointer< vtkIdList > GetCellList (std::string name)
 Request celllist. More...
 

Static Public Member Functions

static std::shared_ptr< DataSetNew ()
 Static constructor of smart pointer to an empty instance. More...
 

Detailed Description

Class to collect all data given on one grid.

Wrapper of vtkRectilinearGrid resp. vtkUnstructuredGrid Using duck typing, this class internally builds copies of the data in the form of a vtkDataSet

Member Enumeration Documentation

◆ DataType

Enum describing different possible data types.

Member Function Documentation

◆ New()

static std::shared_ptr<DataSet> vtkfig::DataSet::New ( )
static

◆ SetSimplexGrid()

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])
Examples
example-gridview2d.cxx, and example-gridview3d.cxx.

◆ SetRectilinearGrid() [1/2]

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
Examples
example-gridview2d.cxx, example-gridview3d.cxx, example-multifig.cxx, and example-multiframe.cxx.

◆ SetRectilinearGrid() [2/2]

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

◆ SetCellRegions()

template<class V >
void vtkfig::DataSet::SetCellRegions ( const V &  cr)
inline

Set cell region data (for multiregion grids)

Template Parameters
VVector class counting from zero with member functions size() and operator[]. std::vector will work.
Parameters
crCell region numbers
nameName of function
+ Here is the call graph for this function:

◆ SetBoundaryCellRegions()

template<class V >
void vtkfig::DataSet::SetBoundaryCellRegions ( const V &  cr)
inline

Set boundary cell region data (for multiregion grids)

Template Parameters
VVector class counting from zero with member functions size() and operator[]. std::vector will work.
Parameters
crCell region numbers
nameName of function
+ Here is the call graph for this function:

◆ SetSimplexGridBoundaryCells()

template<class IV >
void vtkfig::DataSet::SetSimplexGridBoundaryCells ( const IV &  cells)
inline

Set boundary cells for simplex grids.

Template Parameters
IVVector class counting from zero with member functions size() and operator[]. std::vector will work.
Parameters
crCell region numbers
nameName of function

◆ SetCellMaskByMask()

template<class IV >
void vtkfig::DataSet::SetCellMaskByMask ( const IV &  mask_in,
const std::string  name 
)
inline

Mask cells which shall be shown.

Template Parameters
IVVector class counting from zero with member functions size() and operator[]. std::vector will work.
Parameters
mask_inFlag array indicating if cell is in
nameName of the mask for later reference

◆ SetCellMaskByList()

template<class IV >
void vtkfig::DataSet::SetCellMaskByList ( const IV &  cells_in,
const std::string  name 
)
inline

Mask cells which shall be shown.

Template Parameters
IVVector class counting from zero with member functions size() and operator[]. std::vector will work.
Parameters
cells_inList of cells which shall be shown
nameName of the mask for later reference

◆ SetCellMaskByRegionsOmitted()

template<class IV >
void vtkfig::DataSet::SetCellMaskByRegionsOmitted ( const IV &  regions_omitted,
const std::string  name 
)
inline

Mask cells which shall be shown.

Assumes that DataSet::SetCellRegions() was called.

Template Parameters
IVVector class counting from zero with member functions size() and operator[]. std::vector will work.
Parameters
regions_omittedList of regions which shall be omitted
nameName of the mask for later reference

◆ SetPointScalar()

template<class V >
void vtkfig::DataSet::SetPointScalar ( const V &  f,
const std::string  name 
)
inline

Set data of a scalar function defined on the points of the grid.

Template Parameters
VVector class counting from zero with member functions size() and operator[]. std::vector will work.
Parameters
fVector of function values
nameName of function
Examples
example-multifig.cxx, and example-multiframe.cxx.

◆ SetPointVector() [1/3]

template<class V >
void vtkfig::DataSet::SetPointVector ( const V &  u,
const V &  v,
const std::string  name 
)
inline

Set data of a vector function defined on the points on a 2D grid.

Template Parameters
VVector class counting from zero with member functions size() and operator[]. std::vector will work.
Parameters
uVector of x component values
vVector of y component values
nameName of function

◆ SetPointVector() [2/3]

template<class V >
void vtkfig::DataSet::SetPointVector ( const V &  u,
const V &  v,
const V &  w,
const std::string  name 
)
inline

Set data of a vector function defined on the points on a 2D grid.

Template Parameters
VVector class counting from zero with member functions size() and operator[]. std::vector will work.
Parameters
uVector of x component values
vVector of y component values
wVector of z component values
nameName of function

◆ SetPointVector() [3/3]

template<class V >
void vtkfig::DataSet::SetPointVector ( const V &  uvw,
const int  dim,
const std::string  name 
)
inline

Set data of a vector function defined on the points on a 2D grid.

Template Parameters
VVector class counting from zero with member functions size() and operator[]. std::vector will work.
Parameters
uvwVector of row-wise stored component values ( u1 v1, w1, u2,v2,w2 ...)
nameName of function

◆ SetCellScalar()

template<class V >
void vtkfig::DataSet::SetCellScalar ( const V &  f,
const std::string  name 
)
inline

Set data of a scalar function defined on the cells of the grid.

Template Parameters
VVector class counting from zero with member functions size() and operator[]. std::vector will work.
Parameters
fVector of function values
nameName of function

◆ SetBoundaryCellScalar()

template<class V >
void vtkfig::DataSet::SetBoundaryCellScalar ( const V &  f,
const std::string  name 
)
inline

Set data of a scalar function defined on the cells of the grid.

Template Parameters
VVector class counting from zero with member functions size() and operator[]. std::vector will work.
Parameters
fVector of function values
nameName of function

◆ SetCoordinateScaleFactor()

void vtkfig::DataSet::SetCoordinateScaleFactor ( double  factor)
inline

Set scale factor for coordinates.

This has to be set before the grid is added to the dataset. Some vtk algorithms (notably vtkProbeFilter) struggle with small scale grids, and so one can try to increase the magnitude of coordinates.

◆ SetCoordinateScaleFactorX()

void vtkfig::DataSet::SetCoordinateScaleFactorX ( double  factor)
inline

Set scale factor for x coordinates.

This has to be set before the grid is added to the dataset.

◆ SetCoordinateScaleFactorY()

void vtkfig::DataSet::SetCoordinateScaleFactorY ( double  factor)
inline

Set scale factor for y coordinates.

This has to be set before the grid is added to the dataset.

◆ SetCoordinateScaleFactorZ()

void vtkfig::DataSet::SetCoordinateScaleFactorZ ( double  factor)
inline

Set scale factor for z coordinates.

This has to be set before the grid is added to the dataset.

◆ WriteVTK()

void vtkfig::DataSet::WriteVTK ( const std::string  fname,
const std::string  filetype 
)

Write dataset to disk in VTK format.

Parameters
fnameFile name
filetypeFile type: "A" for ASCII, "B" for binary

◆ GetSpaceDimension()

int vtkfig::DataSet::GetSpaceDimension ( )
inline

Request the space dimension of the dataset.

Returns
Space dimension

◆ GetDataType()

DataType vtkfig::DataSet::GetDataType ( )

Request the data type of the dataset.

Returns
Data type

◆ GetVTKDataSet()

vtkSmartPointer<vtkDataSet> vtkfig::DataSet::GetVTKDataSet ( )
inline

Request vtkDataset which contains all the data.

Returns
vtkDataset

◆ GetVTKBoundaryDataSet()

vtkSmartPointer<vtkDataSet> vtkfig::DataSet::GetVTKBoundaryDataSet ( )
inline

Request boundary dataset.

Returns
vtkDataset

◆ GetCellList()

vtkSmartPointer<vtkIdList> vtkfig::DataSet::GetCellList ( std::string  name)
inline

Request celllist.

Returns
Cell List

The documentation for this class was generated from the following file: