GEOS
|
#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... | |
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.
|
strong |
Choice of advanced mesh partitioner.
Enumerator | |
---|---|
parmetis | Use ParMETIS library. |
ptscotch | Use PTScotch library. |
Definition at line 42 of file VTKUtilities.hpp.
string geos::vtk::buildCellBlockName | ( | ElementType const | type, |
int const | regionId | ||
) |
build cell block name from regionId and cellType
[in] | type | The type of element in the region |
[in] | regionId | The region considered |
CellMapType geos::vtk::buildCellMap | ( | vtkDataSet & | mesh, |
string const & | attributeName | ||
) |
Collect lists of VTK cell indices organized by type and attribute value.
[in] | mesh | the vtkUnstructuredGrid or vtkStructuredGrid that is loaded |
[in] | attributeName | name of the VTK data array containing the attribute, if any |
vtkDataArray* geos::vtk::findArrayForImport | ( | vtkDataSet & | mesh, |
string const & | sourceName | ||
) |
Collect the data to be imported.
[in] | mesh | an input mesh |
[in] | sourceName | a field name |
std::vector< int > geos::vtk::findNeighborRanks | ( | std::vector< vtkBoundingBox > | boundingBoxes | ) |
Compute the rank neighbor candidate list.
[in] | boundingBoxes | the bounding boxes used by the VTK partitioner for all ranks |
vtkSmartPointer< vtkMultiProcessController > geos::vtk::getController | ( | ) |
Return a VTK controller for multiprocessing.
bool geos::vtk::hasArray | ( | vtkDataSet & | mesh, |
string const & | sourceName | ||
) |
Check if the vtk mesh as a cell data field of name sourceName
.
[in] | mesh | an input mesh |
[in] | sourceName | a field name |
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
.
cellIds | The cells for which we should copy the data. |
vtkArray | The source. |
wrapper | The destination. |
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
.
cellIds | The cells for which we should copy the data. |
vtkArray | The source. |
wrapper | The destination. |
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.
vtkArray | The source. |
wrapper | The destination. |
AllMeshes geos::vtk::loadAllMeshes | ( | Path const & | filePath, |
string const & | mainBlockName, | ||
array1d< string > const & | faceBlockNames | ||
) |
Load the VTK file into the VTK data structure.
[in] | filePath | the Path of the file to load |
[in] | mainBlockName | The name of the block to import (will be considered for multi-block files only). |
[in] | faceBlockNames | The names of the face blocks to import (will be considered for multi-block files only). |
void geos::vtk::printMeshStatistics | ( | vtkDataSet & | mesh, |
CellMapType const & | cellMap, | ||
MPI_Comm const | comm | ||
) |
Print statistics for a vtk mesh.
[in] | mesh | an input mesh |
cellMap | a map of cell lists grouped by type | |
comm | the MPI communicator |
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.
[in] | logLevel | the log level |
[in] | loadedMesh | the mesh that was loaded on one or several MPI ranks |
[in] | namesToFractures | the fracture meshes |
[in] | comm | the MPI communicator |
[in] | method | the partitionning method |
[in] | partitionRefinement | number of graph partitioning refinement cycles |
[in] | useGlobalIds | controls whether global id arrays from the vtk input should be used |