20 #ifndef GEOS_MESH_GENERATORS_VTKUTILITIES_HPP
21 #define GEOS_MESH_GENERATORS_VTKUTILITIES_HPP
27 #include <vtkDataSet.h>
28 #include <vtkMultiProcessController.h>
29 #include <vtkSmartPointer.h>
32 #include <unordered_set>
59 using CellMapType = std::map< ElementType, std::unordered_map< int, std::vector< vtkIdType > > >;
80 AllMeshes( vtkSmartPointer< vtkDataSet >
const & main,
81 std::map<
string, vtkSmartPointer< vtkDataSet > >
const & faceBlocks )
83 m_faceBlocks( faceBlocks )
117 m_faceBlocks = faceBlocks;
122 vtkSmartPointer< vtkDataSet > m_main;
125 std::map< string, vtkSmartPointer< vtkDataSet > > m_faceBlocks;
136 string const & mainBlockName,
160 vtkSmartPointer< vtkDataSet > loadedMesh,
161 std::map<
string, vtkSmartPointer< vtkDataSet > > & namesToFractures,
164 int const partitionRefinement,
165 int const useGlobalIds );
176 string const & attributeName );
187 MPI_Comm
const comm );
204 bool hasArray( vtkDataSet & mesh,
string const & sourceName );
221 vtkDataArray * vtkArray,
231 vtkDataArray * vtkArray,
vtkDataArray * findArrayForImport(vtkDataSet &mesh, string const &sourceName)
Collect the data to be imported.
vtkSmartPointer< vtkMultiProcessController > getController()
Return a VTK controller for multiprocessing.
void printMeshStatistics(vtkDataSet &mesh, CellMapType const &cellMap, MPI_Comm const comm)
Print statistics for a vtk mesh.
std::map< ElementType, std::unordered_map< int, std::vector< vtkIdType > > > CellMapType
Type of map used to store cell lists.
std::vector< int > findNeighborRanks(std::vector< vtkBoundingBox > boundingBoxes)
Compute the rank neighbor candidate list.
AllMeshes 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.
string buildCellBlockName(ElementType const type, int const regionId)
build cell block name from regionId and cellType
void importRegularField(std::vector< vtkIdType > const &cellIds, vtkDataArray *vtkArray, dataRepository::WrapperBase &wrapper)
Imports 1d and 2d arrays from vtkArray to wrapper, only for cellIds.
AllMeshes loadAllMeshes(Path const &filePath, string const &mainBlockName, array1d< string > const &faceBlockNames)
Load the VTK file into the VTK data structure.
CellMapType buildCellMap(vtkDataSet &mesh, string const &attributeName)
Collect lists of VTK cell indices organized by type and attribute value.
bool hasArray(vtkDataSet &mesh, string const &sourceName)
Check if the vtk mesh as a cell data field of name sourceName.
PartitionMethod
Choice of advanced mesh partitioner.
@ ptscotch
Use PTScotch library.
@ parmetis
Use ParMETIS library.
void importMaterialField(std::vector< vtkIdType > const &cellIds, vtkDataArray *vtkArray, dataRepository::WrapperBase &wrapper)
Imports 2d and 3d arrays from vtkArray to wrapper, only for cellIds.
The CellBlockManager class provides an interface to ObjectManagerBase in order to manage CellBlock da...
Class describing a file Path.
Base class for all wrappers containing common operations.
Gathers all the vtk meshes together.
std::map< string, vtkSmartPointer< vtkDataSet > > & getFaceBlocks()
void setMainMesh(vtkSmartPointer< vtkDataSet > main)
Defines the main 3d mesh for the simulation.
vtkSmartPointer< vtkDataSet > getMainMesh()
void setFaceBlocks(std::map< string, vtkSmartPointer< vtkDataSet > > const &faceBlocks)
Defines the face blocks/fractures.
AllMeshes(vtkSmartPointer< vtkDataSet > const &main, std::map< string, vtkSmartPointer< vtkDataSet > > const &faceBlocks)
Builds the compound from values.
ENUM_STRINGS(LinearSolverParameters::SolverType, "direct", "cg", "gmres", "fgmres", "bicgstab", "preconditioner")
Declare strings associated with enumeration values.
array1d< string > string_array
A 1-dimensional array of geos::string types.
double real64
64-bit floating point type.
std::int32_t integer
Signed integer type.
mapBase< TKEY, TVAL, std::integral_constant< bool, true > > map
Ordered map type.
void writeSurfaces(integer const logLevel, vtkDataSet &mesh, const geos::vtk::CellMapType &cellMap, CellBlockManager &cellBlockManager)
Build the "surface" node sets from the surface information.
ElementType
Denotes type of cell/element shape.
void writeCells(integer const logLevel, vtkDataSet &mesh, const geos::vtk::CellMapType &cellMap, CellBlockManager &cellBlockManager)
Build all the cell blocks.
Array< T, 1 > array1d
Alias for 1D array.
real64 writeNodes(integer const logLevel, vtkDataSet &mesh, string_array &nodesetNames, CellBlockManager &cellBlockManager, const geos::R1Tensor &translate, const geos::R1Tensor &scale)
Build all the vertex blocks.