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 Class Reference

Class to store HDF5 data from Proteus files. More...

Detailed Description

Definition at line 110 of file proteusHdf5.H.

Public Member Functions

 proteusHdf5 (std::string fname, std::string outputFname, std::string edgeSidesetName, std::string exoMeshFName, bool lowOrder=false, bool bndryConst=true)
 Construct Proteus HDF5 object. More...
 
 ~proteusHdf5 ()
 
void getVectorNames (std::string name, std::vector< std::string > &stringVector, int numVectors)
 Get names of Proteus vector fields. More...
 
void getBlockInfo (Group &group, proteusBlock &myBlock)
 Read information for each Proteus block. More...
 
void getBlockXYZ (Group &group, proteusBlock &myBlock)
 Get coordinates of Proteus block. More...
 
void getBlockGlobalID (Group &group, proteusBlock &myBlock)
 Get global coordinates of Proteus block. More...
 
void getBlockElementData (Group &group, proteusBlock &myBlock)
 Get element data for block. More...
 
void getBlockVertexData (Group &group, proteusBlock &myBlock)
 Get vertex data for block. More...
 
int getNumBlocks () const
 Get number of blocks in Proteus file. More...
 
int getNumDimensions () const
 Get number of spatial dimensions in Proteus file. More...
 
int getNumVertexVectors () const
 Get number of vertex vectors in Proteus file. More...
 
int getNumElementVectors () const
 Get number of element vectors in Proteus file. More...
 
int getCharStringLength () const
 Get character string length of names in Proteus file. More...
 
void getBoundaryEdgePts (vtkSmartPointer< vtkPoints > startPts, vtkSmartPointer< vtkPoints > endPts)
 Write boundary condition information to json file. More...
 
void getBoundarySideSets (meshBase *myMeshBase, vtkSmartPointer< vtkPoints > startPts, vtkSmartPointer< vtkPoints > endPts, vtkSmartPointer< vtkIdList > sidesetElementIdList, vtkSmartPointer< vtkIdList > sidesetSideIdList)
 Get boundary sidesets of mesh. More...
 
bool get2dCellNodeOrder (meshBase *myMeshBase)
 Get normal vector for 2d meshes. More...
 
int openFile ()
 Open HDF5 file. More...
 
int closeFile ()
 Close HDF5 file. More...
 
int readTopDataSet (std::string dsetName, std::string &buffer)
 Read string data from top level HDF5 DataSet. More...
 
template<class T >
int readTopDataSet (std::string dsetName, std::vector< T > &buffer)
 Read numeric data from top level HDF5 DataSet. More...
 
int readDataSet (DataSet dset, std::string &buffer)
 Read string data from generic HDF5 DataSet. More...
 
template<class T >
int readDataSet (DataSet dset, std::vector< T > &buffer)
 Read numeric data from generic HDF5 DataSet. More...
 
int readGroup (std::string groupName, Group &group)
 Read existing HDF5 Group. More...
 
template<class T >
int readGroupDataSet (std::string dsetName, Group &group, std::vector< T > &buffer)
 Read numeric data from HDF5 group. More...
 
std::string getFileName ()
 Get HDF5 filename. More...
 
template<class T >
void flattenBuffer2d (std::vector< std::vector< T >> originalBuffer, std::vector< T > &flatBuffer)
 Flatten 2d data buffer. More...
 
template<class T >
void unflattenBuffer2d (std::vector< T > flatBuffer, std::vector< std::vector< T >> &originalBuffer)
 Unflatten 2d data buffer. More...
 
void exportToVTKMesh ()
 Export to VTK format. More...
 
vtkSmartPointer< vtkDataSet > getVTKMesh ()
 Get VTK mesh object. More...
 
void exportToMeshBase ()
 Export HDF5 mesh to meshBase object. More...
 
meshBasegetMeshBase ()
 Get meshBase object. More...
 
void writeVTK ()
 Write mesh to VTK file. More...
 
void setFields (meshBase *myMeshBase, std::vector< std::string > dataNames, std::vector< std::vector< double >> data, int pointOrCell)
 Set meshBase field data. More...
 
void setNumberOfVertices (int numVertices)
 Set number of vertices in meshBase. More...
 
void setNumberOfElements (int numElements)
 Set number of elements in meshBase. More...
 
void setNumberOfDimensions (int numDimensions)
 Set number of spatial dimensions in meshBase. More...
 
void setCoordinates (std::vector< double > xCrd)
 Set meshBase 1d coordinate data. More...
 
void setCoordinates (std::vector< double > xCrd, std::vector< double > yCrd)
 Set meshBase 2d coordinate data. More...
 
void setCoordinates (std::vector< double > xCrd, std::vector< double > yCrd, std::vector< double > zCrd)
 Set meshBase 3d coordinate data. More...
 
void setConnectivities (std::vector< std::map< int, std::vector< int >>> elementConnectivity)
 Set meshBase element connectivity table. More...
 
void setElementTypes (std::vector< int > vtkElementTypes)
 Set VTK element types for each element. More...
 
void setElementTypesList (std::vector< int > vtkElementTypesList)
 Set list of each unique type of VTK element in mesh. More...
 
std::vector< int > getElementConnectivity (int elementId)
 Get element connectivity for global element ID. More...
 

Protected Attributes

std::string h5FileName
 HDF5 filename. More...
 
std::string baseH5FileName
 base of HDF5 filename More...
 
std::string outputFname
 VTK output filename. More...
 
int verb
 verbosity flag More...
 
H5File file
 HDF5 file object. More...
 
vtkSmartPointer< vtkDataSet > vtkMesh
 VTK DataSet object. More...
 
meshBasemyMeshBase
 meshBase object More...
 
std::vector< double > xCrd
 
std::vector< double > yCrd
 
std::vector< double > zCrd
 vector of coordinates More...
 
std::vector< std::map< int, std::vector< int > > > elementConnectivity
 element connectivity table More...
 
int numVertices
 number of vertices in mesh More...
 
int numElements
 number of elements in mesh More...
 
std::vector< int > vtkElementTypes
 vector of VTK element type for each element More...
 
std::vector< int > vtkElementTypesList
 vector storing list of each unique VTK element type More...
 

Private Member Functions

void getBlocks ()
 Read Proteus blocks. More...
 
void getControlInformation ()
 Read Proteus control information. More...
 
void mergeBlocks ()
 Merge Proteus block together. More...
 

Private Attributes

bool lowOrder
 Boolean converting high order cells to low order (useful for visualization) More...
 
int numBlocks
 number of Proteus blocks More...
 
int numDimensions
 number of spatial dimensions More...
 
int numVertexVectors
 number of vertex vector fields in Proteus file More...
 
int numElementVectors
 number of element vector fields in Proteus file More...
 
int charStringLength
 string length of each Proteus variable name More...
 
std::vector< std::string > vertexVectorNames
 vector containing names of vertex vectors More...
 
std::vector< std::string > elementVectorNames
 vector containing names of element vectors More...
 
proteusSuperBlock mySuperBlock
 global block structure for entire Proteus mesh More...
 
std::vector< proteusBlockproteusBlocks
 vector of all individual Proteus-style blocks More...
 

Inherits hdf5Reader.

Constructor & Destructor Documentation

◆ proteusHdf5()

proteusHdf5::proteusHdf5 ( std::string  fname,
std::string  outputFname,
std::string  edgeSidesetName,
std::string  exoMeshFName,
bool  lowOrder = false,
bool  bndryConst = true 
)
Parameters
fieldFNameProteus format HDF5 field filename
meshFNameOutput VTK mesh filename
edgeSidesetNameName of sideset written to output Exodus file
exoMeshFNameOutput Exodus mesh filename
lowOrderBoolean converting high order cells to low order
bndryConstBoolean to employ boundary constraint during refinement

Definition at line 60 of file proteusHdf5.C.

References hdf5Reader::closeFile(), meshBase::convertQuads(), meshBase::Create(), NEM::DRV::TransferDriver::CreateTransferObject(), proteusSuperBlock::elementData, proteusSuperBlock::elementTypes, proteusSuperBlock::elementTypesList, elementVectorNames, hdf5Reader::exportToMeshBase(), hdf5Reader::exportToVTKMesh(), NEM::DRV::ConversionDriver::genExo(), get2dCellNodeOrder(), getBlocks(), getBoundaryEdgePts(), getBoundarySideSets(), getControlInformation(), hdf5Reader::getMeshBase(), getNumDimensions(), getVectorNames(), lowOrder, mergeBlocks(), hdf5Reader::myMeshBase, mySuperBlock, NEM::MSH::New(), numDimensions, proteusSuperBlock::numElements, numElementVectors, numVertexVectors, proteusSuperBlock::numVertices, hdf5Reader::openFile(), meshBase::refineMesh(), hdf5Reader::setConnectivities(), hdf5Reader::setCoordinates(), hdf5Reader::setElementTypes(), hdf5Reader::setElementTypesList(), hdf5Reader::setFields(), meshBase::setMetadata(), hdf5Reader::setNumberOfDimensions(), hdf5Reader::setNumberOfElements(), hdf5Reader::setNumberOfVertices(), proteusSuperBlock::vertexData, vertexVectorNames, proteusSuperBlock::vtkConnectivity, meshBase::write(), proteusSuperBlock::xCrd, proteusSuperBlock::yCrd, and proteusSuperBlock::zCrd.

63  : hdf5Reader(fname, meshFname) {
64  // Set boolean determining whether to output higher order elements
65  // as lower order elements; useful for visualization
66  lowOrder = _lowOrder;
67 
68  // Open HDF5 file
69  openFile();
70 
71  // Get file system metadata
73  getVectorNames("ELEMENT_VECTOR_NAMES", elementVectorNames, numElementVectors);
74  getVectorNames("VERTEX_VECTOR_NAMES", vertexVectorNames, numVertexVectors);
75 
76  // Read in each block
77  getBlocks();
78 
79  // Merge block data
80  mergeBlocks();
81 
82  // Set HDF5 data
90 
91  // Export to VTK file
93 
94  // Export to MeshBase
96 
97  // Get VTK mesh
99 
100  // Add vertex data to meshBase
102 
103  // Add element data to meshBase
105 
106  // If mesh is 2d, get node ordering convention to follow when splitting cells
107  bool rhr;
108  // true = right-hand rule (projection of normal vector onto z-axis is +)
109  // false = left-hand rule (projection of normal vector onto z-axis is -)
110  if (getNumDimensions() == 2) {
111  std::cout << "Checking 2d cell node order" << std::endl;
113  if (rhr) {
114  std::cout << "2d node ordering right-hand rule: +z" << std::endl;
115  } else {
116  std::cout << "2d node ordering right-hand rule: -z" << std::endl;
117  }
118  }
119 
120  // Create a meshbase with quads removed (converted to tris)
121  std::cout << "Removing quad elements (converting to tris) for refinement"
122  << std::endl;
123  meshBase *myMeshBaseNoQuads = myMeshBase->convertQuads();
124 
125  // Check node ordering again after element conversion (all tris)
126  if (getNumDimensions() == 2) {
127  std::cout << "Checking 2d cell node order" << std::endl;
128  rhr = get2dCellNodeOrder(myMeshBaseNoQuads);
129  if (rhr) {
130  std::cout << "2d node ordering right-hand rule: +z" << std::endl;
131  } else {
132  std::cout << "2d node ordering right-hand rule: -z" << std::endl;
133  }
134  }
135 
136  std::cout << "Writing meshBase to file" << std::endl;
137  myMeshBase->write();
138 
139  // Transfer solution onto new mesh
140  std::cout << "Performing solution transfer" << std::endl;
141  // myMeshBase->transfer(myMeshBaseNoQuads, "Consistent Interpolation");
143  myMeshBase, myMeshBaseNoQuads, "Consistent Interpolation");
144  transfer->run();
145 
146  // debug
147  // myMeshBase->write();
148  // myGmshMeshBase->write();
149 
150  // Refine mesh
151  myMeshBaseNoQuads->refineMesh("uniform", 0.5, "refined.vtu", true,
152  bndryConst);
153 
154  // Read in refined mesh
155  meshBase *refinedMeshBase = meshBase::Create("refined.vtu");
156 
157  // Write refined mesh to file
158  myMeshBaseNoQuads->write();
159  // myGmshMeshNoQuads->write();
160 
161  // Write to file
162  // writeVTK();
163 
164  // Check node ordering again after refinement
165  if (getNumDimensions() == 2) {
166  std::cout << "Checking 2d cell node order" << std::endl;
167  rhr = get2dCellNodeOrder(refinedMeshBase);
168  if (rhr) {
169  std::cout << "2d node ordering right-hand rule: +z" << std::endl;
170  } else {
171  std::cout << "2d node ordering right-hand rule: -z" << std::endl;
172  }
173  }
174 
175  // Sideset implementation
176 
177  // Got coordinates for each boundary edge (two pts)
178  vtkSmartPointer<vtkPoints> startPts = vtkSmartPointer<vtkPoints>::New();
179  vtkSmartPointer<vtkPoints> endPts = vtkSmartPointer<vtkPoints>::New();
180  getBoundaryEdgePts(startPts, endPts);
181 
182  // Initialize lists for storing sideset quantities
183  vtkSmartPointer<vtkIdList> sidesetElementIdList =
185  vtkSmartPointer<vtkIdList> sidesetSideIdList =
187  // Find elements and sides in meshbase
188  // Create two arrays:
189  // -one with element ID for sideset cells
190  // -one with side ID for each sideset cell
191  getBoundarySideSets(refinedMeshBase, startPts, endPts, sidesetElementIdList,
192  sidesetSideIdList);
193 
194  // Setting meshbase metadata
195  // Set meshbase sideset metadata
196  vtkSmartPointer<vtkModelMetadata> sidesetData =
198  vtkSmartPointer<vtkStringArray> sidesetNames =
200  // currently supports one sideset
201  sidesetNames->InsertNextValue(edgeSidesetName);
202  sidesetData->SetSideSetNames(sidesetNames);
203 
204  // Convert sideset lists into int*
205  std::vector<int> _sidesetElementIdList;
206  for (int iId = 0; iId < sidesetElementIdList->GetNumberOfIds(); iId++) {
207  _sidesetElementIdList.push_back(sidesetElementIdList->GetId(iId));
208  }
209 
210  std::vector<int> _sidesetSideIdList;
211  for (int iId = 0; iId < sidesetSideIdList->GetNumberOfIds(); iId++) {
212  _sidesetSideIdList.push_back(sidesetSideIdList->GetId(iId));
213  }
214 
215  // Set sidesets to meshbase
216  sidesetData->SetSideSetElementList(_sidesetElementIdList.data());
217  sidesetData->SetSideSetSideList(_sidesetSideIdList.data());
218  int sidesetSizes[1] = {
219  static_cast<int>(sidesetElementIdList->GetNumberOfIds())};
220 
221  sidesetData->SetSideSetSize(sidesetSizes);
222  sidesetData->SetNumberOfSideSets(1);
223  refinedMeshBase->setMetadata(sidesetData);
224 
225  refinedMeshBase->write();
226 
227  std::vector<meshBase *> mbVec;
228  mbVec.push_back(refinedMeshBase);
229  std::string exoName = exoMeshFName;
230 
231  NEM::DRV::ConversionDriver::genExo(mbVec, exoName);
232 
233  // Close HDF5 file
234  closeFile();
235 }
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
void setNumberOfVertices(int numVertices)
Set number of vertices in meshBase.
Definition: hdf5Reader.C:154
std::vector< std::vector< double > > elementData
elementData[field no][element id] = data
Definition: proteusHdf5.H:101
std::vector< double > yCrd
Definition: proteusHdf5.H:76
void refineMesh(const std::string &method, int arrayID, double dev_mult, bool maxIsmin, double edge_scale, const std::string &ofname, bool transferData, double sizeFactor=1., bool constrainBoundary=false)
perform sizefield-based h-refinement.
Definition: meshBase.C:1565
std::vector< int > elementTypes
Definition: proteusHdf5.H:83
void mergeBlocks()
Merge Proteus block together.
Definition: proteusHdf5.C:351
void setElementTypesList(std::vector< int > vtkElementTypesList)
Set list of each unique type of VTK element in mesh.
Definition: hdf5Reader.C:186
std::vector< std::string > elementVectorNames
vector containing names of element vectors
Definition: proteusHdf5.H:255
void setElementTypes(std::vector< int > vtkElementTypes)
Set VTK element types for each element.
Definition: hdf5Reader.C:182
A brief description of meshBase.
Definition: meshBase.H:64
void exportToVTKMesh()
Export to VTK format.
Definition: hdf5Reader.C:241
void getBlocks()
Read Proteus blocks.
Definition: proteusHdf5.C:308
geoMeshBase * New(MeshType meshType)
Create a new mesh object.
void getVectorNames(std::string name, std::vector< std::string > &stringVector, int numVectors)
Get names of Proteus vector fields.
Definition: proteusHdf5.C:267
int getNumDimensions() const
Get number of spatial dimensions in Proteus file.
Definition: proteusHdf5.C:1047
int closeFile()
Close HDF5 file.
Definition: hdf5Reader.C:63
meshBase * getMeshBase()
Get meshBase object.
Definition: hdf5Reader.C:321
static meshBase * Create(const std::string &fname)
Construct vtkMesh from filename.
Definition: meshBase.C:78
bool get2dCellNodeOrder(meshBase *myMeshBase)
Get normal vector for 2d meshes.
Definition: proteusHdf5.C:1055
std::vector< double > zCrd
coordinates separated out into x,y,z coordinates
Definition: proteusHdf5.H:76
void getControlInformation()
Read Proteus control information.
Definition: proteusHdf5.C:239
proteusSuperBlock mySuperBlock
global block structure for entire Proteus mesh
Definition: proteusHdf5.H:258
void getBoundarySideSets(meshBase *myMeshBase, vtkSmartPointer< vtkPoints > startPts, vtkSmartPointer< vtkPoints > endPts, vtkSmartPointer< vtkIdList > sidesetElementIdList, vtkSmartPointer< vtkIdList > sidesetSideIdList)
Get boundary sidesets of mesh.
Definition: proteusHdf5.C:963
int numDimensions
number of spatial dimensions
Definition: proteusHdf5.H:245
int numElements
global number of elements
Definition: proteusHdf5.H:68
void setConnectivities(std::vector< std::map< int, std::vector< int >>> elementConnectivity)
Set meshBase element connectivity table.
Definition: hdf5Reader.C:190
virtual void write() const
write the mesh to file named after the private var &#39;filename&#39;.
Definition: meshBase.H:598
meshBase * convertQuads()
Definition: meshBase.C:1769
void setFields(meshBase *myMeshBase, std::vector< std::string > dataNames, std::vector< std::vector< double >> data, int pointOrCell)
Set meshBase field data.
Definition: hdf5Reader.C:323
hdf5Reader(std::string fname, std::string _outputFname, int _verb=0)
Construct HDF5 reader object.
Definition: hdf5Reader.H:67
void getBoundaryEdgePts(vtkSmartPointer< vtkPoints > startPts, vtkSmartPointer< vtkPoints > endPts)
Write boundary condition information to json file.
Definition: proteusHdf5.C:924
void exportToMeshBase()
Export HDF5 mesh to meshBase object.
Definition: hdf5Reader.C:313
meshBase * myMeshBase
meshBase object
Definition: hdf5Reader.H:378
int numVertexVectors
number of vertex vector fields in Proteus file
Definition: proteusHdf5.H:246
int numVertices
global number of vertices
Definition: proteusHdf5.H:69
void setNumberOfDimensions(int numDimensions)
Set number of spatial dimensions in meshBase.
Definition: hdf5Reader.C:162
std::vector< std::vector< double > > vertexData
cellData[field no][vertex id] = data
Definition: proteusHdf5.H:103
std::vector< double > xCrd
Definition: proteusHdf5.H:76
static std::shared_ptr< TransferBase > CreateTransferObject(meshBase *srcmsh, meshBase *trgmsh, const std::string &method)
void setCoordinates(std::vector< double > xCrd)
Set meshBase 1d coordinate data.
Definition: hdf5Reader.C:166
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
void setNumberOfElements(int numElements)
Set number of elements in meshBase.
Definition: hdf5Reader.C:158
void setMetadata(vtkSmartPointer< vtkModelMetadata > _metadata)
Definition: meshBase.H:440
std::vector< int > elementTypesList
list of element types, ordered by type according to elements vector
Definition: proteusHdf5.H:86
int openFile()
Open HDF5 file.
Definition: hdf5Reader.C:41
static void genExo(std::vector< meshBase *> meshes, const std::string &fname)

◆ ~proteusHdf5()

proteusHdf5::~proteusHdf5 ( )

Definition at line 237 of file proteusHdf5.C.

237 {}

Member Function Documentation

◆ closeFile()

int hdf5Reader::closeFile ( )
inherited
Returns
Integer flag: 0 (failure) or 1 (success)

Definition at line 63 of file hdf5Reader.C.

References hdf5Reader::file, and hdf5Reader::h5FileName.

Referenced by proteusHdf5().

63  {
64  std::cout << "Closing HDF5 file " << h5FileName << std::endl;
65 
66  try {
67  // Close HDF5 file
68  file.close();
69  } catch (FileIException error) {
70 #if H5_VERS_MAJOR >= 1 && H5_VERS_MINOR >= 8 && H5_VERS_RELEASE >= 20
71  error.printErrorStack();
72 #elif H5_VERS_MAJOR >= 1 && H5_VERS_MINOR >= 10 && H5_VERS_RELEASE >= 2
73  error.printErrorStack();
74 #elif H5_VERS_MAJOR >= 1 && H5_VERS_MINOR >= 11 && H5_VERS_RELEASE >= 0
75  error.printErrorStack();
76 #else
77  error.printError();
78 #endif
79  return -1;
80  }
81  return 0;
82 }
H5File file
HDF5 file object.
Definition: hdf5Reader.H:376
std::string h5FileName
HDF5 filename.
Definition: hdf5Reader.H:372

◆ exportToMeshBase()

void hdf5Reader::exportToMeshBase ( )
inherited

Definition at line 313 of file hdf5Reader.C.

References meshBase::Create(), hdf5Reader::exportToVTKMesh(), hdf5Reader::myMeshBase, and hdf5Reader::outputFname.

Referenced by proteusHdf5(), and hdf5Reader::writeVTK().

313  {
314  std::cout << "Exporting hdf5Reader object to meshBase format..." << std::endl;
315  if (!vtkMesh) {
316  exportToVTKMesh();
317  }
319 }
void exportToVTKMesh()
Export to VTK format.
Definition: hdf5Reader.C:241
static meshBase * Create(const std::string &fname)
Construct vtkMesh from filename.
Definition: meshBase.C:78
meshBase * myMeshBase
meshBase object
Definition: hdf5Reader.H:378
std::string outputFname
VTK output filename.
Definition: hdf5Reader.H:374

◆ exportToVTKMesh()

void hdf5Reader::exportToVTKMesh ( )
inherited

Definition at line 241 of file hdf5Reader.C.

References hdf5Reader::getElementConnectivity(), NEM::MSH::New(), hdf5Reader::numDimensions, hdf5Reader::numElements, hdf5Reader::numVertices, points, hdf5Reader::vtkElementTypes, hdf5Reader::xCrd, hdf5Reader::yCrd, and hdf5Reader::zCrd.

Referenced by hdf5Reader::exportToMeshBase(), hdf5Reader::getVTKMesh(), and proteusHdf5().

241  {
242  std::cout << "Exporting hdf5Reader object to VTK format..." << std::endl;
243 
244  if (!vtkMesh) {
245  // declare vtk dataset
246  auto dataSet_tmp = vtkSmartPointer<vtkUnstructuredGrid>::New();
247 
248  // points to be pushed into dataSet
250  // set double precision points
251  points->SetDataTypeToDouble();
252 
253  // allocate size for vtk point container
254  points->SetNumberOfPoints(numVertices);
255 
256  if (numDimensions == 2) {
257  for (int i = 0; i < numVertices; ++i) {
258  points->SetPoint(i, xCrd[i], yCrd[i], 0.0);
259  }
260  } else {
261  for (int i = 0; i < numVertices; ++i) {
262  points->SetPoint(i, xCrd[i], yCrd[i], zCrd[i]);
263  }
264  }
265  // add points to vtk mesh data structure
266  std::cout << " setting VTK point data..." << std::endl;
267  dataSet_tmp->SetPoints(points);
268  // allocate space for elements
269  std::cout << " allocating memory for elements..." << std::endl;
270  dataSet_tmp->Allocate(numElements);
271  // add the elements
272  std::cout << " adding VTK elements..." << std::endl;
273  for (int i = 0; i < numElements; ++i) {
274  auto vtkElmIds = vtkSmartPointer<vtkIdList>::New();
275  std::vector<int> elmIds(getElementConnectivity(i));
276  vtkElmIds->SetNumberOfIds(elmIds.size());
277  for (int j = 0; j < elmIds.size(); ++j) {
278  vtkElmIds->SetId(j, elmIds[j] - 0);
279  }
280  switch (vtkElementTypes[i]) {
281 #if VTK_MAJOR_VERSION > 8 || (VTK_MAJOR_VERSION == 8 && VTK_MINOR_VERSION >= 1)
282  case VTK_LAGRANGE_QUADRILATERAL:
283  dataSet_tmp->InsertNextCell(VTK_LAGRANGE_QUADRILATERAL, vtkElmIds);
284  break;
285 #endif
286  case VTK_QUAD:
287  dataSet_tmp->InsertNextCell(VTK_QUAD, vtkElmIds);
288  break;
289  case VTK_TRIANGLE:
290  dataSet_tmp->InsertNextCell(VTK_TRIANGLE, vtkElmIds);
291  break;
292  default:
293  std::cerr << "Unknown element type " << vtkElementTypes[i]
294  << std::endl;
295  break;
296  }
297  }
298  // build links - element incidences
299  dataSet_tmp->BuildLinks();
300  vtkMesh = dataSet_tmp;
301  }
302 }
geoMeshBase * New(MeshType meshType)
Create a new mesh object.
std::vector< double > yCrd
Definition: hdf5Reader.H:381
std::vector< double > xCrd
Definition: hdf5Reader.H:381
int numVertices
number of vertices in mesh
Definition: hdf5Reader.H:384
std::vector< vtkIdType > points
points given by id in .inp file
Definition: inpGeoMesh.C:133
int numDimensions
number of spatial dimensions
Definition: hdf5Reader.H:386
std::vector< double > zCrd
vector of coordinates
Definition: hdf5Reader.H:381
std::vector< int > getElementConnectivity(int elementId)
Get element connectivity for global element ID.
Definition: hdf5Reader.C:195
int numElements
number of elements in mesh
Definition: hdf5Reader.H:385
std::vector< int > vtkElementTypes
vector of VTK element type for each element
Definition: hdf5Reader.H:388

◆ flattenBuffer2d()

template<class T >
void hdf5Reader::flattenBuffer2d ( std::vector< std::vector< T >>  originalBuffer,
std::vector< T > &  flatBuffer 
)
inlineinherited
Parameters
originalBufferOriginal unflattened data
flatBufferOutput flattened data

Definition at line 241 of file hdf5Reader.H.

242  {
243  // Populate flat buffer
244  for (auto bufItr1 = originalBuffer.begin(); bufItr1 != originalBuffer.end();
245  bufItr1++) {
246  for (auto bufItr2 = (*bufItr1).begin(); bufItr2 != (*bufItr1).end();
247  bufItr2++) {
248  flatBuffer.push_back(*bufItr2);
249  }
250  }
251  }

◆ get2dCellNodeOrder()

bool proteusHdf5::get2dCellNodeOrder ( meshBase myMeshBase)
Parameters
myMeshBasemeshBase object of Proteus mesh
Returns
normal vector convention, projection onto positive z-axis: positive (true) or negative (false)

Definition at line 1055 of file proteusHdf5.C.

References meshBase::getDataSet().

Referenced by proteusHdf5().

1055  {
1056  std::map<int, bool> rhr;
1057  for (int iCell = 0; iCell < myMeshBase->getDataSet()->GetNumberOfCells();
1058  iCell++) {
1059  // Get first cell
1060  vtkSmartPointer<vtkCell> testCell =
1061  myMeshBase->getDataSet()->GetCell(iCell);
1062  vtkSmartPointer<vtkPoints> cellPts = testCell->GetPoints();
1063  double pt0[3];
1064  double pt1[3];
1065  double pt2[3];
1066  double vc0[3];
1067  double vc1[3];
1068  double vc2[3];
1069  double zNorm[3] = {0, 0, 1};
1070  double res;
1071  if (testCell->GetCellType() == VTK_TRIANGLE) {
1072  cellPts->GetPoint(0, pt0);
1073  cellPts->GetPoint(1, pt1);
1074  cellPts->GetPoint(2, pt2);
1075  vtkMath::Subtract(pt1, pt0, vc0);
1076  vtkMath::Subtract(pt2, pt0, vc1);
1077  } else if (testCell->GetCellType() == VTK_QUAD) {
1078  cellPts->GetPoint(0, pt0);
1079  cellPts->GetPoint(1, pt1);
1080  cellPts->GetPoint(2, pt2);
1081  vtkMath::Subtract(pt1, pt0, vc0);
1082  vtkMath::Subtract(pt2, pt0, vc1);
1083  } else if (testCell->GetCellType() != VTK_QUAD &&
1084  testCell->GetCellType() != VTK_TRIANGLE) {
1085  std::cerr << "Only triangular and quadrilateral elements supported."
1086  << std::endl;
1087  exit(-1);
1088  } else {
1089  continue;
1090  }
1091  vtkMath::Cross(vc0, vc1, vc2);
1092  res = vtkMath::Dot(vc2, zNorm);
1093 
1094  bool ans;
1095  if (res > 0) {
1096  ans = true;
1097  } else {
1098  ans = false;
1099  }
1100  if ((testCell->GetCellType() == VTK_QUAD &&
1101  rhr.find(VTK_QUAD) == rhr.end()) ||
1102  (testCell->GetCellType() == VTK_TRIANGLE &&
1103  rhr.find(VTK_TRIANGLE) == rhr.end())) {
1104  rhr[testCell->GetCellType()] = ans;
1105  } else {
1106  if (rhr[testCell->GetCellType()] != ans) {
1107  std::cout << "Mismatched ordering for type " << testCell->GetCellType()
1108  << ", cell number " << iCell << std::endl;
1109  }
1110  }
1111  }
1112  bool conv = false;
1113  for (auto itr = rhr.begin(); itr != rhr.end(); itr++) {
1114  if (itr == rhr.begin()) {
1115  conv = itr->second;
1116  } else {
1117  if (conv != itr->second) {
1118  std::cerr << "Error: Mesh contains multiple element types with "
1119  "different node ordering."
1120  << std::endl;
1121  exit(-1);
1122  }
1123  }
1124  }
1125  return conv;
1126 }
vtkSmartPointer< vtkDataSet > getDataSet() const
get this meshes&#39; dataSet
Definition: meshBase.H:308

◆ getBlockElementData()

void proteusHdf5::getBlockElementData ( Group &  group,
proteusBlock myBlock 
)
Parameters
groupHDF5 Group object containing data
myBlockProteus block object to read

Definition at line 782 of file proteusHdf5.C.

References proteusBlock::elementData, elementVectorNames, proteusBlock::numElements, numElementVectors, hdf5Reader::readGroupDataSet(), and hdf5Reader::unflattenBuffer2d().

Referenced by getBlocks().

782  {
783  // Define GlobalID buffer
784  int bufDim1 = numElementVectors;
785  int bufDim2 = myBlock.numElements;
786  int bufDimFlat = bufDim1 * bufDim2;
787 
788  // Buffer for HDF read (float)
789  std::vector<float> flatBuffer(bufDimFlat, 0.0);
790 
791  // Buffer for proteusHDF5 object (2d double)
792  std::vector<std::vector<double>> buffer(bufDim1,
793  std::vector<double>(bufDim2, 0.0));
794 
795  // Read dataset from group
796  readGroupDataSet("ELEMENTDATA", group, flatBuffer);
797 
798  // Convert to double
799  std::vector<double> doubleFlatBuffer(flatBuffer.begin(), flatBuffer.end());
800 
801  // Unflatten buffer to 2d
802  unflattenBuffer2d(doubleFlatBuffer, buffer);
803  // Todo: Make this a generic templated function like flattenBuffer() above
804  // auto resizedBuffer = unflattenBuffer(flatBuffer);
805 
806  // Output info to screen
807  std::cout << " Element Data info:" << std::endl;
808  std::cout << " Read in " << numElementVectors
809  << " element data:" << std::endl;
810  for (auto edItr = elementVectorNames.begin();
811  edItr != elementVectorNames.end(); edItr++) {
812  std::cout << " " << *edItr << std::endl;
813  }
814 
815  // Assign element data to block
816  myBlock.elementData = buffer;
817 }
int numElementVectors
number of element vector fields in Proteus file
Definition: proteusHdf5.H:247
void unflattenBuffer2d(std::vector< T > flatBuffer, std::vector< std::vector< T >> &originalBuffer)
Unflatten 2d data buffer.
Definition: hdf5Reader.H:258
std::vector< std::string > elementVectorNames
vector containing names of element vectors
Definition: proteusHdf5.H:255
std::vector< std::vector< double > > elementData
element field data
Definition: proteusHdf5.H:59
int numElements
number of elements in block
Definition: proteusHdf5.H:42
int readGroupDataSet(std::string dsetName, Group &group, std::vector< T > &buffer)
Read numeric data from HDF5 group.
Definition: hdf5Reader.H:190

◆ getBlockGlobalID()

void proteusHdf5::getBlockGlobalID ( Group &  group,
proteusBlock myBlock 
)
Parameters
groupHDF5 Group object containing data
myBlockProteus block object to read

Definition at line 758 of file proteusHdf5.C.

References proteusBlock::loc2Glob, proteusBlock::numElements, proteusBlock::numVerticesPerElement, and hdf5Reader::readGroupDataSet().

Referenced by getBlocks().

758  {
759  // Define GlobalID buffer
760  int bufSize = myBlock.numElements * myBlock.numVerticesPerElement;
761  std::vector<int> buffer(bufSize, 0);
762 
763  // Read dataset from group
764  readGroupDataSet("GLOBALID", group, buffer);
765 
766  // Output info to screen
767  std::cout << " Global ID info:" << std::endl;
768  std::cout << " Read in "
769  << myBlock.numElements * myBlock.numVerticesPerElement
770  << " global IDs" << std::endl;
771 
772  // Create a map from local to global ID
773  std::vector<int> tmpVec;
774  for (int iBuf = 0; iBuf < bufSize; iBuf++) {
775  tmpVec.push_back(buffer[iBuf]);
776  }
777 
778  // Assign local to global ID map to block
779  myBlock.loc2Glob = tmpVec;
780 }
int numVerticesPerElement
number of vertices per element (only one element type allowed per block)
Definition: proteusHdf5.H:43
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 readGroupDataSet(std::string dsetName, Group &group, std::vector< T > &buffer)
Read numeric data from HDF5 group.
Definition: hdf5Reader.H:190

◆ getBlockInfo()

void proteusHdf5::getBlockInfo ( Group &  group,
proteusBlock myBlock 
)
Parameters
groupHDF5 Group object containing data
myBlockProteus block object to read

Definition at line 657 of file proteusHdf5.C.

References lowOrder, proteusBlock::numElements, proteusBlock::numVerticesPerElement, proteusBlock::originalVtkElementType, hdf5Reader::readGroupDataSet(), and proteusBlock::vtkElementType.

Referenced by getBlocks().

657  {
658  int bufSize = 3;
659  std::vector<int> buffer(bufSize, 0);
660  readGroupDataSet("INFO", group, buffer);
661 
662  std::vector<int> tmpVec;
663  for (int iBuf = 0; iBuf < bufSize; iBuf++) {
664  tmpVec.push_back(buffer[iBuf]);
665  }
666  myBlock.numElements = tmpVec[0];
667  myBlock.numVerticesPerElement = tmpVec[1];
668 
669  // parse element type
670  // Note: the element ID numbers below don't match up with the PROTEUS
671  // documentation, but have been implemented based on PROTEUS HDF5 test files
672  // provided
673  if (tmpVec[2] == 162) {
674  if (lowOrder) {
675 #if VTK_MAJOR_VERSION > 8 || (VTK_MAJOR_VERSION == 8 && VTK_MINOR_VERSION >= 1)
676  myBlock.originalVtkElementType = VTK_LAGRANGE_QUADRILATERAL;
677 #else
678  myBlock.originalVtkElementType = -1; // older VTK versions don't support
679  // LAGRANGE_QUADRILATERAL types
680 #endif
681  myBlock.vtkElementType = VTK_QUAD;
682  } else {
683 #if VTK_MAJOR_VERSION > 8 || (VTK_MAJOR_VERSION == 8 && VTK_MINOR_VERSION >= 1)
684  myBlock.originalVtkElementType = VTK_LAGRANGE_QUADRILATERAL;
685  myBlock.vtkElementType = VTK_LAGRANGE_QUADRILATERAL;
686 #else
687  myBlock.originalVtkElementType = VTK_QUAD;
688  myBlock.vtkElementType = VTK_QUAD;
689  std::cout << "Warning: VTK version " << vtkVersion::GetVTKSourceVersion()
690  << " does not support Lagrange Quadrilateral cells."
691  << std::endl;
692 #endif
693  }
694  } else if (tmpVec[2] == 160) {
695  myBlock.originalVtkElementType = VTK_QUAD;
696  myBlock.vtkElementType = VTK_QUAD;
697  }
698  // else if (tmpVec[2] == 110)
699  // It looks like some of the triangle elements in Proteus test files have an
700  // element type of 101 instead of 110, as documented in the user manual. I
701  // have just added this as an additional element type to catch this.
702  else if (tmpVec[2] == 110) {
703  myBlock.originalVtkElementType = VTK_TRIANGLE;
704  myBlock.vtkElementType = VTK_TRIANGLE;
705  } else {
706  std::cerr << "Element type not supported: " << tmpVec[2] << std::endl;
707  exit(-1);
708  }
709 
710  // Output info to screen
711  std::cout << " Info:" << std::endl;
712  std::cout << " Number of elements: " << tmpVec[0] << std::endl;
713  std::cout << " Number of vertices per element: "
714  << myBlock.numVerticesPerElement << std::endl;
715  std::cout << " VTK element type: " << myBlock.vtkElementType
716  << std::endl;
717 }
bool lowOrder
Boolean converting high order cells to low order (useful for visualization)
Definition: proteusHdf5.H:240
int originalVtkElementType
original VTK element type.
Definition: proteusHdf5.H:49
int vtkElementType
VTK element type.
Definition: proteusHdf5.H:48
int numVerticesPerElement
number of vertices per element (only one element type allowed per block)
Definition: proteusHdf5.H:43
int numElements
number of elements in block
Definition: proteusHdf5.H:42
int readGroupDataSet(std::string dsetName, Group &group, std::vector< T > &buffer)
Read numeric data from HDF5 group.
Definition: hdf5Reader.H:190

◆ getBlocks()

void proteusHdf5::getBlocks ( )
private

Definition at line 308 of file proteusHdf5.C.

References proteusBlock::blockName, getBlockElementData(), getBlockGlobalID(), getBlockInfo(), getBlockVertexData(), getBlockXYZ(), numBlocks, proteusBlocks, and hdf5Reader::readGroup().

Referenced by proteusHdf5().

308  {
309  std::cout << "Reading in PROTEUS block data..." << std::endl;
310 
311  for (int iBlock = 0; iBlock < numBlocks; iBlock++) {
312  // Set PROTEUS convention block name
313  std::string blockName = "BLOCK";
314  std::string numStr = std::to_string(iBlock + 1); // Proteus blocks 1-indexed
315  while (numStr.length() < 12) {
316  numStr.insert(0, "0");
317  }
318  blockName += numStr;
319  std::cout << "Reading in " + blockName + "..." << std::endl;
320 
321  // Create HDF5 Group object for block
322  Group group;
323  readGroup(blockName, group);
324 
325  // Create block object
326  proteusBlock myBlock = proteusBlock();
327  myBlock.blockName = blockName;
328 
329  // Read block info
330  getBlockInfo(group, myBlock);
331 
332  // Read coordinates
333  getBlockXYZ(group, myBlock);
334 
335  // Read global IDs
336  getBlockGlobalID(group, myBlock);
337 
338  // Read element data
339  getBlockElementData(group, myBlock);
340 
341  // Read vertex data
342  getBlockVertexData(group, myBlock);
343 
344  // Push back proteus block into data structure
345  proteusBlocks.push_back(myBlock);
346  }
347 
348  std::cout << "Finished reading PROTEUS block data..." << std::endl;
349 }
void getBlockVertexData(Group &group, proteusBlock &myBlock)
Get vertex data for block.
Definition: proteusHdf5.C:819
int readGroup(std::string groupName, Group &group)
Read existing HDF5 Group.
Definition: hdf5Reader.C:135
std::string blockName
block name
Definition: proteusHdf5.H:45
void getBlockInfo(Group &group, proteusBlock &myBlock)
Read information for each Proteus block.
Definition: proteusHdf5.C:657
int numBlocks
number of Proteus blocks
Definition: proteusHdf5.H:244
void getBlockGlobalID(Group &group, proteusBlock &myBlock)
Get global coordinates of Proteus block.
Definition: proteusHdf5.C:758
void getBlockElementData(Group &group, proteusBlock &myBlock)
Get element data for block.
Definition: proteusHdf5.C:782
std::vector< proteusBlock > proteusBlocks
vector of all individual Proteus-style blocks
Definition: proteusHdf5.H:261
void getBlockXYZ(Group &group, proteusBlock &myBlock)
Get coordinates of Proteus block.
Definition: proteusHdf5.C:719
stores information for each block of Proteus data
Definition: proteusHdf5.H:40

◆ getBlockVertexData()

void proteusHdf5::getBlockVertexData ( Group &  group,
proteusBlock myBlock 
)
Parameters
groupVector of Proteus vector field names
myBlockNumber of vector fields

Definition at line 819 of file proteusHdf5.C.

References proteusBlock::numElements, numVertexVectors, proteusBlock::numVerticesPerElement, hdf5Reader::readGroupDataSet(), hdf5Reader::unflattenBuffer2d(), proteusBlock::vertexData, and vertexVectorNames.

Referenced by getBlocks().

819  {
820  // Define GlobalID buffer
821  int bufDim1 = numVertexVectors;
822  int bufDim2 = myBlock.numElements * myBlock.numVerticesPerElement;
823  int bufDimFlat = bufDim1 * bufDim2;
824 
825  // Buffer for HDF read (float)
826  std::vector<float> flatBuffer(bufDimFlat, 0.0);
827 
828  // Buffer for proteusHDF5 object (2d double)
829  std::vector<std::vector<double>> buffer(bufDim1,
830  std::vector<double>(bufDim2, 0.0));
831 
832  // Read dataset from group
833  readGroupDataSet("VERTEXDATA", group, flatBuffer);
834 
835  // Convert to double
836  std::vector<double> doubleFlatBuffer(flatBuffer.begin(), flatBuffer.end());
837 
838  // Unflatten buffer to 2d
839  unflattenBuffer2d(doubleFlatBuffer, buffer);
840 
841  // Output info to screen
842  std::cout << " Vertex Data info:" << std::endl;
843  std::cout << " Read in " << numVertexVectors
844  << " vertex data:" << std::endl;
845  for (auto vdItr = vertexVectorNames.begin(); vdItr != vertexVectorNames.end();
846  vdItr++) {
847  std::cout << " " << *vdItr << std::endl;
848  }
849 
850  // Assign vertex data to block
851  myBlock.vertexData = buffer;
852 }
void unflattenBuffer2d(std::vector< T > flatBuffer, std::vector< std::vector< T >> &originalBuffer)
Unflatten 2d data buffer.
Definition: hdf5Reader.H:258
std::vector< std::vector< double > > vertexData
vertex field data
Definition: proteusHdf5.H:60
int numVerticesPerElement
number of vertices per element (only one element type allowed per block)
Definition: proteusHdf5.H:43
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 readGroupDataSet(std::string dsetName, Group &group, std::vector< T > &buffer)
Read numeric data from HDF5 group.
Definition: hdf5Reader.H:190
std::vector< std::string > vertexVectorNames
vector containing names of vertex vectors
Definition: proteusHdf5.H:253

◆ getBlockXYZ()

void proteusHdf5::getBlockXYZ ( Group &  group,
proteusBlock myBlock 
)
Parameters
groupHDF5 Group object containing data
myBlockProteus block object to read

Definition at line 719 of file proteusHdf5.C.

References numDimensions, proteusBlock::numElements, proteusBlock::numVerticesPerElement, hdf5Reader::readGroupDataSet(), and proteusBlock::vertices.

Referenced by getBlocks().

719  {
720  // Define xyz buffer
721  int bufSize =
722  numDimensions * myBlock.numElements * myBlock.numVerticesPerElement;
723  std::vector<float> buffer(bufSize, 0);
724 
725  // Read dataset from group
726  readGroupDataSet("XYZ", group, buffer);
727 
728  // Output info to screen
729  std::cout << " Coordinate info:" << std::endl;
730  std::cout << " Read in " << bufSize << " vertex coordinates"
731  << std::endl;
732  for (int iDim = 0; iDim < numDimensions; iDim++) {
733  std::cout << " " << myBlock.numElements * myBlock.numVerticesPerElement
734  << " in dimension " << iDim << std::endl;
735  }
736 
737  // Convert vertex coordinates into an un-flattened array using row-major
738  // ordering
739  std::vector<std::vector<std::vector<double>>> tmpVec(
740  myBlock.numElements,
741  std::vector<std::vector<double>>(
742  numDimensions,
743  std::vector<double>(myBlock.numVerticesPerElement, 0.0)));
744  int tmpInd = 0;
745  for (int iElem = 0; iElem < myBlock.numElements; iElem++) {
746  for (int iDim = 0; iDim < numDimensions; iDim++) {
747  for (int iVert = 0; iVert < myBlock.numVerticesPerElement; iVert++) {
748  tmpVec[iElem][iDim][iVert] = double(buffer[tmpInd]);
749  tmpInd++;
750  }
751  }
752  }
753 
754  // Assign vertex coordinates to block
755  myBlock.vertices = tmpVec;
756 }
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
number of elements in block
Definition: proteusHdf5.H:42
std::vector< std::vector< std::vector< double > > > vertices
vertices
Definition: proteusHdf5.H:55
int readGroupDataSet(std::string dsetName, Group &group, std::vector< T > &buffer)
Read numeric data from HDF5 group.
Definition: hdf5Reader.H:190

◆ getBoundaryEdgePts()

void proteusHdf5::getBoundaryEdgePts ( vtkSmartPointer< vtkPoints >  startPts,
vtkSmartPointer< vtkPoints >  endPts 
)
Parameters
myMeshBasemeshBase object of Proteus mesh
edgeSidesetNameName of sidesetGet boundary edge points of mesh using MAdLib
startPtslist containing first vertex of each boundary edge
endPtslist containing second vertex of each boundary edge

Definition at line 924 of file proteusHdf5.C.

References numDimensions.

Referenced by proteusHdf5().

925  {
926  if (numDimensions != 2) {
927  std::cerr << "Error: Sidesets only supports for 2D meshes" << std::endl;
928  exit(-1);
929  }
930 
931  MAd::pGModel tmpMdl = NULL;
932  MAd::GM_create(&tmpMdl, "");
933  MAd::pMesh refinedMadMesh = MAd::M_new(tmpMdl);
934  MAd::M_load(refinedMadMesh, "refined.msh");
935 
936  MAd::GM_create(&tmpMdl, "");
937  MAd::pMesh skinnedMadMesh = MAd::M_new(tmpMdl);
938 
939  std::vector<int> regIds;
940  refinedMadMesh->skin_me(skinnedMadMesh, regIds);
941 
942  MAd::EIter eit = MAd::M_edgeIter(skinnedMadMesh);
943  MAd::pEdge pe;
944  int edgeId = 0;
945  double coord[3];
946  while ((pe = MAd::EIter_next(eit))) {
947  MAd::pVertex v = MAd::E_vertex(pe, 0);
948  coord[0] = v->X;
949  coord[1] = v->Y;
950  coord[2] = v->Z;
951  startPts->InsertPoint(edgeId, coord);
952 
953  v = MAd::E_vertex(pe, 1);
954  coord[0] = v->X;
955  coord[1] = v->Y;
956  coord[2] = v->Z;
957  endPts->InsertPoint(edgeId, coord);
958 
959  edgeId++;
960  }
961 }
int numDimensions
number of spatial dimensions
Definition: proteusHdf5.H:245

◆ getBoundarySideSets()

void proteusHdf5::getBoundarySideSets ( meshBase myMeshBase,
vtkSmartPointer< vtkPoints >  startPts,
vtkSmartPointer< vtkPoints >  endPts,
vtkSmartPointer< vtkIdList >  sidesetElementIdList,
vtkSmartPointer< vtkIdList >  sidesetSideIdList 
)
Parameters
myMeshBasemeshBase object of Proteus mesh
startPtslist containing first vertex of each boundary edge
endPtslist containing second vertex of each boundary edge
sidesetElementIdListlist containing element IDs for each edge in boundary sideset
sidesetSideIdListlist containing side IDs for each edge in boundary sideset

Definition at line 963 of file proteusHdf5.C.

References meshBase::getDataSet(), and NEM::MSH::New().

Referenced by proteusHdf5().

967  {
968  // Build up locator
969  vtkSmartPointer<vtkPointLocator> pointLocator =
971  pointLocator->SetDataSet(myMeshBase->getDataSet());
972  pointLocator->BuildLocator();
973 
974  // edge Pt iterators
975  vtkIdType id1, id2;
976  vtkSmartPointer<vtkIdList> cellList1 = vtkSmartPointer<vtkIdList>::New();
977  vtkSmartPointer<vtkIdList> cellList2 = vtkSmartPointer<vtkIdList>::New();
978  vtkCell *foundCell;
979  vtkSmartPointer<vtkIdList> edgeIdListCompare =
981  vtkSmartPointer<vtkIdList> edgeIdListOriginal =
983  vtkSmartPointer<vtkIdList> edgeIdListIntersection =
985  for (int iBndEdge = 0; iBndEdge < startPts->GetNumberOfPoints(); iBndEdge++) {
986  // std::cout << "boundary edge: " << iBndEdge << std::endl;
987  edgeIdListOriginal->Reset();
988 
989  id1 = pointLocator->FindClosestPoint(startPts->GetPoint(iBndEdge));
990  id2 = pointLocator->FindClosestPoint(endPts->GetPoint(iBndEdge));
991 
992  edgeIdListOriginal->InsertNextId(id1);
993  edgeIdListOriginal->InsertNextId(id2);
994 
995  myMeshBase->getDataSet()->GetPointCells(id1, cellList1);
996  myMeshBase->getDataSet()->GetPointCells(id2, cellList2);
997 
998  cellList1->IntersectWith(cellList2);
999  if (cellList1->GetNumberOfIds() != 1) {
1000  std::cerr << "Error: found multiple cells for boundary edge in a 2d mesh"
1001  << std::endl;
1002  exit(-1);
1003  } else {
1004  sidesetElementIdList->InsertNextId(cellList1->GetId(0));
1005  // std::cout << " elem ID = " << cellList1->GetId(0) << std::endl;
1006  foundCell = myMeshBase->getDataSet()->GetCell(cellList1->GetId(0));
1007  vtkSmartPointer<vtkIdList> cellPts = vtkSmartPointer<vtkIdList>::New();
1008  myMeshBase->getDataSet()->GetCellPoints(cellList1->GetId(0), cellPts);
1009  for (int ipt = 0; ipt < cellPts->GetNumberOfIds(); ipt++) {
1010  double position[3];
1011  myMeshBase->getDataSet()->GetPoint(cellPts->GetId(ipt), position);
1012  // std::cout << " ipt " << ipt << " position = " << position[0] << ",
1013  // "
1014  // << position[1] << ", " << position[2] << std::endl;
1015  }
1016  }
1017 
1018  for (int iCellEdge = 0; iCellEdge < foundCell->GetNumberOfEdges();
1019  iCellEdge++) {
1020  edgeIdListCompare = foundCell->GetEdge(iCellEdge)->GetPointIds();
1021 
1022  // Alternate 'intersection' implementation due to pointer issues
1023  edgeIdListIntersection->Reset();
1024  for (int iId = 0; iId < edgeIdListOriginal->GetNumberOfIds(); iId++) {
1025  for (int jId = 0; jId < edgeIdListCompare->GetNumberOfIds(); jId++) {
1026  if (edgeIdListOriginal->GetId(iId) == edgeIdListCompare->GetId(jId)) {
1027  edgeIdListIntersection->InsertNextId(
1028  edgeIdListOriginal->GetId(iId));
1029  }
1030  }
1031  }
1032 
1033  if (edgeIdListIntersection->GetNumberOfIds() == 2) {
1034  sidesetSideIdList->InsertNextId(iCellEdge);
1035  // std::cout << " side ID = " << iCellEdge << std::endl;
1036  break;
1037  }
1038  if (iCellEdge == foundCell->GetNumberOfEdges() - 1) {
1039  std::cout << "Couldn't find a matching side!" << std::endl;
1040  }
1041  }
1042  }
1043 }
geoMeshBase * New(MeshType meshType)
Create a new mesh object.
vtkSmartPointer< vtkDataSet > getDataSet() const
get this meshes&#39; dataSet
Definition: meshBase.H:308

◆ getCharStringLength()

int proteusHdf5::getCharStringLength ( ) const
Returns
Length of string

Definition at line 1053 of file proteusHdf5.C.

References charStringLength.

1053 { return (charStringLength); }
int charStringLength
string length of each Proteus variable name
Definition: proteusHdf5.H:249

◆ getControlInformation()

void proteusHdf5::getControlInformation ( )
private

Definition at line 239 of file proteusHdf5.C.

References charStringLength, numBlocks, numDimensions, numElementVectors, numVertexVectors, and hdf5Reader::readTopDataSet().

Referenced by proteusHdf5().

239  {
240  std::cout << "Reading in PROTEUS control data..." << std::endl;
241 
242  // Set control buffer information
243  H5std_string dsetName("CONTROL");
244  int bufSize = 5;
245  std::vector<int> buffer(bufSize, 0);
246 
247  // Read dataset
248  readTopDataSet(dsetName, buffer);
249 
250  // Set control data
251  numBlocks = buffer[0];
252  numDimensions = buffer[1];
253  numVertexVectors = buffer[2];
254  numElementVectors = buffer[3];
255  charStringLength = buffer[4];
256 
257  // Output to info to screen
258  std::cout << " Number of Blocks: " << numBlocks << std::endl;
259  std::cout << " Number of Dimensions: " << numDimensions << std::endl;
260  std::cout << " Number of Vertex Vectors: " << numVertexVectors << std::endl;
261  std::cout << " Number of Element Vectors: " << numElementVectors
262  << std::endl;
263  std::cout << " Length of Character String: " << charStringLength
264  << std::endl;
265 }
int numElementVectors
number of element vector fields in Proteus file
Definition: proteusHdf5.H:247
int numDimensions
number of spatial dimensions
Definition: proteusHdf5.H:245
int numBlocks
number of Proteus blocks
Definition: proteusHdf5.H:244
int readTopDataSet(std::string dsetName, std::string &buffer)
Read string data from top level HDF5 DataSet.
Definition: hdf5Reader.C:85
int numVertexVectors
number of vertex vector fields in Proteus file
Definition: proteusHdf5.H:246
int charStringLength
string length of each Proteus variable name
Definition: proteusHdf5.H:249

◆ getElementConnectivity()

std::vector< int > hdf5Reader::getElementConnectivity ( int  elementId)
inherited
Parameters
elementIdElement ID
Returns
Element connectivity table

Definition at line 195 of file hdf5Reader.C.

References hdf5Reader::elementConnectivity, hdf5Reader::vtkElementTypes, and hdf5Reader::vtkElementTypesList.

Referenced by hdf5Reader::exportToVTKMesh().

195  {
196  int nVerts;
197  int connId = -1;
198 #if VTK_MAJOR_VERSION > 8 || (VTK_MAJOR_VERSION == 8 && VTK_MINOR_VERSION >= 1)
199  if (vtkElementTypes[elementId] == VTK_LAGRANGE_QUADRILATERAL) {
200  nVerts = 16;
201  auto itr = std::find(vtkElementTypesList.begin(), vtkElementTypesList.end(),
202  VTK_LAGRANGE_QUADRILATERAL);
203  if (itr != vtkElementTypesList.end()) {
204  connId = distance(vtkElementTypesList.begin(), itr);
205  }
206  } else if (vtkElementTypes[elementId] == VTK_QUAD)
207 #else
208  if (vtkElementTypes[elementId] == VTK_QUAD)
209 #endif
210  {
211  nVerts = 4;
212  auto itr = std::find(vtkElementTypesList.begin(), vtkElementTypesList.end(),
213  VTK_QUAD);
214  if (itr != vtkElementTypesList.end()) {
215  connId = distance(vtkElementTypesList.begin(), itr);
216  }
217  } else if (vtkElementTypes[elementId] == VTK_TRIANGLE) {
218  nVerts = 3;
219  auto itr = std::find(vtkElementTypesList.begin(), vtkElementTypesList.end(),
220  VTK_TRIANGLE);
221  if (itr != vtkElementTypesList.end()) {
222  connId = distance(vtkElementTypesList.begin(), itr);
223  }
224  } else {
225  std::cerr << "Element type not supported: " << vtkElementTypes[elementId]
226  << std::endl;
227  exit(-1);
228  }
229  if (connId != -1) {
230  std::vector<int> subConnectivity(
231  elementConnectivity[connId][elementId].begin(),
232  elementConnectivity[connId][elementId].begin() + nVerts);
233  return subConnectivity;
234  } else {
235  std::cerr << "Element type " << vtkElementTypes[elementId]
236  << " not found in connectivity table" << std::endl;
237  exit(-1);
238  }
239 }
std::vector< int > vtkElementTypesList
vector storing list of each unique VTK element type
Definition: hdf5Reader.H:389
std::vector< std::map< int, std::vector< int > > > elementConnectivity
element connectivity table
Definition: hdf5Reader.H:383
std::vector< int > vtkElementTypes
vector of VTK element type for each element
Definition: hdf5Reader.H:388

◆ getFileName()

std::string hdf5Reader::getFileName ( )
inherited
Returns
HDF5 filename

Definition at line 39 of file hdf5Reader.C.

References hdf5Reader::h5FileName.

39 { return h5FileName; }
std::string h5FileName
HDF5 filename.
Definition: hdf5Reader.H:372

◆ getMeshBase()

meshBase * hdf5Reader::getMeshBase ( )
inherited
Returns
meshBase object

Definition at line 321 of file hdf5Reader.C.

References hdf5Reader::myMeshBase.

Referenced by proteusHdf5().

321 { return myMeshBase; }
meshBase * myMeshBase
meshBase object
Definition: hdf5Reader.H:378

◆ getNumBlocks()

int proteusHdf5::getNumBlocks ( ) const
Returns
Number of blocks

Definition at line 1045 of file proteusHdf5.C.

References numBlocks.

1045 { return (numBlocks); }
int numBlocks
number of Proteus blocks
Definition: proteusHdf5.H:244

◆ getNumDimensions()

int proteusHdf5::getNumDimensions ( ) const
Returns
Number of dimensions

Definition at line 1047 of file proteusHdf5.C.

References numDimensions.

Referenced by proteusHdf5().

1047 { return (numDimensions); }
int numDimensions
number of spatial dimensions
Definition: proteusHdf5.H:245

◆ getNumElementVectors()

int proteusHdf5::getNumElementVectors ( ) const
Returns
Number of element vectors

Definition at line 1051 of file proteusHdf5.C.

References numElementVectors.

1051 { return (numElementVectors); }
int numElementVectors
number of element vector fields in Proteus file
Definition: proteusHdf5.H:247

◆ getNumVertexVectors()

int proteusHdf5::getNumVertexVectors ( ) const
Returns
Number of vertex vectors

Definition at line 1049 of file proteusHdf5.C.

References numVertexVectors.

1049 { return (numVertexVectors); }
int numVertexVectors
number of vertex vector fields in Proteus file
Definition: proteusHdf5.H:246

◆ getVectorNames()

void proteusHdf5::getVectorNames ( std::string  name,
std::vector< std::string > &  stringVector,
int  numVectors 
)
Parameters
nameString description of Proteus vector field names to read
stringVectorVector of Proteus vector field names
numVectorsNumber of vector fields

Definition at line 267 of file proteusHdf5.C.

References charStringLength, and hdf5Reader::readTopDataSet().

Referenced by proteusHdf5().

269  {
270  std::cout << "Reading in PROTEUS " << name << "..." << std::endl;
271 
272  // Set vector buffer information
273  H5std_string dsetName(name);
274  std::string buffer;
275 
276  // Read dataset
277  readTopDataSet(dsetName, buffer);
278 
279  // Set element vector names
280  int bufSize = numVectors * charStringLength;
281  std::string tmp = "";
282  std::string buf = "";
283  int iCharCnt = 0;
284  int iChar = 0;
285  while (iChar <= bufSize) {
286  buf = "";
287  buf += buffer[iChar];
288  if (!(iCharCnt < charStringLength)) {
289  stringVector.push_back(std::string(tmp));
290  tmp = "";
291  iCharCnt = 0;
292  }
293  if (buf.compare(" ") != 0) {
294  tmp += buffer[iChar];
295  }
296  iChar++;
297  iCharCnt++;
298  }
299 
300  // Output info to screen
301  std::cout << "Found " << name << ":" << std::endl;
302  for (auto nameItr = stringVector.begin(); nameItr != stringVector.end();
303  nameItr++) {
304  std::cout << " " << *nameItr << std::endl;
305  }
306 }
int readTopDataSet(std::string dsetName, std::string &buffer)
Read string data from top level HDF5 DataSet.
Definition: hdf5Reader.C:85
int charStringLength
string length of each Proteus variable name
Definition: proteusHdf5.H:249

◆ getVTKMesh()

vtkSmartPointer< vtkDataSet > hdf5Reader::getVTKMesh ( )
inherited
Returns
Pointer to vtkDataSet

Definition at line 304 of file hdf5Reader.C.

References hdf5Reader::exportToVTKMesh(), and hdf5Reader::vtkMesh.

304  {
305  if (vtkMesh)
306  return vtkMesh;
307  else {
308  exportToVTKMesh();
309  return vtkMesh;
310  }
311 }
void exportToVTKMesh()
Export to VTK format.
Definition: hdf5Reader.C:241
vtkSmartPointer< vtkDataSet > vtkMesh
VTK DataSet object.
Definition: hdf5Reader.H:377

◆ mergeBlocks()

void proteusHdf5::mergeBlocks ( )
private

Definition at line 351 of file proteusHdf5.C.

References proteusSuperBlock::coordinates, hdf5Reader::elementConnectivity, proteusSuperBlock::elementData, proteusSuperBlock::elements, proteusSuperBlock::elementTypes, proteusSuperBlock::elementTypesList, proteusSuperBlock::maxNumVerticesPerElement, mySuperBlock, numDimensions, proteusSuperBlock::numElements, numElementVectors, numVertexVectors, proteusSuperBlock::numVertices, proteusSuperBlock::originalElementTypesList, proteusBlocks, proteusSuperBlock::vertexData, proteusSuperBlock::vtkConnectivity, proteusSuperBlock::xCrd, proteusSuperBlock::yCrd, and proteusSuperBlock::zCrd.

Referenced by proteusHdf5().

351  {
352  std::cout << "Merging PROTEUS block data..." << std::endl;
353 
354  // Determine total number of vertices based on maximum global ID
355  int numGlobalVertices = 0;
356 
357  // Determine total number of elements
358  int numGlobalElements = 0;
359  int maxNumVerticesPerElement = 0;
360 
361  // Determine number of element types
362  std::vector<int> elementTypesList;
363  std::vector<int> originalElementTypesList;
364 
365  for (auto blockItr = proteusBlocks.begin(); blockItr != proteusBlocks.end();
366  blockItr++) {
367  // Get number of global vertex Ids
368  int localMaxVertexId = *std::max_element((*blockItr).loc2Glob.begin(),
369  (*blockItr).loc2Glob.end());
370  if (localMaxVertexId > numGlobalVertices) {
371  numGlobalVertices = localMaxVertexId;
372  }
373 
374  // Accumulate number of global element Ids
375  numGlobalElements += (*blockItr).numElements;
376 
377  // Get maximum number of vertices per element
378  if ((*blockItr).numVerticesPerElement > maxNumVerticesPerElement) {
379  maxNumVerticesPerElement = (*blockItr).numVerticesPerElement;
380  }
381 
382  if (std::find(originalElementTypesList.begin(),
383  originalElementTypesList.end(),
384  (*blockItr).originalVtkElementType) ==
385  originalElementTypesList.end()) {
386  originalElementTypesList.push_back((*blockItr).originalVtkElementType);
387  elementTypesList.push_back((*blockItr).vtkElementType);
388  }
389  }
390  numGlobalVertices++; // global vertices are 0-indexed, so total number is
391  // maxID+1
392 
393  // Set superblock info
394  mySuperBlock.numVertices = numGlobalVertices;
395  mySuperBlock.numElements = numGlobalElements;
396  mySuperBlock.maxNumVerticesPerElement = maxNumVerticesPerElement;
397 
398  std::cout << " Found " << numGlobalVertices << " global vertices"
399  << std::endl;
400  std::cout << " Found " << numGlobalElements << " global elements"
401  << std::endl;
402 
403  // Vector to store global vertices
404  std::vector<std::vector<double>> coordinates(
405  numGlobalVertices, std::vector<double>(numDimensions, 0.0));
406 
407  // Vector to store global elements
408  std::vector<std::vector<int>> elements(
409  numGlobalElements, std::vector<int>(maxNumVerticesPerElement, 0));
410 
411  // Map of element to VTK element type
412  std::vector<int> elementTypes(numGlobalElements);
413 
414  // Vector to store vertex data
415  std::vector<std::vector<double>> vertexData(
416  numVertexVectors, std::vector<double>(numGlobalVertices, 0.0));
417 
418  // Vector to store element data
419  std::vector<std::vector<double>> elementData(
420  numElementVectors, std::vector<double>(numGlobalElements, 0.0));
421 
422  // Instantiate variables for loops below
423  int locVertId, globElemId, locElemId;
424  int lId, gId, vertNo, dimNo, origId;
425  globElemId = 0;
426 
427  // Loop over blocks
428  std::cout << " Accumulating mesh information from each block..."
429  << std::endl;
430  for (auto blockItr = proteusBlocks.begin(); blockItr != proteusBlocks.end();
431  blockItr++) {
432  std::cout << " " << (*blockItr).blockName << ":" << std::endl;
433  locVertId = 0;
434  locElemId = 0;
435  // Loop over element coordinates
436  for (auto elemItr = (*blockItr).vertices.begin();
437  elemItr != (*blockItr).vertices.end(); elemItr++) {
438  // Loop over dimension coordinates
439  origId = locVertId;
440  for (auto dimItr = (*elemItr).begin(); dimItr != (*elemItr).end();
441  dimItr++) {
442  locVertId = origId;
443  dimNo = distance((*elemItr).begin(), dimItr);
444  // Loop over vertex coordinates
445  for (auto vertItr = (*dimItr).begin(); vertItr != (*dimItr).end();
446  vertItr++) {
447  vertNo = distance((*dimItr).begin(), vertItr);
448  elements[globElemId][vertNo] = (*blockItr).loc2Glob[locVertId];
449  coordinates[(*blockItr).loc2Glob[locVertId]][dimNo] = *vertItr;
450  locVertId++;
451  }
452  elementTypes[globElemId] = (*blockItr).vtkElementType;
453  }
454 
455  // Accumulate element field data
456  for (int iElemData = 0; iElemData < numElementVectors; iElemData++) {
457  elementData[iElemData][globElemId] =
458  (*blockItr).elementData[iElemData][locElemId];
459  }
460 
461  locElemId++;
462  globElemId++;
463  }
464 
465  // Loop over vertices
466  for (auto vertItr = (*blockItr).loc2Glob.begin();
467  vertItr != (*blockItr).loc2Glob.end(); vertItr++) {
468  lId = distance((*blockItr).loc2Glob.begin(), vertItr);
469  gId = *vertItr;
470  for (int iVertData = 0; iVertData < numVertexVectors; iVertData++) {
471  vertexData[iVertData][gId] = (*blockItr).vertexData[iVertData][lId];
472  }
473  }
474 
475  std::cout << " " << (*blockItr).vertices.size() << " vertices"
476  << std::endl;
477  std::cout << " " << (*blockItr).numElements << " elements"
478  << std::endl;
479  std::cout << " " << numVertexVectors << " vertex fields"
480  << std::endl;
481  std::cout << " " << numElementVectors << " element fields"
482  << std::endl;
483  }
484 
485  // Re-structure data into separate x,y,z coordinates (todo: could be done in
486  // loops above)
487  std::cout << " Separating coordinate arrays into individual dimension "
488  "arrays..."
489  << std::endl;
490  for (auto vertItr = coordinates.begin(); vertItr != coordinates.end();
491  vertItr++) {
492  mySuperBlock.xCrd.push_back((*vertItr)[0]);
493  if (numDimensions > 1) {
494  mySuperBlock.yCrd.push_back((*vertItr)[1]);
495  }
496  if (numDimensions > 2) {
497  mySuperBlock.zCrd.push_back((*vertItr)[2]);
498  }
499  }
500 
501  // Convert element connectivity information into a format that is easier to
502  // parse for hdf5Reader class. (todo: chould be done in loops above) Main
503  // structure shown below: elementConnectivity[element type no][global
504  // element id][vertex no] = [global vertex id] elementTypesList[element type
505  // no] = vtk element type integer
506  std::cout << " Sorting element connectivities by element type..."
507  << std::endl;
508  std::vector<std::map<int, std::vector<int>>> elementConnectivity;
509 
510  // Loop over original element type list (contains all element types read in
511  // across all blocks)
512  auto origTypeItr = originalElementTypesList.begin();
513  while (origTypeItr != originalElementTypesList.end()) {
514  // temporary map for storing connectivity for each element type
515  std::map<int, std::vector<int>> tmpElementConnectivity;
516 
517  int vertNum, elemId;
518  // Loop over elements and their types
519  auto elemItr = elements.begin();
520  auto elemTypeItr = elementTypes.begin();
521  while (elemItr != elements.end() || elemTypeItr != elementTypes.end()) {
522  // Get global element id
523  elemId = distance(elements.begin(), elemItr);
524  // Loop over all vertices of each element
525  for (auto vertItr = (*elemItr).begin(); vertItr != (*elemItr).end();
526  vertItr++) {
527  // Get vertex number
528  vertNum = distance((*elemItr).begin(), vertItr);
529  // Check element type against each supported type
530  if (*elemTypeItr == VTK_QUAD) {
531 // Remove extra vertices from high-order Lagrange quads if output
532 // to lower-order quad is desired
533 #if VTK_MAJOR_VERSION > 8 || (VTK_MAJOR_VERSION == 8 && VTK_MINOR_VERSION >= 1)
534  if (*origTypeItr == VTK_LAGRANGE_QUADRILATERAL) {
535  if (vertNum == 0 || vertNum == 3 || vertNum == 12 ||
536  vertNum == 15) {
537  tmpElementConnectivity[elemId].push_back(*vertItr);
538  }
539  } else {
540  tmpElementConnectivity[elemId].push_back(*vertItr);
541  }
542 #else
543  if (*origTypeItr == -1) // lower VTK versions don't have
544  // LAGRANGE_QUADRILATERAL types
545  {
546  if (vertNum == 0 || vertNum == 3 || vertNum == 12 ||
547  vertNum == 15) {
548  tmpElementConnectivity[elemId].push_back(*vertItr);
549  }
550  } else {
551  tmpElementConnectivity[elemId].push_back(*vertItr);
552  }
553 #endif
554  }
555 #if VTK_MAJOR_VERSION > 8 || (VTK_MAJOR_VERSION == 8 && VTK_MINOR_VERSION >= 1)
556  else if (*elemTypeItr == VTK_LAGRANGE_QUADRILATERAL) {
557  tmpElementConnectivity[elemId].push_back(*vertItr);
558  }
559 #endif
560  else if (*elemTypeItr == VTK_TRIANGLE) {
561  tmpElementConnectivity[elemId].push_back(*vertItr);
562  }
563  }
564 
565  // Increment element and element type iterators
566  if (elemItr != elements.end()) {
567  ++elemItr;
568  }
569  if (elemTypeItr != elementTypes.end()) {
570  ++elemTypeItr;
571  }
572  }
573 
574  // Push back each inidividual element type connectivity
575  // into vector containing all types
576  elementConnectivity.push_back(tmpElementConnectivity);
577 
578  if (origTypeItr != originalElementTypesList.end()) {
579  ++origTypeItr;
580  }
581  }
582 
583  // Re-order vertices per VTK convention
584  std::cout << " Re-ordering element vertices per VTK convention..."
585  << std::endl;
586  int vertCount;
587  int swapVert1 = -1, swapVert2 = -1;
588  int connId;
589  // If VTK_QUAD, re-order vertices
590  auto testIt =
591  std::find(elementTypesList.begin(), elementTypesList.end(), VTK_QUAD);
592  if (testIt != elementTypesList.end()) {
593  std::cout << " Re-ordering VTK_QUAD element vertices..."
594  << std::endl;
595  connId = distance(elementTypesList.begin(), testIt);
596  for (auto elemIt = elementConnectivity[connId].begin();
597  elemIt != elementConnectivity[connId].end(); elemIt++) {
598  vertCount = 0;
599  for (int i = 0; i < (elemIt->second).size(); i++) {
600  if (vertCount == 2) {
601  swapVert1 = (elemIt->second)[i];
602  }
603  if (vertCount == 3) {
604  swapVert2 = (elemIt->second)[i];
605  (elemIt->second)[i] = swapVert1;
606  (elemIt->second)[i - 1] = swapVert2;
607  vertCount = -1;
608  }
609  vertCount++;
610  }
611  }
612  }
613  // If VTK_TRI, re-order vertices
614  /*
615  testIt = std::find(elementTypesList.begin(), elementTypesList.end(),
616  VTK_TRIANGLE); if (testIt != elementTypesList.end())
617  {
618  std::cout << " Re-ordering VTK_TRIANGLE element vertices..." <<
619  std::endl; connId = distance(elementTypesList.begin(), testIt); for (auto
620  elemIt = elementConnectivity[connId].begin(); elemIt !=
621  elementConnectivity[connId].end(); elemIt++)
622  {
623  vertCount = 0;
624  for (int i = 0; i < (elemIt->second).size(); i++)
625  {
626  if (vertCount == 1)
627  {
628  swapVert1 = (elemIt->second)[i];
629  }
630  if (vertCount == 2)
631  {
632  swapVert2 = (elemIt->second)[i];
633  (elemIt->second)[i] = swapVert1;
634  (elemIt->second)[i-1] = swapVert2;
635  vertCount = -1;
636  }
637  vertCount++;
638  }
639  }
640  }
641  */
642 
643  // Set data in superBlock struct
644  std::cout << " Setting mesh data in PROTEUS superblock..." << std::endl;
645  mySuperBlock.coordinates = coordinates;
646  mySuperBlock.elements = elements;
647  mySuperBlock.vertexData = vertexData;
648  mySuperBlock.elementData = elementData;
650  mySuperBlock.elementTypesList = elementTypesList;
651  mySuperBlock.originalElementTypesList = originalElementTypesList;
652  mySuperBlock.elementTypes = elementTypes;
653 
654  std::cout << "Finished merging PROTEUS block data..." << std::endl;
655 }
int numElementVectors
number of element vector fields in Proteus file
Definition: proteusHdf5.H:247
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
std::vector< double > yCrd
Definition: proteusHdf5.H:76
std::vector< int > elementTypes
Definition: proteusHdf5.H:83
std::vector< int > originalElementTypesList
list of original element types, ordered by type according to elements vector
Definition: proteusHdf5.H:90
std::vector< std::map< int, std::vector< int > > > elementConnectivity
element connectivity table
Definition: hdf5Reader.H:383
std::vector< double > zCrd
coordinates separated out into x,y,z coordinates
Definition: proteusHdf5.H:76
proteusSuperBlock mySuperBlock
global block structure for entire Proteus mesh
Definition: proteusHdf5.H:258
int numDimensions
number of spatial dimensions
Definition: proteusHdf5.H:245
int numElements
global number of elements
Definition: proteusHdf5.H:68
std::vector< proteusBlock > proteusBlocks
vector of all individual Proteus-style blocks
Definition: proteusHdf5.H:261
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
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< double > xCrd
Definition: proteusHdf5.H:76
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< 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

◆ openFile()

int hdf5Reader::openFile ( )
inherited
Returns
Integer flag: 0 (failure) or 1 (success)

Definition at line 41 of file hdf5Reader.C.

References hdf5Reader::file, and hdf5Reader::h5FileName.

Referenced by proteusHdf5().

41  {
42  std::cout << "Opening HDF5 file " << h5FileName << std::endl;
43 
44  try {
45  // Open HDF5 file
46  H5std_string h5FileName_(h5FileName);
47  file.openFile(h5FileName_, H5F_ACC_RDWR);
48  } catch (FileIException error) {
49 #if H5_VERS_MAJOR >= 1 && H5_VERS_MINOR >= 8 && H5_VERS_RELEASE >= 20
50  error.printErrorStack();
51 #elif H5_VERS_MAJOR >= 1 && H5_VERS_MINOR >= 10 && H5_VERS_RELEASE >= 2
52  error.printErrorStack();
53 #elif H5_VERS_MAJOR >= 1 && H5_VERS_MINOR >= 11 && H5_VERS_RELEASE >= 0
54  error.printErrorStack();
55 #else
56  error.printError();
57 #endif
58  return -1;
59  }
60  return 0;
61 }
H5File file
HDF5 file object.
Definition: hdf5Reader.H:376
std::string h5FileName
HDF5 filename.
Definition: hdf5Reader.H:372

◆ readDataSet() [1/2]

int hdf5Reader::readDataSet ( DataSet  dset,
std::string &  buffer 
)
inherited
Parameters
dsetHDF5 dataset
bufferBuffer into which data is read
Returns
Integer flag: 0 (failure) or 1 (success)

Definition at line 109 of file hdf5Reader.C.

Referenced by hdf5Reader::readTopDataSet().

109  {
110  try {
111  // Get type of data
112  DataType dtype = dset.getDataType();
113  // Read data
114  dset.read(buffer, dtype);
115  // Close dataset
116  dset.close();
117  }
118  // catch failure caused by the H5File operations
119  catch (DataSetIException error) {
120 #if H5_VERS_MAJOR >= 1 && H5_VERS_MINOR >= 8 && H5_VERS_RELEASE >= 20
121  error.printErrorStack();
122 #elif H5_VERS_MAJOR >= 1 && H5_VERS_MINOR >= 10 && H5_VERS_RELEASE >= 2
123  error.printErrorStack();
124 #elif H5_VERS_MAJOR >= 1 && H5_VERS_MINOR >= 11 && H5_VERS_RELEASE >= 0
125  error.printErrorStack();
126 #else
127  error.printError();
128 #endif
129  return -1;
130  }
131  return 0;
132 }

◆ readDataSet() [2/2]

template<class T >
int hdf5Reader::readDataSet ( DataSet  dset,
std::vector< T > &  buffer 
)
inlineinherited
Parameters
dsetHDF5 dataset
bufferBuffer into which numeric data is read
Returns
Integer flag: 0 (failure) or 1 (success)

Definition at line 148 of file hdf5Reader.H.

148  {
149  try {
150  // Get type of data
151  DataType dtype = dset.getDataType();
152  // Read data
153  dset.read(&buffer[0], dtype);
154  // Close dataset
155  dset.close();
156  }
157  // catch failure caused by the H5File operations
158  catch (DataSetIException error) {
159 #if H5_VERS_MAJOR >= 1 && H5_VERS_MINOR >= 8 && H5_VERS_RELEASE >= 20
160  error.printErrorStack();
161 #elif H5_VERS_MAJOR >= 1 && H5_VERS_MINOR >= 10 && H5_VERS_RELEASE >= 2
162  error.printErrorStack();
163 #elif H5_VERS_MAJOR >= 1 && H5_VERS_MINOR >= 11 && H5_VERS_RELEASE >= 0
164  error.printErrorStack();
165 #else
166  error.printError();
167 #endif
168  return -1;
169  }
170 
171  return 0;
172  }

◆ readGroup()

int hdf5Reader::readGroup ( std::string  groupName,
Group &  group 
)
inherited
Parameters
groupNameHDF5 group object name
groupHDF5 group object
Returns
Integer flag: 0 (failure) or 1 (success)

Definition at line 135 of file hdf5Reader.C.

References hdf5Reader::file.

Referenced by getBlocks().

135  {
136  try {
137  // Open group
138  group = Group(file.openGroup(groupName));
139  } catch (FileIException error) {
140 #if H5_VERS_MAJOR >= 1 && H5_VERS_MINOR >= 8 && H5_VERS_RELEASE >= 20
141  error.printErrorStack();
142 #elif H5_VERS_MAJOR >= 1 && H5_VERS_MINOR >= 10 && H5_VERS_RELEASE >= 2
143  error.printErrorStack();
144 #elif H5_VERS_MAJOR >= 1 && H5_VERS_MINOR >= 11 && H5_VERS_RELEASE >= 0
145  error.printErrorStack();
146 #else
147  error.printError();
148 #endif
149  return -1;
150  }
151  return 0;
152 }
H5File file
HDF5 file object.
Definition: hdf5Reader.H:376

◆ readGroupDataSet()

template<class T >
int hdf5Reader::readGroupDataSet ( std::string  dsetName,
Group &  group,
std::vector< T > &  buffer 
)
inlineinherited
Parameters
dsetNameDataSet name
groupHDF5 group object
bufferBuffer into which numeric data is read
Returns
Integer flag: 0 (failure) or 1 (success)

Definition at line 190 of file hdf5Reader.H.

Referenced by getBlockElementData(), getBlockGlobalID(), getBlockInfo(), getBlockVertexData(), and getBlockXYZ().

191  {
192  try {
193  // Open dataset
194  DataSet dset = group.openDataSet(dsetName);
195  // Get type of data
196  DataType dtype = dset.getDataType();
197  // Read data
198  dset.read(&buffer[0], dtype);
199  // Close dataset
200  dset.close();
201  }
202  // catch failure caused by the H5File operations
203  catch (GroupIException error) {
204 #if H5_VERS_MAJOR >= 1 && H5_VERS_MINOR >= 8 && H5_VERS_RELEASE >= 20
205  error.printErrorStack();
206 #elif H5_VERS_MAJOR >= 1 && H5_VERS_MINOR >= 10 && H5_VERS_RELEASE >= 2
207  error.printErrorStack();
208 #elif H5_VERS_MAJOR >= 1 && H5_VERS_MINOR >= 11 && H5_VERS_RELEASE >= 0
209  error.printErrorStack();
210 #else
211  error.printError();
212 #endif
213  return -1;
214  } catch (DataSetIException error) {
215 #if H5_VERS_MAJOR >= 1 && H5_VERS_MINOR >= 8 && H5_VERS_RELEASE >= 20
216  error.printErrorStack();
217 #elif H5_VERS_MAJOR >= 1 && H5_VERS_MINOR >= 10 && H5_VERS_RELEASE >= 2
218  error.printErrorStack();
219 #elif H5_VERS_MAJOR >= 1 && H5_VERS_MINOR >= 11 && H5_VERS_RELEASE >= 0
220  error.printErrorStack();
221 #else
222  error.printError();
223 #endif
224  return -1;
225  }
226 
227  return 0;
228  }

◆ readTopDataSet() [1/2]

int hdf5Reader::readTopDataSet ( std::string  dsetName,
std::string &  buffer 
)
inherited
Parameters
dsetNameDataset name
bufferBuffer into which string data is read
Returns
Integer flag: 0 (failure) or 1 (success)

Definition at line 85 of file hdf5Reader.C.

References hdf5Reader::file, and hdf5Reader::readDataSet().

Referenced by getControlInformation(), and getVectorNames().

85  {
86  try {
87  // Open dataset
88  DataSet dset = file.openDataSet(dsetName);
89  // Read dataset
90  readDataSet(dset, buffer);
91  }
92  // catch failure caused by the H5File operations
93  catch (FileIException error) {
94 #if H5_VERS_MAJOR >= 1 && H5_VERS_MINOR >= 8 && H5_VERS_RELEASE >= 20
95  error.printErrorStack();
96 #elif H5_VERS_MAJOR >= 1 && H5_VERS_MINOR >= 10 && H5_VERS_RELEASE >= 2
97  error.printErrorStack();
98 #elif H5_VERS_MAJOR >= 1 && H5_VERS_MINOR >= 11 && H5_VERS_RELEASE >= 0
99  error.printErrorStack();
100 #else
101  error.printError();
102 #endif
103  return -1;
104  }
105  return 0;
106 }
H5File file
HDF5 file object.
Definition: hdf5Reader.H:376
int readDataSet(DataSet dset, std::string &buffer)
Read string data from generic HDF5 DataSet.
Definition: hdf5Reader.C:109

◆ readTopDataSet() [2/2]

template<class T >
int hdf5Reader::readTopDataSet ( std::string  dsetName,
std::vector< T > &  buffer 
)
inlineinherited
Parameters
dsetNameDataset name
bufferBuffer into which numeric data is read
Returns
Integer flag: 0 (failure) or 1 (success)

Definition at line 115 of file hdf5Reader.H.

115  {
116  try {
117  std::cout << "H5_VERS_MAJOR = " << H5_VERS_MAJOR << std::endl;
118  std::cout << "H5_VERS_MINOR = " << H5_VERS_MINOR << std::endl;
119  std::cout << "H5_VERS_RELEASE = " << H5_VERS_RELEASE << std::endl;
120  std::cout << "H5_VERS_SUBRELEASE = " << H5_VERS_SUBRELEASE << std::endl;
121  // Open dataset
122  DataSet dset = file.openDataSet(dsetName);
123  readDataSet(dset, buffer);
124  }
125  // catch failure caused by the H5File operations
126  catch (FileIException error) {
127 #if H5_VERS_MAJOR >= 1 && H5_VERS_MINOR >= 8 && H5_VERS_RELEASE >= 20
128  error.printErrorStack();
129 #elif H5_VERS_MAJOR >= 1 && H5_VERS_MINOR >= 10 && H5_VERS_RELEASE >= 2
130  error.printErrorStack();
131 #elif H5_VERS_MAJOR >= 1 && H5_VERS_MINOR >= 11 && H5_VERS_RELEASE >= 0
132  error.printErrorStack();
133 #else
134  error.printError();
135 #endif
136  return -1;
137  }
138 
139  return 0;
140  }
H5File file
HDF5 file object.
Definition: hdf5Reader.H:376
int readDataSet(DataSet dset, std::string &buffer)
Read string data from generic HDF5 DataSet.
Definition: hdf5Reader.C:109

◆ setConnectivities()

void hdf5Reader::setConnectivities ( std::vector< std::map< int, std::vector< int >>>  elementConnectivity)
inherited
Parameters
elementConnectivityElement connectivity table

Definition at line 190 of file hdf5Reader.C.

References hdf5Reader::elementConnectivity.

Referenced by proteusHdf5().

191  {
192  elementConnectivity = _elementConnectivity;
193 }

◆ setCoordinates() [1/3]

void hdf5Reader::setCoordinates ( std::vector< double >  xCrd)
inherited
Parameters
xCrdvector of x-coordinates

Definition at line 166 of file hdf5Reader.C.

References hdf5Reader::xCrd.

Referenced by proteusHdf5().

166 { xCrd = _xCrd; }

◆ setCoordinates() [2/3]

void hdf5Reader::setCoordinates ( std::vector< double >  xCrd,
std::vector< double >  yCrd 
)
inherited
Parameters
xCrdvector of x-coordinates
yCrdvector of y-coordinates

Definition at line 168 of file hdf5Reader.C.

References hdf5Reader::xCrd, and hdf5Reader::yCrd.

169  {
170  xCrd = _xCrd;
171  yCrd = _yCrd;
172 }

◆ setCoordinates() [3/3]

void hdf5Reader::setCoordinates ( std::vector< double >  xCrd,
std::vector< double >  yCrd,
std::vector< double >  zCrd 
)
inherited
Parameters
xCrdvector of x-coordinates
yCrdvector of y-coordinates
zCrdvector of z-coordinates

Definition at line 174 of file hdf5Reader.C.

References hdf5Reader::xCrd, hdf5Reader::yCrd, and hdf5Reader::zCrd.

176  {
177  xCrd = _xCrd;
178  yCrd = _yCrd;
179  zCrd = _zCrd;
180 }

◆ setElementTypes()

void hdf5Reader::setElementTypes ( std::vector< int >  vtkElementTypes)
inherited
Parameters
vtkElementTypesVector of VTK element types for each element

Definition at line 182 of file hdf5Reader.C.

References hdf5Reader::vtkElementTypes.

Referenced by proteusHdf5().

182  {
183  vtkElementTypes = _vtkElementTypes;
184 }

◆ setElementTypesList()

void hdf5Reader::setElementTypesList ( std::vector< int >  vtkElementTypesList)
inherited
Parameters
vtkElementTypesListVector of each unique VTK element type in mesh

Definition at line 186 of file hdf5Reader.C.

References hdf5Reader::vtkElementTypesList.

Referenced by proteusHdf5().

186  {
187  vtkElementTypesList = _vtkElementTypesList;
188 }

◆ setFields()

void hdf5Reader::setFields ( meshBase myMeshBase,
std::vector< std::string >  dataNames,
std::vector< std::vector< double >>  data,
int  pointOrCell 
)
inherited
Parameters
myMeshBasemeshBase object
dataNamesVector of names of data fields
dataVector storing data
pointOrCellBoolean indicating point or cell data

Definition at line 323 of file hdf5Reader.C.

References data, meshBase::setCellDataArray(), and meshBase::setPointDataArray().

Referenced by proteusHdf5().

326  {
327  if (pointOrCell == 0) {
328  std::cout << "Setting hdf5Reader object point fields..." << std::endl;
329  } else if (pointOrCell == 1) {
330  std::cout << "Setting hdf5Reader object cell fields..." << std::endl;
331  }
332 
333  // Set element data
334  auto dataNameItr = dataNames.begin();
335  auto dataItr = data.begin();
336  while (dataNameItr != dataNames.end() || dataItr != data.end()) {
337  std::cout << " " << *dataNameItr << std::endl;
338 
339  if (pointOrCell == 0) {
340  myMeshBase->setPointDataArray((*dataNameItr).c_str(), *dataItr);
341  } else if (pointOrCell == 1) {
342  myMeshBase->setCellDataArray((*dataNameItr).c_str(), *dataItr);
343  }
344 
345  if (dataNameItr != dataNames.end()) {
346  dataNameItr++;
347  }
348  if (dataItr != data.end()) {
349  dataItr++;
350  }
351  }
352 }
virtual void setPointDataArray(const std::string &name, const std::vector< std::vector< double >> &data)
register data to dataSet&#39;s point data
Definition: meshBase.H:319
virtual void setCellDataArray(const std::string &name, const std::vector< std::vector< double >> &data)
register data to dataSet&#39;s cell data
Definition: meshBase.H:334

◆ setNumberOfDimensions()

void hdf5Reader::setNumberOfDimensions ( int  numDimensions)
inherited
Parameters
numDimensionsNumber of spatial dimensions

Definition at line 162 of file hdf5Reader.C.

References hdf5Reader::numDimensions.

Referenced by proteusHdf5().

162  {
163  numDimensions = _numDimensions;
164 }
int numDimensions
number of spatial dimensions
Definition: hdf5Reader.H:386

◆ setNumberOfElements()

void hdf5Reader::setNumberOfElements ( int  numElements)
inherited
Parameters
numElementsNumber of elements

Definition at line 158 of file hdf5Reader.C.

References hdf5Reader::numElements.

Referenced by proteusHdf5().

158  {
159  numElements = _numElements;
160 }
int numElements
number of elements in mesh
Definition: hdf5Reader.H:385

◆ setNumberOfVertices()

void hdf5Reader::setNumberOfVertices ( int  numVertices)
inherited
Parameters
numVerticesNumber of vertices

Definition at line 154 of file hdf5Reader.C.

References hdf5Reader::numVertices.

Referenced by proteusHdf5().

154  {
155  numVertices = _numVertices;
156 }
int numVertices
number of vertices in mesh
Definition: hdf5Reader.H:384

◆ unflattenBuffer2d()

template<class T >
void hdf5Reader::unflattenBuffer2d ( std::vector< T >  flatBuffer,
std::vector< std::vector< T >> &  originalBuffer 
)
inlineinherited
Parameters
flatBufferOriginal flattened data
originalBufferOutput unflattened data

Definition at line 258 of file hdf5Reader.H.

References data.

Referenced by getBlockElementData(), and getBlockVertexData().

259  {
260  // Populate original buffer
261  auto flatItr = flatBuffer.begin();
262  for (auto bufItr1 = originalBuffer.begin(); bufItr1 != originalBuffer.end();
263  bufItr1++) {
264  for (auto bufItr2 = (*bufItr1).begin(); bufItr2 != (*bufItr1).end();
265  bufItr2++) {
266  *bufItr2 = *flatItr;
267  flatItr++;
268  }
269  }
270  }

◆ writeVTK()

void hdf5Reader::writeVTK ( )
inherited

Definition at line 354 of file hdf5Reader.C.

References hdf5Reader::exportToMeshBase(), hdf5Reader::myMeshBase, and meshBase::write().

354  {
355  std::cout << "Writing hdf5Reader object to VTK file" << std::endl;
356  if (!myMeshBase) {
358  }
359  myMeshBase->write();
360 }
virtual void write() const
write the mesh to file named after the private var &#39;filename&#39;.
Definition: meshBase.H:598
void exportToMeshBase()
Export HDF5 mesh to meshBase object.
Definition: hdf5Reader.C:313
meshBase * myMeshBase
meshBase object
Definition: hdf5Reader.H:378

Member Data Documentation

◆ baseH5FileName

std::string hdf5Reader::baseH5FileName
protectedinherited

Definition at line 373 of file hdf5Reader.H.

◆ charStringLength

int proteusHdf5::charStringLength
private

Definition at line 249 of file proteusHdf5.H.

Referenced by getCharStringLength(), getControlInformation(), and getVectorNames().

◆ elementConnectivity

std::vector<std::map<int, std::vector<int> > > hdf5Reader::elementConnectivity
protectedinherited

◆ elementVectorNames

std::vector<std::string> proteusHdf5::elementVectorNames
private

Definition at line 255 of file proteusHdf5.H.

Referenced by getBlockElementData(), and proteusHdf5().

◆ file

H5File hdf5Reader::file
protectedinherited

◆ h5FileName

std::string hdf5Reader::h5FileName
protectedinherited

◆ lowOrder

bool proteusHdf5::lowOrder
private

Definition at line 240 of file proteusHdf5.H.

Referenced by getBlockInfo(), and proteusHdf5().

◆ myMeshBase

meshBase* hdf5Reader::myMeshBase
protectedinherited

◆ mySuperBlock

proteusSuperBlock proteusHdf5::mySuperBlock
private

Definition at line 258 of file proteusHdf5.H.

Referenced by mergeBlocks(), and proteusHdf5().

◆ numBlocks

int proteusHdf5::numBlocks
private

Definition at line 244 of file proteusHdf5.H.

Referenced by getBlocks(), getControlInformation(), and getNumBlocks().

◆ numDimensions

int proteusHdf5::numDimensions
private

◆ numElements

int hdf5Reader::numElements
protectedinherited

Definition at line 385 of file hdf5Reader.H.

Referenced by hdf5Reader::exportToVTKMesh(), and hdf5Reader::setNumberOfElements().

◆ numElementVectors

int proteusHdf5::numElementVectors
private

◆ numVertexVectors

int proteusHdf5::numVertexVectors
private

◆ numVertices

int hdf5Reader::numVertices
protectedinherited

Definition at line 384 of file hdf5Reader.H.

Referenced by hdf5Reader::exportToVTKMesh(), and hdf5Reader::setNumberOfVertices().

◆ outputFname

std::string hdf5Reader::outputFname
protectedinherited

Definition at line 374 of file hdf5Reader.H.

Referenced by hdf5Reader::exportToMeshBase().

◆ proteusBlocks

std::vector<proteusBlock> proteusHdf5::proteusBlocks
private

Definition at line 261 of file proteusHdf5.H.

Referenced by getBlocks(), and mergeBlocks().

◆ verb

int hdf5Reader::verb
protectedinherited

Definition at line 375 of file hdf5Reader.H.

◆ vertexVectorNames

std::vector<std::string> proteusHdf5::vertexVectorNames
private

Definition at line 253 of file proteusHdf5.H.

Referenced by getBlockVertexData(), and proteusHdf5().

◆ vtkElementTypes

std::vector<int> hdf5Reader::vtkElementTypes
protectedinherited

◆ vtkElementTypesList

std::vector<int> hdf5Reader::vtkElementTypesList
protectedinherited

◆ vtkMesh

vtkSmartPointer<vtkDataSet> hdf5Reader::vtkMesh
protectedinherited

Definition at line 377 of file hdf5Reader.H.

Referenced by hdf5Reader::getVTKMesh().

◆ xCrd

std::vector<double> hdf5Reader::xCrd
protectedinherited

Definition at line 381 of file hdf5Reader.H.

Referenced by hdf5Reader::exportToVTKMesh(), and hdf5Reader::setCoordinates().

◆ yCrd

std::vector<double> hdf5Reader::yCrd
protectedinherited

Definition at line 381 of file hdf5Reader.H.

Referenced by hdf5Reader::exportToVTKMesh(), and hdf5Reader::setCoordinates().

◆ zCrd

std::vector<double> hdf5Reader::zCrd
protectedinherited

Definition at line 381 of file hdf5Reader.H.

Referenced by hdf5Reader::exportToVTKMesh(), and hdf5Reader::setCoordinates().


The documentation for this class was generated from the following files: