GEOS
Classes | Namespaces | Typedefs | Enumerations | Functions
VTKUtilities.hpp File Reference
#include "common/DataTypes.hpp"
#include "common/MpiWrapper.hpp"
#include "mesh/generators/CellBlockManager.hpp"
#include <vtkDataSet.h>
#include <vtkMultiProcessController.h>
#include <vtkSmartPointer.h>
#include <numeric>
#include <unordered_set>

Go to the source code of this file.

Classes

class  geos::vtk::AllMeshes
 Gathers all the vtk meshes together. More...
 

Namespaces

 geos
 

Typedefs

using geos::vtk::CellMapType = std::map< ElementType, std::unordered_map< int, std::vector< vtkIdType > > >
 Type of map used to store cell lists. More...
 

Enumerations

enum class  geos::vtk::PartitionMethod : integer { parmetis , ptscotch }
 Choice of advanced mesh partitioner. More...
 

Functions

 geos::vtk::ENUM_STRINGS (PartitionMethod, "parmetis", "ptscotch")
 Strings for VTKMeshGenerator::PartitionMethod enumeration.
 
vtkSmartPointer< vtkMultiProcessController > geos::vtk::getController ()
 Return a VTK controller for multiprocessing. More...
 
AllMeshes geos::vtk::loadAllMeshes (Path const &filePath, string const &mainBlockName, array1d< string > const &faceBlockNames)
 Load the VTK file into the VTK data structure. More...
 
std::vector< int > geos::vtk::findNeighborRanks (std::vector< vtkBoundingBox > boundingBoxes)
 Compute the rank neighbor candidate list. More...
 
AllMeshes geos::vtk::redistributeMeshes (integer const logLevel, vtkSmartPointer< vtkDataSet > loadedMesh, std::map< string, vtkSmartPointer< vtkDataSet > > &namesToFractures, MPI_Comm const comm, PartitionMethod const method, int const partitionRefinement, int const useGlobalIds)
 Generate global point/cell IDs and redistribute the mesh among MPI ranks. More...
 
CellMapType geos::vtk::buildCellMap (vtkDataSet &mesh, string const &attributeName)
 Collect lists of VTK cell indices organized by type and attribute value. More...
 
void geos::vtk::printMeshStatistics (vtkDataSet &mesh, CellMapType const &cellMap, MPI_Comm const comm)
 Print statistics for a vtk mesh. More...
 
vtkDataArray * geos::vtk::findArrayForImport (vtkDataSet &mesh, string const &sourceName)
 Collect the data to be imported. More...
 
bool geos::vtk::hasArray (vtkDataSet &mesh, string const &sourceName)
 Check if the vtk mesh as a cell data field of name sourceName. More...
 
string geos::vtk::buildCellBlockName (ElementType const type, int const regionId)
 build cell block name from regionId and cellType More...
 
void geos::vtk::importMaterialField (std::vector< vtkIdType > const &cellIds, vtkDataArray *vtkArray, dataRepository::WrapperBase &wrapper)
 Imports 2d and 3d arrays from vtkArray to wrapper, only for cellIds. More...
 
void geos::vtk::importRegularField (std::vector< vtkIdType > const &cellIds, vtkDataArray *vtkArray, dataRepository::WrapperBase &wrapper)
 Imports 1d and 2d arrays from vtkArray to wrapper, only for cellIds. More...
 
void geos::vtk::importRegularField (vtkDataArray *vtkArray, dataRepository::WrapperBase &wrapper)
 Imports 1d and 2d arrays from vtkArray to wrapper, for all the elements/cells of the provided wrapper. More...
 
real64 geos::writeNodes (integer const logLevel, vtkDataSet &mesh, string_array &nodesetNames, CellBlockManager &cellBlockManager, const geos::R1Tensor &translate, const geos::R1Tensor &scale)
 Build all the vertex blocks. More...
 
void geos::writeCells (integer const logLevel, vtkDataSet &mesh, const geos::vtk::CellMapType &cellMap, CellBlockManager &cellBlockManager)
 Build all the cell blocks. More...
 
void geos::writeSurfaces (integer const logLevel, vtkDataSet &mesh, const geos::vtk::CellMapType &cellMap, CellBlockManager &cellBlockManager)
 Build the "surface" node sets from the surface information. More...
 

Typedef Documentation

◆ CellMapType

using geos::vtk::CellMapType = typedef std::map< ElementType, std::unordered_map< int, std::vector< vtkIdType > > >

Type of map used to store cell lists.

This should be an unordered_map, but some outdated standard libraries on some systems do not provide std::hash specialization for enums. This is not performance critical though.

Definition at line 59 of file VTKUtilities.hpp.

Enumeration Type Documentation

◆ PartitionMethod

enum geos::vtk::PartitionMethod : integer
strong

Choice of advanced mesh partitioner.

Enumerator
parmetis 

Use ParMETIS library.

ptscotch 

Use PTScotch library.

Definition at line 42 of file VTKUtilities.hpp.

Function Documentation

◆ buildCellBlockName()

string geos::vtk::buildCellBlockName ( ElementType const  type,
int const  regionId 
)

build cell block name from regionId and cellType

Parameters
[in]typeThe type of element in the region
[in]regionIdThe region considered
Returns
The cell block name

◆ buildCellMap()

CellMapType geos::vtk::buildCellMap ( vtkDataSet &  mesh,
string const &  attributeName 
)

Collect lists of VTK cell indices organized by type and attribute value.

Parameters
[in]meshthe vtkUnstructuredGrid or vtkStructuredGrid that is loaded
[in]attributeNamename of the VTK data array containing the attribute, if any
Returns
A map from element type to a map of attribute to the associated cell ids for the current rank. The map contains entries for all types and attribute values across all MPI ranks, even if there are no cells on current rank (then the list will be empty).

◆ findArrayForImport()

vtkDataArray* geos::vtk::findArrayForImport ( vtkDataSet &  mesh,
string const &  sourceName 
)

Collect the data to be imported.

Parameters
[in]meshan input mesh
[in]sourceNamea field name
Returns
A VTK data array pointer.

◆ findNeighborRanks()

std::vector< int > geos::vtk::findNeighborRanks ( std::vector< vtkBoundingBox >  boundingBoxes)

Compute the rank neighbor candidate list.

Parameters
[in]boundingBoxesthe bounding boxes used by the VTK partitioner for all ranks
Returns
the list of neighboring MPI ranks, will be updated

◆ getController()

vtkSmartPointer< vtkMultiProcessController > geos::vtk::getController ( )

Return a VTK controller for multiprocessing.

Returns
Return a VTK controller for multiprocessing.

◆ hasArray()

bool geos::vtk::hasArray ( vtkDataSet &  mesh,
string const &  sourceName 
)

Check if the vtk mesh as a cell data field of name sourceName.

Parameters
[in]meshan input mesh
[in]sourceNamea field name
Returns
The boolean result.
Note
No check is performed to see if the type of the vtk array matches any requirement.

◆ importMaterialField()

void geos::vtk::importMaterialField ( std::vector< vtkIdType > const &  cellIds,
vtkDataArray *  vtkArray,
dataRepository::WrapperBase wrapper 
)

Imports 2d and 3d arrays from vtkArray to wrapper, only for cellIds.

Parameters
cellIdsThe cells for which we should copy the data.
vtkArrayThe source.
wrapperThe destination.

◆ importRegularField() [1/2]

void geos::vtk::importRegularField ( std::vector< vtkIdType > const &  cellIds,
vtkDataArray *  vtkArray,
dataRepository::WrapperBase wrapper 
)

Imports 1d and 2d arrays from vtkArray to wrapper, only for cellIds.

Parameters
cellIdsThe cells for which we should copy the data.
vtkArrayThe source.
wrapperThe destination.

◆ importRegularField() [2/2]

void geos::vtk::importRegularField ( vtkDataArray *  vtkArray,
dataRepository::WrapperBase wrapper 
)

Imports 1d and 2d arrays from vtkArray to wrapper, for all the elements/cells of the provided wrapper.

Parameters
vtkArrayThe source.
wrapperThe destination.

◆ loadAllMeshes()

AllMeshes geos::vtk::loadAllMeshes ( Path const &  filePath,
string const &  mainBlockName,
array1d< string > const &  faceBlockNames 
)

Load the VTK file into the VTK data structure.

Parameters
[in]filePaththe Path of the file to load
[in]mainBlockNameThe name of the block to import (will be considered for multi-block files only).
[in]faceBlockNamesThe names of the face blocks to import (will be considered for multi-block files only).
Returns
The compound of the main mesh and the face block meshes.

◆ printMeshStatistics()

void geos::vtk::printMeshStatistics ( vtkDataSet &  mesh,
CellMapType const &  cellMap,
MPI_Comm const  comm 
)

Print statistics for a vtk mesh.

Parameters
[in]meshan input mesh
cellMapa map of cell lists grouped by type
commthe MPI communicator

◆ redistributeMeshes()

AllMeshes geos::vtk::redistributeMeshes ( integer const  logLevel,
vtkSmartPointer< vtkDataSet >  loadedMesh,
std::map< string, vtkSmartPointer< vtkDataSet > > &  namesToFractures,
MPI_Comm const  comm,
PartitionMethod const  method,
int const  partitionRefinement,
int const  useGlobalIds 
)

Generate global point/cell IDs and redistribute the mesh among MPI ranks.

Parameters
[in]logLevelthe log level
[in]loadedMeshthe mesh that was loaded on one or several MPI ranks
[in]namesToFracturesthe fracture meshes
[in]commthe MPI communicator
[in]methodthe partitionning method
[in]partitionRefinementnumber of graph partitioning refinement cycles
[in]useGlobalIdscontrols whether global id arrays from the vtk input should be used
Returns
the vtk grid redistributed