GEOS
VTKSuperCellPartitioning.hpp
1 
16 #ifndef GEOS_MESH_GENERATORS_VTKSUPERCELLPARTITIONING_HPP
17 #define GEOS_MESH_GENERATORS_VTKSUPERCELLPARTITIONING_HPP
18 
19 #include "common/DataTypes.hpp"
20 #include "common/MpiWrapper.hpp"
21 #include "common/TimingMacros.hpp"
24 #include "LvArray/src/ArrayOfArrays.hpp"
25 
26 #include <vtkSmartPointer.h>
27 
28 class vtkUnstructuredGrid;
29 class vtkDataSet;
30 
31 namespace geos
32 {
33 namespace vtk
34 {
35 
39 enum class InitialDistributionStrategy
40 {
41  MORTON,
42  BLOCK
43 };
44 
51 {
54 
57 
59  std::set< vtkIdType > atomicSuperCells;
60 };
61 
73 SuperCellInfo tagCellsWithSuperCellIds(
74  vtkSmartPointer< vtkUnstructuredGrid > cells3D,
75  stdMap< string, ArrayOfArrays< localIndex, int64_t > > const & fractureNeighbors,
76  integer fractureWeight );
77 
88 SuperCellInfo reconstructSuperCellInfo( vtkSmartPointer< vtkUnstructuredGrid > mesh,
89  integer fractureWeight );
90 
102 vtkSmartPointer< vtkDataSet >
103 redistributeBySuperCellBlocks( vtkSmartPointer< vtkUnstructuredGrid > cells3D,
104  MPI_Comm comm,
105  InitialDistributionStrategy strategy = InitialDistributionStrategy::MORTON );
106 
120 std::pair< ArrayOfArrays< pmet_idx_t, pmet_idx_t >, array1d< pmet_idx_t > >
121 buildSuperCellGraph(
122  vtkSmartPointer< vtkUnstructuredGrid > cells3D,
123  ArrayOfArrays< pmet_idx_t, pmet_idx_t > const & baseGraph,
124  arrayView1d< pmet_idx_t const > const & baseElemDist,
125  SuperCellInfo const & info,
126  MPI_Comm comm );
127 
139 void validateSuperCellGraph(
140  ArrayOfArrays< pmet_idx_t, pmet_idx_t > const & superGraph,
141  arrayView1d< pmet_idx_t const > const & superElemDist,
142  arrayView1d< pmet_idx_t const > const & vertexWeights,
143  MPI_Comm comm );
144 
158 expandSuperCellPartitioningToCells(
159  vtkSmartPointer< vtkUnstructuredGrid > cells3D,
160  array1d< int64_t > const & superPartitioning,
161  stdMap< vtkIdType, localIndex > const & superCellIdToLocalIdx,
162  MPI_Comm comm );
163 
164 } // namespace vtk
165 } // namespace geos
166 
167 #endif /* GEOS_MESH_GENERATORS_VTKSUPERCELLPARTITIONING_HPP */
ArrayView< T, 1 > arrayView1d
Alias for 1D array view.
Definition: DataTypes.hpp:179
internal::StdMapWrapper< std::map< Key, T, Compare, Allocator >, USE_STD_CONTAINER_BOUNDS_CHECKING > stdMap
int integer
Signed integer type.
Definition: DataTypes.hpp:81
Array< T, 1 > array1d
Alias for 1D array.
Definition: DataTypes.hpp:175
LvArray::ArrayOfArrays< T, INDEX_TYPE, LvArray::ChaiBuffer > ArrayOfArrays
Array of variable-sized arrays. See LvArray::ArrayOfArrays for details.
Definition: DataTypes.hpp:281
Super-cell metadata for constrained partitioning.
stdMap< vtkIdType, vtkIdType > vertexWeights
SuperCellId to vertex weight.
stdMap< vtkIdType, stdVector< vtkIdType > > superCellToOriginalCells
SuperCellId to global cell IDs.
std::set< vtkIdType > atomicSuperCells
SuperCellIds containing multiple cells.