A collection of static methods to assist converting OpenFOAM data structures into VTK internal data structures. More...
Classes | |
struct | Patch |
Convert OpenFOAM patch to vtkPolyData. More... | |
Functions | |
UList< uint8_t > | asUList (vtkUnsignedCharArray *array, const label size) |
Wrap vtkUnsignedCharArray as a UList. More... | |
UList< vtkIdType > | asUList (vtkIdTypeArray *array, const label size) |
Wrap vtkIdTypeArray as a UList. More... | |
vtkSmartPointer< vtkPoints > | Points (const UList< point > &pts) |
Return a list of points as vtkPoints. More... | |
vtkSmartPointer< vtkPoints > | Points (const UList< point > &pts, const labelUList &addr) |
Return an indirect list of points as vtkPoints. More... | |
template<class Face > | |
vtkSmartPointer< vtkCellArray > | Faces (const UList< Face > &faces) |
Convert a list of faces (or triFaces) to vtk polygon cells. More... | |
vtkSmartPointer< vtkPolyData > | Vertices (const UList< point > &pts) |
Return vtkPolyData of vertices for each point. More... | |
vtkSmartPointer< vtkPolyData > | Vertices (const UList< point > &pts, const labelUList &addr) |
Return vtkPolyData of vertices for each point. More... | |
scalarMinMax | rangeOf (vtkDataArray *data) |
Min/Max of scalar, or mag() of non-scalars. Includes nullptr check. More... | |
template<class Type > | |
void | remapTuple (float data[]) |
Remapping for some OpenFOAM data types (eg, symmTensor) More... | |
template<> | |
void | remapTuple< symmTensor > (float data[]) |
Template specialization for symmTensor ordering. More... | |
template<class Type > | |
void | remapTuple (double data[]) |
Remapping for some OpenFOAM data types (eg, symmTensor) More... | |
template<> | |
void | remapTuple< symmTensor > (double data[]) |
Template specialization for symmTensor ordering. More... | |
template<class Type > | |
void | foamToVtkTuple (float output[], const Type &val) |
Copy/transcribe OpenFOAM data types to VTK format. More... | |
template<class Type > | |
void | foamToVtkTuple (double output[], const Type &val) |
Copy/transcribe OpenFOAM data types to VTK format. More... | |
template<class Type > | |
label | transcribeFloatData (vtkFloatArray *array, const UList< Type > &input, vtkIdType start=0) |
Copy list to pre-allocated vtk array. More... | |
template<class Type > | |
vtkSmartPointer< vtkFloatArray > | zeroField (const word &name, const label size) |
Create named field initialized to zero. More... | |
template<class Type > | |
vtkSmartPointer< vtkFloatArray > | convertFieldToVTK (const word &name, const UList< Type > &fld) |
Convert field data to a vtkFloatArray. More... | |
vtkSmartPointer< vtkCellArray > | identityVertices (const label size) |
An identity list of VTK_VERTEX. More... | |
A collection of static methods to assist converting OpenFOAM data structures into VTK internal data structures.
Remapping of the symmTensor order is required in input or output directions. OpenFOAM uses (XX, XY, XZ, YY, YZ, ZZ) order, VTK uses (XX, YY, ZZ, XY, YZ, XZ) order.
|
inline |
Wrap vtkUnsignedCharArray as a UList.
Definition at line 32 of file foamVtkToolsI.H.
Referenced by vtuAdaptor::internal().
|
inline |
Wrap vtkIdTypeArray as a UList.
Definition at line 45 of file foamVtkToolsI.H.
Return a list of points as vtkPoints.
Definition at line 59 of file foamVtkToolsI.H.
References UList< T >::size().
Referenced by Patch::mesh(), and Patch::points().
|
inline |
Return an indirect list of points as vtkPoints.
Definition at line 76 of file foamVtkToolsI.H.
References UList< T >::size().
vtkSmartPointer< vtkCellArray > Faces | ( | const UList< Face > & | faces | ) |
Convert a list of faces (or triFaces) to vtk polygon cells.
Definition at line 44 of file foamVtkToolsTemplates.C.
References cells, f(), and UList< T >::size().
Referenced by Patch::faces(), and Patch::mesh().
Return vtkPolyData of vertices for each point.
Definition at line 93 of file foamVtkToolsI.H.
References UList< T >::size().
|
inline |
Return vtkPolyData of vertices for each point.
Definition at line 105 of file foamVtkToolsI.H.
References UList< T >::size().
|
inline |
Min/Max of scalar, or mag() of non-scalars. Includes nullptr check.
Definition at line 116 of file foamVtkToolsI.H.
References range.
|
inline |
Remapping for some OpenFOAM data types (eg, symmTensor)
data[in,out] | The data to be remapped in-place |
Definition at line 256 of file foamVtkTools.H.
|
inline |
Template specialization for symmTensor ordering.
|
inline |
Remapping for some OpenFOAM data types (eg, symmTensor)
data[in,out] | The data to be remapped in-place |
Definition at line 265 of file foamVtkTools.H.
|
inline |
Template specialization for symmTensor ordering.
|
inline |
Copy/transcribe OpenFOAM data types to VTK format.
This allows a change of data type (float vs double) as well as addressing any swapping issues (eg, symmTensor)
output[out] | The output scratch space. Must be long enough to hold the result. |
val[in] | The input data to copy/transcribe |
Definition at line 213 of file foamVtkToolsI.H.
References Foam::component().
Referenced by vtuAdaptor::convertField().
|
inline |
Copy/transcribe OpenFOAM data types to VTK format.
This allows a change of data type (float vs double) as well as addressing any swapping issues (eg, symmTensor)
output[out] | The output scratch space. Must be long enough to hold the result. |
val[in] | The input data to copy/transcribe |
Definition at line 228 of file foamVtkToolsI.H.
References Foam::component().
label transcribeFloatData | ( | vtkFloatArray * | array, |
const UList< Type > & | input, | ||
vtkIdType | start = 0 |
||
) |
Copy list to pre-allocated vtk array.
start | The write offset into output array |
Referenced by convertFieldToVTK().
Create named field initialized to zero.
Definition at line 326 of file foamVtkToolsTemplates.C.
References Foam::name().
Convert field data to a vtkFloatArray.
Definition at line 354 of file foamVtkToolsTemplates.C.
References fld(), Foam::name(), and transcribeFloatData().
|
inline |