GEOS
VTKUtilities.hpp
Go to the documentation of this file.
1 /*
2  * ------------------------------------------------------------------------------------------------------------
3  * SPDX-License-Identifier: LGPL-2.1-only
4  *
5  * Copyright (c) 2016-2024 Lawrence Livermore National Security LLC
6  * Copyright (c) 2018-2024 TotalEnergies
7  * Copyright (c) 2018-2024 The Board of Trustees of the Leland Stanford Junior University
8  * Copyright (c) 2023-2024 Chevron
9  * Copyright (c) 2019- GEOS/GEOSX Contributors
10  * All rights reserved
11  *
12  * See top level LICENSE, COPYRIGHT, CONTRIBUTORS, NOTICE, and ACKNOWLEDGEMENTS files for details.
13  * ------------------------------------------------------------------------------------------------------------
14  */
15 
20 #ifndef GEOS_MESH_GENERATORS_VTKUTILITIES_HPP
21 #define GEOS_MESH_GENERATORS_VTKUTILITIES_HPP
22 
23 #include "common/DataTypes.hpp"
24 #include "common/MpiWrapper.hpp"
26 
27 #include <vtkDataSet.h>
28 #include <vtkMultiProcessController.h>
29 #include <vtkSmartPointer.h>
30 
31 #include <numeric>
32 #include <unordered_set>
33 
34 namespace geos
35 {
36 namespace vtk
37 {
38 class CollocatedNodes;
39 
44 {
45  parmetis,
46  ptscotch,
47 };
48 
51  "parmetis",
52  "ptscotch" );
53 
61 
66 vtkSmartPointer< vtkMultiProcessController > getController();
67 
71 class AllMeshes
72 {
73 public:
74  AllMeshes() = default;
75 
81  AllMeshes( vtkSmartPointer< vtkDataSet > const & main,
82  stdMap< string, vtkSmartPointer< vtkDataSet > > const & faceBlocks )
83  : m_main( main ),
84  m_faceBlocks( faceBlocks )
85  { }
86 
90  vtkSmartPointer< vtkDataSet > getMainMesh()
91  {
92  return m_main;
93  }
94 
99  {
100  return m_faceBlocks;
101  }
102 
107  void setMainMesh( vtkSmartPointer< vtkDataSet > main )
108  {
109  m_main = main;
110  }
111 
116  void setFaceBlocks( stdMap< string, vtkSmartPointer< vtkDataSet > > const & faceBlocks )
117  {
118  m_faceBlocks = faceBlocks;
119  }
120 
121 private:
123  vtkSmartPointer< vtkDataSet > m_main;
124 
127 };
128 
136 AllMeshes loadAllMeshes( Path const & filePath,
137  string const & mainBlockName,
138  string_array const & faceBlockNames );
139 
147 
162 AllMeshes
163 redistributeMeshes( integer const logLevel,
164  vtkSmartPointer< vtkDataSet > loadedMesh,
165  stdMap< string, vtkSmartPointer< vtkDataSet > > & namesToFractures,
166  MPI_Comm const comm,
167  PartitionMethod const method,
168  int const partitionRefinement,
169  int const partitionFractureWeight,
170  int const useGlobalIds,
171  string const & structuredIndexAttributeName,
172  int const numPartZ );
173 
182 CellMapType buildCellMap( vtkDataSet & mesh,
183  string const & attributeName );
184 
192 void printMeshStatistics( vtkDataSet & mesh,
193  CellMapType const & cellMap,
194  MPI_Comm const comm );
195 
202 vtkDataArray * findArrayForImport( vtkDataSet & mesh, string const & sourceName );
203 
211 bool hasArray( vtkDataSet & mesh, string const & sourceName );
212 
219 string buildCellBlockName( ElementType const type, int const regionId );
220 
228  vtkDataArray * vtkArray,
229  dataRepository::WrapperBase & wrapper );
230 
238  vtkDataArray * vtkArray,
239  dataRepository::WrapperBase & wrapper );
240 
246 void importRegularField( vtkDataArray * vtkArray,
247  dataRepository::WrapperBase & wrapper );
248 
261  vtkIdList * fractureNodeIds,
262  CollocatedNodes const & collocatedNodes,
263  stdMap< vtkIdType, std::set< vtkIdType > > const & nodesToCells );
264 
265 
266 } // namespace vtk
267 
277 void writeNodes( integer const logLevel,
278  vtkDataSet & mesh,
279  string_array & nodesetNames,
280  CellBlockManager & cellBlockManager,
281  const geos::R1Tensor & translate,
282  const geos::R1Tensor & scale );
283 
292 void writeCells( integer const logLevel,
293  vtkDataSet & mesh,
294  vtk::CellMapType const & cellMap,
295  string const & structuredIndexAttributeName,
296  CellBlockManager & cellBlockManager );
297 
307 void writeSurfaces( integer const logLevel,
308  vtkDataSet & mesh,
309  const geos::vtk::CellMapType & cellMap,
310  CellBlockManager & cellBlockManager );
311 
317 std::pair< real64, real64 > getGlobalLengthAndOffset( vtkDataSet & mesh );
318 
319 } // namespace geos
320 
321 #endif /* GEOS_MESH_GENERATORS_VTKUTILITIES_HPP */
stdMap< ElementType, stdUnorderedMap< int, stdVector< vtkIdType > > > CellMapType
Type of map used to store cell lists.
AllMeshes loadAllMeshes(Path const &filePath, string const &mainBlockName, string_array const &faceBlockNames)
Load the VTK file into the VTK data structure.
void importMaterialField(stdVector< vtkIdType > const &cellIds, vtkDataArray *vtkArray, dataRepository::WrapperBase &wrapper)
Imports 2d and 3d arrays from vtkArray to wrapper, only for cellIds.
vtkDataArray * findArrayForImport(vtkDataSet &mesh, string const &sourceName)
Collect the data to be imported.
AllMeshes redistributeMeshes(integer const logLevel, vtkSmartPointer< vtkDataSet > loadedMesh, stdMap< string, vtkSmartPointer< vtkDataSet > > &namesToFractures, MPI_Comm const comm, PartitionMethod const method, int const partitionRefinement, int const partitionFractureWeight, int const useGlobalIds, string const &structuredIndexAttributeName, int const numPartZ)
Generate global point/cell IDs and redistribute the mesh among MPI ranks.
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.
stdVector< vtkIdType > findMatchingCellsForFractureElement(vtkIdList *fractureNodeIds, CollocatedNodes const &collocatedNodes, stdMap< vtkIdType, std::set< vtkIdType > > const &nodesToCells)
Find 3D cells whose faces exactly match a fracture element.
void importRegularField(stdVector< vtkIdType > const &cellIds, vtkDataArray *vtkArray, dataRepository::WrapperBase &wrapper)
Imports 1d and 2d arrays from vtkArray to wrapper, only for cellIds.
string buildCellBlockName(ElementType const type, int const regionId)
build cell block name from regionId and cellType
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.
stdVector< int > findNeighborRanks(stdVector< vtkBoundingBox > boundingBoxes)
Compute the rank neighbor candidate list.
The CellBlockManager class provides an interface to ObjectManagerBase in order to manage CellBlock da...
Class describing a file Path.
Definition: Path.hpp:35
Base class for all wrappers containing common operations.
Definition: WrapperBase.hpp:56
Gathers all the vtk meshes together.
AllMeshes(vtkSmartPointer< vtkDataSet > const &main, stdMap< string, vtkSmartPointer< vtkDataSet > > const &faceBlocks)
Builds the compound from values.
void setMainMesh(vtkSmartPointer< vtkDataSet > main)
Defines the main 3d mesh for the simulation.
vtkSmartPointer< vtkDataSet > getMainMesh()
void setFaceBlocks(stdMap< string, vtkSmartPointer< vtkDataSet > > const &faceBlocks)
Defines the face blocks/fractures.
stdMap< string, vtkSmartPointer< vtkDataSet > > & getFaceBlocks()
Convenience wrapper around the raw vtk information.
stdVector< string > string_array
A 1-dimensional array of geos::string types.
Definition: DataTypes.hpp:361
std::pair< real64, real64 > getGlobalLengthAndOffset(vtkDataSet &mesh)
Compute the global length of the mesh and its offset.
void writeNodes(integer const logLevel, vtkDataSet &mesh, string_array &nodesetNames, CellBlockManager &cellBlockManager, const geos::R1Tensor &translate, const geos::R1Tensor &scale)
Build all the vertex blocks.
void writeCells(integer const logLevel, vtkDataSet &mesh, vtk::CellMapType const &cellMap, string const &structuredIndexAttributeName, CellBlockManager &cellBlockManager)
Build all the cell blocks.
internal::StdMapWrapper< std::map< Key, T, Compare, Allocator >, USE_STD_CONTAINER_BOUNDS_CHECKING > stdMap
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.
Definition: ElementType.hpp:32
int integer
Signed integer type.
Definition: DataTypes.hpp:81
ENUM_STRINGS(LinearSolverParameters::SolverType, "direct", "cg", "gmres", "fgmres", "bicgstab", "richardson", "preconditioner")
Declare strings associated with enumeration values.
internal::StdVectorWrapper< T, Allocator, USE_STD_CONTAINER_BOUNDS_CHECKING > stdVector