VTKFIG
0.25.0
Easy VTK based in situ visualization
|
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 ®ions_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< DataSet > | New () |
Static constructor of smart pointer to an empty instance. More... | |
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
|
strong |
Enum describing different possible data types.
|
static |
Static constructor of smart pointer to an empty instance.
|
inline |
Enter data of a simplex grid.
V | Vector class counting from zero with member functions size() and operator[]. std::vector will work. |
IV | Vector class counting from zero with member functions size() and operator[]. std::vector will work. |
dim | Space dimension (2 or 3) |
points | Point 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]) |
cells | Simplex 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]) |
|
inline |
Enter data of a 2D rectilinear grid.
V | Vector class counting from zero with member functions size() and operator[]. std::vector will work. |
x | x-coordinates |
y | y-coordinates |
|
inline |
Enter data of a 3D rectilinear grid.
V | Vector class counting from zero with member functions size() and operator[]. std::vector will work. |
x | x-coordinates |
y | y-coordinates |
z | z-coordinates |
|
inline |
Set cell region data (for multiregion grids)
V | Vector class counting from zero with member functions size() and operator[]. std::vector will work. |
cr | Cell region numbers |
name | Name of function |
|
inline |
Set boundary cell region data (for multiregion grids)
V | Vector class counting from zero with member functions size() and operator[]. std::vector will work. |
cr | Cell region numbers |
name | Name of function |
|
inline |
Set boundary cells for simplex grids.
IV | Vector class counting from zero with member functions size() and operator[]. std::vector will work. |
cr | Cell region numbers |
name | Name of function |
|
inline |
Mask cells which shall be shown.
IV | Vector class counting from zero with member functions size() and operator[]. std::vector will work. |
mask_in | Flag array indicating if cell is in |
name | Name of the mask for later reference |
|
inline |
Mask cells which shall be shown.
IV | Vector class counting from zero with member functions size() and operator[]. std::vector will work. |
cells_in | List of cells which shall be shown |
name | Name of the mask for later reference |
|
inline |
Mask cells which shall be shown.
Assumes that DataSet::SetCellRegions() was called.
IV | Vector class counting from zero with member functions size() and operator[]. std::vector will work. |
regions_omitted | List of regions which shall be omitted |
name | Name of the mask for later reference |
|
inline |
Set data of a scalar function defined on the points of the grid.
V | Vector class counting from zero with member functions size() and operator[]. std::vector will work. |
f | Vector of function values |
name | Name of function |
|
inline |
Set data of a vector function defined on the points on a 2D grid.
V | Vector class counting from zero with member functions size() and operator[]. std::vector will work. |
u | Vector of x component values |
v | Vector of y component values |
name | Name of function |
|
inline |
Set data of a vector function defined on the points on a 2D grid.
V | Vector class counting from zero with member functions size() and operator[]. std::vector will work. |
u | Vector of x component values |
v | Vector of y component values |
w | Vector of z component values |
name | Name of function |
|
inline |
Set data of a vector function defined on the points on a 2D grid.
V | Vector class counting from zero with member functions size() and operator[]. std::vector will work. |
uvw | Vector of row-wise stored component values ( u1 v1, w1, u2,v2,w2 ...) |
name | Name of function |
|
inline |
Set data of a scalar function defined on the cells of the grid.
V | Vector class counting from zero with member functions size() and operator[]. std::vector will work. |
f | Vector of function values |
name | Name of function |
|
inline |
Set data of a scalar function defined on the cells of the grid.
V | Vector class counting from zero with member functions size() and operator[]. std::vector will work. |
f | Vector of function values |
name | Name of function |
|
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.
|
inline |
Set scale factor for x coordinates.
This has to be set before the grid is added to the dataset.
|
inline |
Set scale factor for y coordinates.
This has to be set before the grid is added to the dataset.
|
inline |
Set scale factor for z coordinates.
This has to be set before the grid is added to the dataset.
void vtkfig::DataSet::WriteVTK | ( | const std::string | fname, |
const std::string | filetype | ||
) |
Write dataset to disk in VTK format.
fname | File name |
filetype | File type: "A" for ASCII, "B" for binary |
|
inline |
Request the space dimension of the dataset.
DataType vtkfig::DataSet::GetDataType | ( | ) |
Request the data type of the dataset.
|
inline |
Request vtkDataset which contains all the data.
|
inline |
Request boundary dataset.
|
inline |
Request celllist.