GEOS
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
geos::vtk::VTKPolyDataWriterInterface Class Reference

Encapsulate output methods for vtk. More...

#include <VTKPolyDataWriterInterface.hpp>

Public Member Functions

 VTKPolyDataWriterInterface (string outputName)
 Constructor. More...
 
void setWriteGhostCells (bool writeGhostCells)
 Defines if the vtk outputs should contain the ghost cells. More...
 
void setWriteFaceElementsAs3D (bool writeFaceElementsAs3D)
 Defines whether in the vtk output facelements should be 2D or 3D. More...
 
void setPlotLevel (integer plotLevel)
 Sets the plot level. More...
 
void setOutputMode (VTKOutputMode mode)
 Set the binary mode. More...
 
void setOutputRegionType (VTKRegionTypes regionType)
 Set the output region type. More...
 
void setOutputLocation (string outputDir, string outputName)
 Set the output directory name. More...
 
void setOnlyPlotSpecifiedFieldNamesFlag (integer const onlyPlotSpecifiedFieldNames)
 Set the flag to decide whether we only plot the fields specified by fieldNames, or if we also plot fields based on plotLevel. More...
 
void setFieldNames (string_array const &fieldNames)
 Set the names of the fields to output. More...
 
void setLevelNames (string_array const &levelNames)
 Set the names of the mesh levels to output. More...
 
void setNumberOfTargetProcesses (integer const numberOfTargetProcesses)
 Set the Number Of Target Processes. More...
 
void write (real64 time, integer cycle, DomainPartition const &domain)
 Main method of this class. Write all the files for one time step. More...
 
void clearData ()
 Clears the datasets accumulated in the pvd writer.
 

Protected Member Functions

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 &region, 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 &region, vtkUnstructuredGrid *ug)
 Writes an unstructured grid. More...
 

Protected Attributes

string m_outputDir
 Output directory name.
 
string m_outputName
 Name of the PVD file and associated directory containing the outputs.
 
VTKPVDWriter m_pvd
 
bool m_writeGhostCells
 Should the vtk files contain the ghost cells or not.
 
dataRepository::PlotLevel m_plotLevel
 Maximum plot level to be written.
 
integer m_onlyPlotSpecifiedFieldNames
 Flag to decide whether we only plot the fields specified by fieldNames, or if we also plot fields based on plotLevel.
 
bool m_requireFieldRegistrationCheck
 Flag to decide whether we check that the specified fieldNames are actually registered.
 
std::set< stringm_fieldNames
 Names of the fields to output.
 
std::set< stringm_levelNames
 Names of the mesh levels to output (an empty array means all levels are saved)
 
integer m_previousCycle
 The previousCycle.
 
VTKOutputMode m_outputMode
 Output mode, could be ASCII or BINARAY.
 
VTKRegionTypes m_outputRegionType
 Region output type, could be CELL, WELL, SURFACE, or ALL.
 
bool m_writeFaceElementsAs3D
 Defines whether to plot a faceElement as a 3D volumetric element or not.
 
integer m_numberOfTargetProcesses
 Number of target processes to aggregate the data to be written.
 
std::map< string, stdVector< integer > > m_targetProcessesId
 Map a region name to the array of ranks outputed for it.
 

Detailed Description

Encapsulate output methods for vtk.

Definition at line 77 of file VTKPolyDataWriterInterface.hpp.

Constructor & Destructor Documentation

◆ VTKPolyDataWriterInterface()

geos::vtk::VTKPolyDataWriterInterface::VTKPolyDataWriterInterface ( string  outputName)
explicit

Constructor.

Parameters
[in]outputNamefolder name in which all the files will be written

Member Function Documentation

◆ isFieldPlotEnabled()

bool geos::vtk::VTKPolyDataWriterInterface::isFieldPlotEnabled ( dataRepository::WrapperBase const &  wrapper) const
protected

Check if plotting is enabled for this field.

Parameters
[in]wrapperthe wrapper
Returns
true if this wrapper should be plot, false otherwise

◆ setFieldNames()

void geos::vtk::VTKPolyDataWriterInterface::setFieldNames ( string_array const &  fieldNames)
inline

Set the names of the fields to output.

Parameters
[in]fieldNamesthe fields to output

Definition at line 158 of file VTKPolyDataWriterInterface.hpp.

◆ setLevelNames()

void geos::vtk::VTKPolyDataWriterInterface::setLevelNames ( string_array const &  levelNames)
inline

Set the names of the mesh levels to output.

Parameters
[in]levelNamesthe 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

Set the Number Of Target Processes.

Parameters
[in]numberOfTargetProcessesthe number of processes

Definition at line 175 of file VTKPolyDataWriterInterface.hpp.

◆ 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]onlyPlotSpecifiedFieldNamesthe 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]outputDirglobal output directory location
[in]outputNamename 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

Set the binary mode.

Parameters
[in]modeoutput mode to be set

Definition at line 119 of file VTKPolyDataWriterInterface.hpp.

◆ setOutputRegionType()

void geos::vtk::VTKPolyDataWriterInterface::setOutputRegionType ( VTKRegionTypes  regionType)
inline

Set the output region type.

Parameters
[in]regionTypeoutput region type to be set

Definition at line 128 of file VTKPolyDataWriterInterface.hpp.

◆ 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]plotLevelthe 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
writeFaceElementsAs3DThe 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
writeGhostCellsThe boolean flag.

Definition at line 90 of file VTKPolyDataWriterInterface.hpp.

◆ write()

void geos::vtk::VTKPolyDataWriterInterface::write ( real64  time,
integer  cycle,
DomainPartition const &  domain 
)

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
      • rank0
      • rank1
      • rank2
      • ...
    • CellElementRegion2
      • rank0
      • rank1
      • rank2
      • ...
    • ... -WellElementRegion
    • Well1
      • rank0
      • rank1
      • rank2
      • ...
    • Well2
      • rank0
      • rank1
      • rank2
      • ...
        Parameters
        [in]timethe time step to be written
        [in]cyclethe current cycle of event
        [in]domainthe computation domain of this rank

◆ writeCellElementRegions()

void geos::vtk::VTKPolyDataWriterInterface::writeCellElementRegions ( real64  time,
ElementRegionManager const &  elemManager,
NodeManager const &  nodeManager,
string const &  path 
)
protected

Writes the files for all the CellElementRegions.

There will be one file written per CellElementRegion and per rank.

Parameters
[in]timethe time-step
[in]elemManagerthe ElementRegionManager containing the CellElementRegions to be output
[in]nodeManagerthe NodeManager containing the nodes of the domain to be output
[in]paththe 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]subRegionElementRegion being written
[in]cellDataa 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]pointDataa VTK object containing all the fields associated with the nodes
[in]nodeIndiceslist of local node indices to write
[in]nodeManagerthe 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]regionParticleRegion being written
[in]cellDataa 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]timethe time-step
[in]particleManagerthe ParticleManager containing the ParticleRegions to be output
[in]paththe root path where the mesh will be written

◆ writeSurfaceElementRegions()

void geos::vtk::VTKPolyDataWriterInterface::writeSurfaceElementRegions ( real64  time,
ElementRegionManager const &  elemManager,
NodeManager const &  nodeManager,
EmbeddedSurfaceNodeManager const &  embSurfNodeManager,
FaceManager const &  faceManager,
string const &  path 
)
protected

Writes the files containing the faces elements.

There will be one file written per FaceElementRegion and per rank

Parameters
[in]timethe time-step
[in]elemManagerthe ElementRegionManager containing the FaceElementRegions to be output
[in]nodeManagerthe NodeManager containing the nodes of the domain to be output
[in]embSurfNodeManagerthe embedded surface node Manager.
[in]faceManagerthe faceManager.
[in]paththe 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]pathdirectory path for the grid file
[in]regionElementRegionBase beeing written
[in]uga VTK SmartPointer to the VTK unstructured grid.

◆ writeVtmFile()

void geos::vtk::VTKPolyDataWriterInterface::writeVtmFile ( integer const  cycle,
DomainPartition const &  domain,
VTKVTMWriter const &  vtmWriter 
) const
protected

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]cyclethe current cycle number
[in]domainthe DomainPartition containing all the regions to be output and referred to in the VTM file
[in]vtmWritera writer specialized for the VTM file format

◆ writeWellElementRegions()

void geos::vtk::VTKPolyDataWriterInterface::writeWellElementRegions ( real64  time,
ElementRegionManager const &  elemManager,
NodeManager const &  nodeManager,
string const &  path 
)
protected

Writes the files containing the well representation.

There will be one file written per WellElementRegion and per rank

Parameters
[in]timethe time-step
[in]elemManagerthe ElementRegionManager containing the WellElementRegions to be output
[in]nodeManagerthe NodeManager containing the nodes of the domain to be output
[in]pathThe path to the file to output

Member Data Documentation

◆ m_pvd

VTKPVDWriter geos::vtk::VTKPolyDataWriterInterface::m_pvd
protected

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: