GEOS
VTKMeshGenerator.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_VTKMESHGENERATOR_HPP
21 #define GEOS_MESH_GENERATORS_VTKMESHGENERATOR_HPP
22 
26 #include "mesh/mpiCommunications/SpatialPartition.hpp"
27 
28 class vtkDataSet;
29 
30 namespace geos
31 {
32 
38 {
39 public:
40 
46  VTKMeshGenerator( const string & name,
47  Group * const parent );
48 
53  static string catalogName() { return "VTKMesh"; }
54 
92  void fillCellBlockManager( CellBlockManager & cellBlockManager, SpatialPartition & partition ) override;
93 
95  string const & blockName,
96  string const & meshFieldName,
97  bool isMaterialField,
98  dataRepository::WrapperBase & wrapper ) const override;
99 
100  void freeResources() override;
101 
102 protected:
103  void postInputInitialization() override;
104 
105 private:
106 
109  {
110  constexpr static char const * regionAttributeString() { return "regionAttribute"; }
111  constexpr static char const * structuredIndexAttributeString() { return "structuredIndexAttribute"; }
112  constexpr static char const * mainBlockNameString() { return "mainBlockName"; }
113  constexpr static char const * faceBlockNamesString() { return "faceBlocks"; }
114  constexpr static char const * nodesetNamesString() { return "nodesetNames"; }
115  constexpr static char const * partitionRefinementString() { return "partitionRefinement"; }
116  constexpr static char const * partitionMethodString() { return "partitionMethod"; }
117  constexpr static char const * useGlobalIdsString() { return "useGlobalIds"; }
118  constexpr static char const * dataSourceString() { return "dataSourceName"; }
119  };
120 
121  struct groupKeyStruct
122  {
123  constexpr static char const * regionString() { return "VTKRegion"; }
124  };
126 
127  void importVolumicFieldOnArray( string const & cellBlockName,
128  string const & meshFieldName,
129  bool isMaterialField,
130  dataRepository::WrapperBase & wrapper ) const;
131 
132  void importSurfacicFieldOnArray( string const & faceBlockName,
133  string const & meshFieldName,
134  dataRepository::WrapperBase & wrapper ) const;
135 
140  // TODO can we use unique_ptr to hold mesh?
141  vtkSmartPointer< vtkDataSet > m_vtkMesh;
142 
144  string m_regionAttributeName;
145 
147  string m_structuredIndexAttributeName;
148 
150  string m_mainBlockName;
151 
153  string_array m_faceBlockNames;
154 
156  stdMap< string, vtkSmartPointer< vtkDataSet > > m_faceBlockMeshes;
157 
159  string_array m_nodesetNames;
160 
162  integer m_partitionRefinement = 0;
163 
165  integer m_useGlobalIds = 0;
166 
168  vtk::PartitionMethod m_partitionMethod = vtk::PartitionMethod::parmetis;
169 
171  vtk::CellMapType m_cellMap;
172 
174  string m_dataSourceName;
175 
177  VTKHierarchicalDataSource * m_dataSource;
178 
179 };
180 
181 } // namespace geos
182 
183 #endif /* GEOS_MESH_GENERATORS_VTKMESHGENERATOR_HPP */
The CellBlockManager class provides an interface to ObjectManagerBase in order to manage CellBlock da...
Base class for external mesh generators (importers).
Block
Describe which kind of block must be considered.
The VTKMeshGenerator class provides a class implementation of VTK generated meshes.
VTKMeshGenerator(const string &name, Group *const parent)
Main constructor for MeshGenerator base class.
void fillCellBlockManager(CellBlockManager &cellBlockManager, SpatialPartition &partition) override
Generate the mesh using the VTK library.
void postInputInitialization() override
void importFieldOnArray(Block block, string const &blockName, string const &meshFieldName, bool isMaterialField, dataRepository::WrapperBase &wrapper) const override
import field from the mesh on the array accessible via the given wrapper.
static string catalogName()
Return the name of the VTKMeshGenerator in object Catalog.
void freeResources() override
Free internal resources associated with mesh/data import.
Base class for all wrappers containing common operations.
Definition: WrapperBase.hpp:56
stdVector< string > string_array
A 1-dimensional array of geos::string types.
Definition: DataTypes.hpp:361
int integer
Signed integer type.
Definition: DataTypes.hpp:81
Structure to hold scoped key names.
Definition: Group.hpp:1444