Encapsulate output methods for vtk.
More...
#include <VTKPolyDataWriterInterface.hpp>
|
bool | isFieldPlotEnabled (dataRepository::WrapperBase const &wrapper) const |
| Check if plotting is enabled for this field. More...
|
|
void | writeCellElementRegions (real64 time, ElementRegionManager const &elemManager, NodeManager const &nodeManager, string const &path) |
| Writes the files for all the CellElementRegions. More...
|
|
void | writeWellElementRegions (real64 time, ElementRegionManager const &elemManager, NodeManager const &nodeManager, string const &path) |
| Writes the files containing the well representation. More...
|
|
void | writeParticleRegions (real64 const time, ParticleManager const &particleManager, string const &path) |
| Writes the files containing the particle representation. More...
|
|
void | writeSurfaceElementRegions (real64 time, ElementRegionManager const &elemManager, NodeManager const &nodeManager, EmbeddedSurfaceNodeManager const &embSurfNodeManager, FaceManager const &faceManager, string const &path) |
| Writes the files containing the faces elements. More...
|
|
void | writeVtmFile (integer const cycle, DomainPartition const &domain, VTKVTMWriter const &vtmWriter) const |
| Writes a VTM file for the time-step time . More...
|
|
void | writeNodeFields (NodeManager const &nodeManager, arrayView1d< localIndex const > const &nodeIndices, vtkPointData *pointData) const |
| Write all the fields associated to the nodes of nodeManager if their plotlevel is <= m_plotLevel. More...
|
|
void | writeElementFields (ElementRegionBase const &subRegion, vtkCellData *cellData) const |
| Writes all the fields associated to the elements of er if their plotlevel is <= m_plotLevel. More...
|
|
void | writeParticleFields (ParticleRegionBase const ®ion, vtkCellData *cellData) const |
| Writes all the fields associated to the elements of region if their plotlevel is <= m_plotLevel. More...
|
|
void | writeUnstructuredGrid (string const &path, ObjectManagerBase const ®ion, vtkUnstructuredGrid *ug) |
| Writes an unstructured grid. More...
|
|
Encapsulate output methods for vtk.
Definition at line 77 of file VTKPolyDataWriterInterface.hpp.
◆ VTKPolyDataWriterInterface()
geos::vtk::VTKPolyDataWriterInterface::VTKPolyDataWriterInterface |
( |
string |
outputName | ) |
|
|
explicit |
Constructor.
- Parameters
-
[in] | outputName | folder name in which all the files will be written |
◆ isFieldPlotEnabled()
Check if plotting is enabled for this field.
- Parameters
-
- Returns
- true if this wrapper should be plot, false otherwise
◆ setFieldNames()
void geos::vtk::VTKPolyDataWriterInterface::setFieldNames |
( |
string_array const & |
fieldNames | ) |
|
|
inline |
◆ setLevelNames()
void geos::vtk::VTKPolyDataWriterInterface::setLevelNames |
( |
string_array const & |
levelNames | ) |
|
|
inline |
Set the names of the mesh levels to output.
- Parameters
-
[in] | levelNames | the mesh levels to output (an empty array means all levels are saved) |
Definition at line 167 of file VTKPolyDataWriterInterface.hpp.
◆ setNumberOfTargetProcesses()
void geos::vtk::VTKPolyDataWriterInterface::setNumberOfTargetProcesses |
( |
integer const |
numberOfTargetProcesses | ) |
|
|
inline |
◆ setOnlyPlotSpecifiedFieldNamesFlag()
void geos::vtk::VTKPolyDataWriterInterface::setOnlyPlotSpecifiedFieldNamesFlag |
( |
integer const |
onlyPlotSpecifiedFieldNames | ) |
|
|
inline |
Set the flag to decide whether we only plot the fields specified by fieldNames, or if we also plot fields based on plotLevel.
- Parameters
-
[in] | onlyPlotSpecifiedFieldNames | the flag |
Definition at line 149 of file VTKPolyDataWriterInterface.hpp.
◆ setOutputLocation()
void geos::vtk::VTKPolyDataWriterInterface::setOutputLocation |
( |
string |
outputDir, |
|
|
string |
outputName |
|
) |
| |
|
inline |
Set the output directory name.
- Parameters
-
[in] | outputDir | global output directory location |
[in] | outputName | name of the VTK output subdirectory and corresponding PVD file |
Definition at line 138 of file VTKPolyDataWriterInterface.hpp.
◆ setOutputMode()
void geos::vtk::VTKPolyDataWriterInterface::setOutputMode |
( |
VTKOutputMode |
mode | ) |
|
|
inline |
◆ setOutputRegionType()
void geos::vtk::VTKPolyDataWriterInterface::setOutputRegionType |
( |
VTKRegionTypes |
regionType | ) |
|
|
inline |
◆ setPlotLevel()
void geos::vtk::VTKPolyDataWriterInterface::setPlotLevel |
( |
integer |
plotLevel | ) |
|
|
inline |
Sets the plot level.
All fields have an associated plot level. If it is <= to plotLevel
, the field will be output.
- Parameters
-
[in] | plotLevel | the limit plotlevel |
Definition at line 110 of file VTKPolyDataWriterInterface.hpp.
◆ setWriteFaceElementsAs3D()
void geos::vtk::VTKPolyDataWriterInterface::setWriteFaceElementsAs3D |
( |
bool |
writeFaceElementsAs3D | ) |
|
|
inline |
Defines whether in the vtk output facelements should be 2D or 3D.
- Parameters
-
writeFaceElementsAs3D | The boolean flag. |
Definition at line 99 of file VTKPolyDataWriterInterface.hpp.
◆ setWriteGhostCells()
void geos::vtk::VTKPolyDataWriterInterface::setWriteGhostCells |
( |
bool |
writeGhostCells | ) |
|
|
inline |
Defines if the vtk outputs should contain the ghost cells.
- Parameters
-
writeGhostCells | The boolean flag. |
Definition at line 90 of file VTKPolyDataWriterInterface.hpp.
◆ write()
Main method of this class. Write all the files for one time step.
This method writes a .pvd file (if a previous one was created from a precedent time step, it is overwritten). The .pvd file contains relative path to every .vtm files (one vtm file per time step). This method triggers also the writing of a .vtm file. A .vtm file containts relative paths to blocks with the following hierarchy :
- CellElementRegion
- CellElementRegion1
- CellElementRegion2
- ... -WellElementRegion
- Well1
- Well2
- rank0
- rank1
- rank2
- ...
- Parameters
-
[in] | time | the time step to be written |
[in] | cycle | the current cycle of event |
[in] | domain | the computation domain of this rank |
◆ writeCellElementRegions()
Writes the files for all the CellElementRegions.
There will be one file written per CellElementRegion and per rank.
- Parameters
-
[in] | time | the time-step |
[in] | elemManager | the ElementRegionManager containing the CellElementRegions to be output |
[in] | nodeManager | the NodeManager containing the nodes of the domain to be output |
[in] | path | the path to the file to output |
◆ writeElementFields()
void geos::vtk::VTKPolyDataWriterInterface::writeElementFields |
( |
ElementRegionBase const & |
subRegion, |
|
|
vtkCellData * |
cellData |
|
) |
| const |
|
protected |
Writes all the fields associated to the elements of er
if their plotlevel is <= m_plotLevel.
- Parameters
-
[in] | subRegion | ElementRegion being written |
[in] | cellData | a VTK object containing all the fields associated with the elements |
◆ writeNodeFields()
void geos::vtk::VTKPolyDataWriterInterface::writeNodeFields |
( |
NodeManager const & |
nodeManager, |
|
|
arrayView1d< localIndex const > const & |
nodeIndices, |
|
|
vtkPointData * |
pointData |
|
) |
| const |
|
protected |
Write all the fields associated to the nodes of nodeManager
if their plotlevel is <= m_plotLevel.
- Parameters
-
[in] | pointData | a VTK object containing all the fields associated with the nodes |
[in] | nodeIndices | list of local node indices to write |
[in] | nodeManager | the NodeManager associated with the domain being written |
◆ writeParticleFields()
void geos::vtk::VTKPolyDataWriterInterface::writeParticleFields |
( |
ParticleRegionBase const & |
region, |
|
|
vtkCellData * |
cellData |
|
) |
| const |
|
protected |
Writes all the fields associated to the elements of region
if their plotlevel is <= m_plotLevel.
- Parameters
-
[in] | region | ParticleRegion being written |
[in] | cellData | a VTK object containing all the fields associated with the elements |
◆ writeParticleRegions()
void geos::vtk::VTKPolyDataWriterInterface::writeParticleRegions |
( |
real64 const |
time, |
|
|
ParticleManager const & |
particleManager, |
|
|
string const & |
path |
|
) |
| |
|
protected |
Writes the files containing the particle representation.
There will be one file written per ParticleRegion and per rank
- Parameters
-
[in] | time | the time-step |
[in] | particleManager | the ParticleManager containing the ParticleRegions to be output |
[in] | path | the root path where the mesh will be written |
◆ writeSurfaceElementRegions()
Writes the files containing the faces elements.
There will be one file written per FaceElementRegion and per rank
- Parameters
-
[in] | time | the time-step |
[in] | elemManager | the ElementRegionManager containing the FaceElementRegions to be output |
[in] | nodeManager | the NodeManager containing the nodes of the domain to be output |
[in] | embSurfNodeManager | the embedded surface node Manager. |
[in] | faceManager | the faceManager. |
[in] | path | the path to the output file. |
◆ writeUnstructuredGrid()
void geos::vtk::VTKPolyDataWriterInterface::writeUnstructuredGrid |
( |
string const & |
path, |
|
|
ObjectManagerBase const & |
region, |
|
|
vtkUnstructuredGrid * |
ug |
|
) |
| |
|
protected |
Writes an unstructured grid.
The unstructured grid is the last element in the hierarchy of the output, it contains the cells connectivities and the vertices coordinates as long as the data fields associated with it
- Parameters
-
[in] | path | directory path for the grid file |
[in] | region | ElementRegionBase beeing written |
[in] | ug | a VTK SmartPointer to the VTK unstructured grid. |
◆ writeVtmFile()
Writes a VTM file for the time-step time
.
a VTM file is a VTK Multiblock file. It contains relative path to different files organized in blocks.
- Parameters
-
[in] | cycle | the current cycle number |
[in] | domain | the DomainPartition containing all the regions to be output and referred to in the VTM file |
[in] | vtmWriter | a writer specialized for the VTM file format |
◆ writeWellElementRegions()
Writes the files containing the well representation.
There will be one file written per WellElementRegion and per rank
- Parameters
-
[in] | time | the time-step |
[in] | elemManager | the ElementRegionManager containing the WellElementRegions to be output |
[in] | nodeManager | the NodeManager containing the nodes of the domain to be output |
[in] | path | The path to the file to output |
◆ m_pvd
A writter specialized for PVD files. There is one PVD file per simulation. It is the root file containing all the paths to the VTM files.
Definition at line 344 of file VTKPolyDataWriterInterface.hpp.
The documentation for this class was generated from the following file: