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 Total, S.A
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 
43 {
44  parmetis,
45  ptscotch,
46 };
47 
50  "parmetis",
51  "ptscotch" );
52 
59 using CellMapType = std::map< ElementType, std::unordered_map< int, std::vector< vtkIdType > > >;
60 
65 vtkSmartPointer< vtkMultiProcessController > getController();
66 
70 class AllMeshes
71 {
72 public:
73  AllMeshes() = default;
74 
80  AllMeshes( vtkSmartPointer< vtkDataSet > const & main,
81  std::map< string, vtkSmartPointer< vtkDataSet > > const & faceBlocks )
82  : m_main( main ),
83  m_faceBlocks( faceBlocks )
84  { }
85 
89  vtkSmartPointer< vtkDataSet > getMainMesh()
90  {
91  return m_main;
92  }
93 
97  std::map< string, vtkSmartPointer< vtkDataSet > > & getFaceBlocks()
98  {
99  return m_faceBlocks;
100  }
101 
106  void setMainMesh( vtkSmartPointer< vtkDataSet > main )
107  {
108  m_main = main;
109  }
110 
115  void setFaceBlocks( std::map< string, vtkSmartPointer< vtkDataSet > > const & faceBlocks )
116  {
117  m_faceBlocks = faceBlocks;
118  }
119 
120 private:
122  vtkSmartPointer< vtkDataSet > m_main;
123 
125  std::map< string, vtkSmartPointer< vtkDataSet > > m_faceBlocks;
126 };
127 
135 AllMeshes loadAllMeshes( Path const & filePath,
136  string const & mainBlockName,
137  array1d< string > const & faceBlockNames );
138 
144 std::vector< int >
145 findNeighborRanks( std::vector< vtkBoundingBox > boundingBoxes );
146 
158 AllMeshes
159 redistributeMeshes( integer const logLevel,
160  vtkSmartPointer< vtkDataSet > loadedMesh,
161  std::map< string, vtkSmartPointer< vtkDataSet > > & namesToFractures,
162  MPI_Comm const comm,
163  PartitionMethod const method,
164  int const partitionRefinement,
165  int const useGlobalIds );
166 
175 CellMapType buildCellMap( vtkDataSet & mesh,
176  string const & attributeName );
177 
185 void printMeshStatistics( vtkDataSet & mesh,
186  CellMapType const & cellMap,
187  MPI_Comm const comm );
188 
195 vtkDataArray * findArrayForImport( vtkDataSet & mesh, string const & sourceName );
196 
204 bool hasArray( vtkDataSet & mesh, string const & sourceName );
205 
212 string buildCellBlockName( ElementType const type, int const regionId );
213 
220 void importMaterialField( std::vector< vtkIdType > const & cellIds,
221  vtkDataArray * vtkArray,
222  dataRepository::WrapperBase & wrapper );
223 
230 void importRegularField( std::vector< vtkIdType > const & cellIds,
231  vtkDataArray * vtkArray,
232  dataRepository::WrapperBase & wrapper );
233 
239 void importRegularField( vtkDataArray * vtkArray,
240  dataRepository::WrapperBase & wrapper );
241 
242 
243 } // namespace vtk
244 
255 real64 writeNodes( integer const logLevel,
256  vtkDataSet & mesh,
257  string_array & nodesetNames,
258  CellBlockManager & cellBlockManager,
259  const geos::R1Tensor & translate,
260  const geos::R1Tensor & scale );
261 
269 void writeCells( integer const logLevel,
270  vtkDataSet & mesh,
271  const geos::vtk::CellMapType & cellMap,
272  CellBlockManager & cellBlockManager );
273 
283 void writeSurfaces( integer const logLevel,
284  vtkDataSet & mesh,
285  const geos::vtk::CellMapType & cellMap,
286  CellBlockManager & cellBlockManager );
287 
288 } // namespace geos
289 
290 #endif /* GEOS_MESH_GENERATORS_VTKUTILITIES_HPP */
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.
Definition: Path.hpp:33
Base class for all wrappers containing common operations.
Definition: WrapperBase.hpp:56
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.
Definition: DataTypes.hpp:392
double real64
64-bit floating point type.
Definition: DataTypes.hpp:99
std::int32_t integer
Signed integer type.
Definition: DataTypes.hpp:82
mapBase< TKEY, TVAL, std::integral_constant< bool, true > > map
Ordered map type.
Definition: DataTypes.hpp:369
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
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.
Definition: DataTypes.hpp:176
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.