NEMoSys  0.63.0
A modular, extensible resource with robust automated mesh generation, mesh quality analysis, adaptive mesh refinement, and data transfer between arbitrary meshes.
proteusHdf5.H
Go to the documentation of this file.
1 /*******************************************************************************
2 * Promesh *
3 * Copyright (C) 2022, IllinoisRocstar LLC. All rights reserved. *
4 * *
5 * Promesh is the property of IllinoisRocstar LLC. *
6 * *
7 * IllinoisRocstar LLC *
8 * Champaign, IL *
9 * www.illinoisrocstar.com *
10 * promesh@illinoisrocstar.com *
11 *******************************************************************************/
12 /*******************************************************************************
13 * This file is part of Promesh *
14 * *
15 * This version of Promesh is free software: you can redistribute it and/or *
16 * modify it under the terms of the GNU Lesser General Public License as *
17 * published by the Free Software Foundation, either version 3 of the License, *
18 * or (at your option) any later version. *
19 * *
20 * Promesh is distributed in the hope that it will be useful, but WITHOUT ANY *
21 * WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS *
22 * FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for more *
23 * details. *
24 * *
25 * You should have received a copy of the GNU Lesser General Public License *
26 * along with this program. If not, see <https://www.gnu.org/licenses/>. *
27 * *
28 *******************************************************************************/
29 #ifndef NEMOSYS_PROTEUSHDF5_H_
30 #define NEMOSYS_PROTEUSHDF5_H_
31 
32 #include <string>
33 
34 #include "IO/hdf5Reader.H"
35 #include "Mesh/meshBase.H"
36 
37 /* Special purpose class for Proteus HDF5 files */
38 /** @brief stores information for each block of Proteus data
39  **/
40 struct proteusBlock {
41  // block info
42  int numElements; /**< number of elements in block */
43  int numVerticesPerElement; /**< number of vertices per element (only one
44  element type allowed per block) */
45  std::string blockName; /**< block name */
46 
47  // element type
48  int vtkElementType; /**< VTK element type */
49  int originalVtkElementType; /**< original VTK element type. This is different
50  if we are converting a high-order element into
51  a low-order element; often useful for
52  visualization */
53 
54  // vertex info
55  std::vector<std::vector<std::vector<double>>> vertices; /**< vertices */
56  std::vector<int> loc2Glob; /**< map between local and global IDs */
57 
58  // field data
59  std::vector<std::vector<double>> elementData; /**< element field data */
60  std::vector<std::vector<double>> vertexData; /**< vertex field data */
62 };
63 
64 /** @brief stores information for each superblock of Proteus data
65  **/
67  // global mesh info
68  int numElements; /**< global number of elements */
69  int numVertices; /**< global number of vertices */
70  int maxNumVerticesPerElement; /**< maximum number of vertices for all element
71  types */
72 
73  // vertex coordinate info
74  std::vector<std::vector<double>>
75  coordinates; /**< coordinates[vertex id][dim] = coordinate */
76  std::vector<double> xCrd, yCrd,
77  zCrd; /**< coordinates separated out into x,y,z coordinates */
78 
79  // element info
80  std::vector<std::vector<int>>
81  elements; /**< connectivity info, elements[id][vertex no] = vertex id */
82  std::vector<int>
83  elementTypes; // element type info, elementType[id] = vtk element type */
84 
85  // element types list
86  std::vector<int> elementTypesList; /**< list of element types, ordered by type
87  according to elements vector */
88 
89  std::vector<int>
90  originalElementTypesList; /**< list of original element types, ordered by
91  type according to elements vector */
92 
93  // element connectivity sorted by element type and global id
94  std::vector<std::map<int, std::vector<int>>>
95  vtkConnectivity; /**< connectivity array for all element types by VTK
96  convention, vtkConnectivity[element type][global
97  element id][vertex no] = global vertex id */
98 
99  // field data
100  std::vector<std::vector<double>>
101  elementData; /**< elementData[field no][element id] = data */
102  std::vector<std::vector<double>>
103  vertexData; /**< cellData[field no][vertex id] = data */
104 
106 };
107 
108 /** @brief Class to store HDF5 data from Proteus files
109  **/
110 class proteusHdf5 : public hdf5Reader {
111 public:
112  /** @brief Construct Proteus HDF5 object
113  @param fieldFName Proteus format HDF5 field filename
114  @param meshFName Output VTK mesh filename
115  @param edgeSidesetName Name of sideset written to output Exodus file
116  @param exoMeshFName Output Exodus mesh filename
117  @param lowOrder Boolean converting high order cells to low order
118  @param bndryConst Boolean to employ boundary constraint during refinement
119  **/
120  proteusHdf5(std::string fname, std::string outputFname,
121  std::string edgeSidesetName, std::string exoMeshFName,
122  bool lowOrder = false, bool bndryConst = true);
123 
124  ~proteusHdf5();
125 
126  /** @brief Get names of Proteus vector fields
127  @param name String description of Proteus vector field names to read
128  @param stringVector Vector of Proteus vector field names
129  @param numVectors Number of vector fields
130  **/
131  void getVectorNames(std::string name, std::vector<std::string> &stringVector,
132  int numVectors);
133 
134  /** @brief Read information for each Proteus block
135  @param group HDF5 Group object containing data
136  @param myBlock Proteus block object to read
137  **/
138  void getBlockInfo(Group &group, proteusBlock &myBlock);
139 
140  /** @brief Get coordinates of Proteus block
141  @param group HDF5 Group object containing data
142  @param myBlock Proteus block object to read
143  **/
144  void getBlockXYZ(Group &group, proteusBlock &myBlock);
145 
146  /** @brief Get global coordinates of Proteus block
147  @param group HDF5 Group object containing data
148  @param myBlock Proteus block object to read
149  **/
150  void getBlockGlobalID(Group &group, proteusBlock &myBlock);
151 
152  /** @brief Get element data for block
153  @param group HDF5 Group object containing data
154  @param myBlock Proteus block object to read
155  **/
156  void getBlockElementData(Group &group, proteusBlock &myBlock);
157 
158  /** @brief Get vertex data for block
159  @param group Vector of Proteus vector field names
160  @param myBlock Number of vector fields
161  **/
162  void getBlockVertexData(Group &group, proteusBlock &myBlock);
163 
164  // access control data
165  /**@brief Get number of blocks in Proteus file
166  @return Number of blocks
167  **/
168  int getNumBlocks() const;
169 
170  /**@brief Get number of spatial dimensions in Proteus file
171  @return Number of dimensions
172  **/
173  int getNumDimensions() const;
174 
175  /**@brief Get number of vertex vectors in Proteus file
176  @return Number of vertex vectors
177  **/
178  int getNumVertexVectors() const;
179 
180  /**@brief Get number of element vectors in Proteus file
181  @return Number of element vectors
182  **/
183  int getNumElementVectors() const;
184 
185  /**@brief Get character string length of names in Proteus file
186  @return Length of string
187  **/
188  int getCharStringLength() const;
189 
190  ///** @brief Write boundary condition information to json file
191  // @param myMeshBase meshBase object of Proteus mesh
192  // @param edgeSidesetName Name of sideset
193  // **/
194  // void writeBcFile(meshBase *myMeshBase, std::string edgeSidesetName);
195 
196  /** @brief Get boundary edge points of mesh using MAdLib
197  @param startPts list containing first vertex of each boundary edge
198  @param endPts list containing second vertex of each boundary edge
199  **/
200  void getBoundaryEdgePts(vtkSmartPointer<vtkPoints> startPts,
201  vtkSmartPointer<vtkPoints> endPts);
202 
203  /** @brief Get boundary sidesets of mesh
204  @param myMeshBase meshBase object of Proteus mesh
205  @param startPts list containing first vertex of each boundary edge
206  @param endPts list containing second vertex of each boundary edge
207  @param sidesetElementIdList list containing element IDs for each edge in
208  boundary sideset
209  @param sidesetSideIdList list containing side IDs for each edge in
210  boundary sideset
211  **/
212  void getBoundarySideSets(meshBase *myMeshBase,
213  vtkSmartPointer<vtkPoints> startPts,
214  vtkSmartPointer<vtkPoints> endPts,
215  vtkSmartPointer<vtkIdList> sidesetElementIdList,
216  vtkSmartPointer<vtkIdList> sidesetSideIdList);
217 
218  /** @brief Get normal vector for 2d meshes
219  @param myMeshBase meshBase object of Proteus mesh
220  @return normal vector convention, projection onto positive z-axis:
221  positive (true) or negative (false)
222  **/
223  bool get2dCellNodeOrder(meshBase *myMeshBase);
224 
225 private:
226  /** @brief Read Proteus blocks
227  **/
228  void getBlocks();
229 
230  /** @brief Read Proteus control information
231  **/
232  void getControlInformation();
233 
234  /** @brief Merge Proteus block together
235  **/
236  void mergeBlocks();
237 
238  // management data
239 private:
240  bool lowOrder; /**< Boolean converting high order cells to low order (useful
241  for visualization)*/
242 
243  // Proteus control data
244  int numBlocks; /**< number of Proteus blocks */
245  int numDimensions; /**< number of spatial dimensions */
246  int numVertexVectors; /**< number of vertex vector fields in Proteus file */
247  int numElementVectors; /**< number of element vector fields in Proteus file
248  */
249  int charStringLength; /**< string length of each Proteus variable name */
250 
251  // field data names
252  std::vector<std::string>
253  vertexVectorNames; /**< vector containing names of vertex vectors */
254  std::vector<std::string>
255  elementVectorNames; /**< vector containing names of element vectors */
256 
258  mySuperBlock; /**< global block structure for entire Proteus mesh */
259 
260  std::vector<proteusBlock>
261  proteusBlocks; /**< vector of all individual Proteus-style blocks */
262 };
263 
264 #endif // NEMOSYS_PROTEUSHDF5_H_
int numElementVectors
number of element vector fields in Proteus file
Definition: proteusHdf5.H:247
bool lowOrder
Boolean converting high order cells to low order (useful for visualization)
Definition: proteusHdf5.H:240
std::vector< std::vector< double > > coordinates
coordinates[vertex id][dim] = coordinate
Definition: proteusHdf5.H:75
std::vector< std::vector< double > > elementData
elementData[field no][element id] = data
Definition: proteusHdf5.H:101
int originalVtkElementType
original VTK element type.
Definition: proteusHdf5.H:49
std::vector< int > elementTypes
Definition: proteusHdf5.H:83
stores information for each superblock of Proteus data
Definition: proteusHdf5.H:66
std::vector< int > originalElementTypesList
list of original element types, ordered by type according to elements vector
Definition: proteusHdf5.H:90
std::vector< std::string > elementVectorNames
vector containing names of element vectors
Definition: proteusHdf5.H:255
A brief description of meshBase.
Definition: meshBase.H:64
std::vector< std::vector< double > > elementData
element field data
Definition: proteusHdf5.H:59
int vtkElementType
VTK element type.
Definition: proteusHdf5.H:48
std::vector< double > zCrd
coordinates separated out into x,y,z coordinates
Definition: proteusHdf5.H:76
std::vector< std::vector< double > > vertexData
vertex field data
Definition: proteusHdf5.H:60
proteusSuperBlock mySuperBlock
global block structure for entire Proteus mesh
Definition: proteusHdf5.H:258
std::string blockName
block name
Definition: proteusHdf5.H:45
int numVerticesPerElement
number of vertices per element (only one element type allowed per block)
Definition: proteusHdf5.H:43
int numDimensions
number of spatial dimensions
Definition: proteusHdf5.H:245
int numElements
global number of elements
Definition: proteusHdf5.H:68
int numBlocks
number of Proteus blocks
Definition: proteusHdf5.H:244
std::vector< proteusBlock > proteusBlocks
vector of all individual Proteus-style blocks
Definition: proteusHdf5.H:261
std::vector< int > loc2Glob
map between local and global IDs
Definition: proteusHdf5.H:56
int numElements
number of elements in block
Definition: proteusHdf5.H:42
int numVertexVectors
number of vertex vector fields in Proteus file
Definition: proteusHdf5.H:246
int maxNumVerticesPerElement
maximum number of vertices for all element types
Definition: proteusHdf5.H:70
std::vector< std::vector< std::vector< double > > > vertices
vertices
Definition: proteusHdf5.H:55
Class to store HDF5 data from Proteus files.
Definition: proteusHdf5.H:110
int numVertices
global number of vertices
Definition: proteusHdf5.H:69
std::vector< std::vector< double > > vertexData
cellData[field no][vertex id] = data
Definition: proteusHdf5.H:103
std::vector< std::map< int, std::vector< int > > > vtkConnectivity
connectivity array for all element types by VTK convention, vtkConnectivity[element type][global elem...
Definition: proteusHdf5.H:95
std::vector< std::string > vertexVectorNames
vector containing names of vertex vectors
Definition: proteusHdf5.H:253
Class for reading data from HDF5 files.
Definition: hdf5Reader.H:60
std::vector< int > elementTypesList
list of element types, ordered by type according to elements vector
Definition: proteusHdf5.H:86
std::vector< std::vector< int > > elements
connectivity info, elements[id][vertex no] = vertex id
Definition: proteusHdf5.H:81
stores information for each block of Proteus data
Definition: proteusHdf5.H:40
int charStringLength
string length of each Proteus variable name
Definition: proteusHdf5.H:249