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.
NEM::MSH Namespace Reference

Namespaces

 EXOMesh
 

Classes

class  dataSetRegionBoundaryFilter
 Like vtkDataSetRegionSurfaceFilter, this filter extracts boundaries of materials, including both interfaces and external boundaries. More...
 
class  exoGeoMesh
 A concrete implementation of geoMeshBase representing an Exodus II mesh. More...
 
class  foamGeoMesh
 A concrete implementation of geoMeshBase representing a mesh in a fvMesh. More...
 
class  geoMeshBase
 abstract class to specify geometry and mesh data More...
 
class  gmshGeoMesh
 A concrete implementation of geoMeshBase representing a Gmsh mesh. More...
 
class  GmshInterface
 management class for Gmsh interface More...
 
class  inpGeoMesh
 Class representing meshes in CalculiX input deck (similar to ABAQUS) More...
 
class  mergeCells
 
class  OmegaHInterface
 management class for Omega_h::Library More...
 
class  oshGeoMesh
 A concrete implementation of geoMeshBase representing a Omega_h::Mesh. More...
 
struct  SM_StdContainerWrapperFromIter
 
class  smeshGeoMesh
 
class  vtkGeoMesh
 A concrete implementation of geoMeshBase representing a mesh in a vtkUnstructuredGrid. More...
 

Typedefs

using nemId_t = std::int64_t
 

Functions

bool is_close (double val1, double val2, double floor, double rel_tol)
 
int diffMesh (geoMeshBase *gmb1, geoMeshBase *gmb2, double floor=1e-9, double relTol=1e-6, double numCellsTol=-1., double numPointsTol=-1.)
 Compare two geoMeshBase objects. More...
 
 vtkStandardNewMacro (exoGeoMesh)
 
 vtkStandardNewMacro (foamGeoMesh)
 
MeshType MeshTypeFromFilename (const std::string &fileName)
 
geoMeshBaseRead (const std::string &fileName)
 Read a mesh from file. More...
 
geoMeshBaseRead (const std::string &fileName, MeshType type)
 Read a mesh from file. More...
 
geoMeshBaseNew (MeshType type)
 Create a new mesh object. More...
 
geoMeshBaseNew (const std::string &fileName)
 Create a new mesh object. More...
 
 vtkStandardNewMacro (gmshGeoMesh)
 
 vtkStandardNewMacro (inpGeoMesh)
 
 vtkStandardNewMacro (oshGeoMesh)
 
VTKCellType getVTKTypeFromOmega_hFamilyDim (Omega_h_Family family, Omega_h::Int dim)
 
Omega_h_Family getOmega_hFamilyFromVTKType (VTKCellType vtkCellType)
 
Omega_h::Int getOmega_hDimFromVTKType (VTKCellType vtkCellType)
 
template<typename OT , typename VT >
void Otag2Varray (const Omega_h::TagBase *tagBase, VT *vtkArray)
 
void getVtkDataArrayFromOmega_hTag (const Omega_h::TagBase *tagBase, vtkSmartPointer< vtkDataArray > &vtkArray)
 
int oshFace2vtkFace (int oshFace, Omega_h_Family oshFamily, Omega_h::Int oshDim)
 
template<typename VT , typename OT >
void Varray2Otag (Omega_h::Mesh *oshMesh, VT *vtkArray, Omega_h::Int oshDim)
 
template<typename VT >
void Varray2Otag2 (Omega_h::Mesh *oshMesh, VT *vtkArray, Omega_h::Int oshDim)
 
void getOmega_hArrayFromVtkDataArray (Omega_h::Mesh *oshMesh, vtkSmartPointer< vtkDataArray > &vtkArray, Omega_h::Int oshDim)
 
 vtkStandardNewMacro (smeshGeoMesh)
 
 mesh_ (gen_->CreateMesh(false))
 
template<class PtrSMDSIterator , class VALUE = typename std::decay< PtrSMDSIterator>::type::element_type::value_type, class EqualVALUE = std::equal_to<VALUE>>
SM_StdContainerWrapperFromIter< typename std::decay< PtrSMDSIterator >::type, VALUE, EqualVALUE > containerWrapper (PtrSMDSIterator &&iter)
 
 vtkStandardNewMacro (vtkGeoMesh)
 
 vtkStandardNewMacro (dataSetRegionBoundaryFilter)
 
 vtkStandardNewMacro (mergeCells)
 

Typedef Documentation

◆ nemId_t

using NEM::MSH::nemId_t = typedef std::int64_t

Definition at line 56 of file geoMeshBase.H.

Function Documentation

◆ containerWrapper()

template<class PtrSMDSIterator , class VALUE = typename std::decay< PtrSMDSIterator>::type::element_type::value_type, class EqualVALUE = std::equal_to<VALUE>>
SM_StdContainerWrapperFromIter<typename std::decay<PtrSMDSIterator>::type, VALUE, EqualVALUE> NEM::MSH::containerWrapper ( PtrSMDSIterator &&  iter)

◆ diffMesh()

NEMOSYS_EXPORT int NEM::MSH::diffMesh ( geoMeshBase gmb1,
geoMeshBase gmb2,
double  floor = 1e-9,
double  relTol = 1e-6,
double  numCellsTol = -1.,
double  numPointsTol = -1. 
)
Parameters
gmb1Reference geoMeshBase
gmb2Other geoMeshBase
floorIf a floating point value is less than or equal to floor, it is treated as 0. Should be non-negative
relTolRelative tolerance for comparing floating point values
numCellsTolTolerance for number of cells, if positive. If non-positive, exact match required for number of cells and cell data compared at each cell.
numPointsTolTolerance for number of points, if positive. If non-positive, exact match required for number of points and point data compared at each cell.
Returns
0 if gmb1 and gmb2 are similar; 1 otherwise

Definition at line 51 of file diffMesh.C.

References NEM::MSH::geoMeshBase::getCell(), NEM::MSH::geoMeshBase::getCellDataArrayCopy(), NEM::MSH::geoMeshBase::getGeoEntArrayName(), NEM::MSH::geoMeshBase::getNumberOfCellDataArrays(), NEM::MSH::geoMeshBase::getNumberOfCells(), NEM::MSH::geoMeshBase::getNumberOfPointDataArrays(), NEM::MSH::geoMeshBase::getNumberOfPoints(), NEM::MSH::geoMeshBase::getPoint(), NEM::MSH::geoMeshBase::getPointDataArrayCopy(), and is_close().

52  {
53  // Check if equal number of points and cells.
54  {
55  auto points1 = gmb1->getNumberOfPoints();
56  auto points2 = gmb2->getNumberOfPoints();
57  if (numPointsTol > 0) {
58  if (!is_close(points1, points2, 0, numPointsTol)) {
59  std::cerr << "Meshes don't have similar number of points: " << points1
60  << " vs " << points2 << std::endl;
61  return 1;
62  }
63  } else {
64  if (points1 != points2) {
65  std::cerr << "Meshes don't have same number of points: " << points1
66  << " vs " << points2 << std::endl;
67  return 1;
68  }
69  }
70  }
71  {
72  auto cells1 = gmb1->getNumberOfCells();
73  auto cells2 = gmb2->getNumberOfCells();
74  if (numCellsTol > 0) {
75  if (!is_close(cells1, cells2, 0, numCellsTol)) {
76  std::cerr << "Meshes don't have similar number of cells: " << cells1
77  << " vs " << cells2 << std::endl;
78  return 1;
79  }
80  } else {
81  if (cells1 != cells2) {
82  std::cerr << "Meshes don't have same number of cells: " << cells1
83  << " vs " << cells2 << std::endl;
84  return 1;
85  }
86  }
87  }
88 
89  if (numPointsTol <= 0) {
90  // Check if point coordinates match within tolerance
91  for (int i = 0; i < gmb1->getNumberOfPoints(); ++i) {
92  std::array<double, 3> coord1{};
93  std::array<double, 3> coord2{};
94  gmb1->getPoint(i, &coord1);
95  gmb2->getPoint(i, &coord2);
96  for (int j = 0; j < 3; ++j) {
97  if (!is_close(coord1[j], coord2[j], floor, relTol)) {
98  std::cerr << "Meshes differ in point coordinates" << std::endl;
99  std::cerr << "Index " << i << " Component " << j << std::endl;
100  std::cerr << "Coord 1 " << std::setprecision(15) << coord1[j]
101  << " Coord 2 " << std::setprecision(15) << coord2[j]
102  << std::endl;
103  std::cerr << "Meshes differ in point coordinates" << std::endl;
104  return 1;
105  }
106  }
107  }
108  }
109 
110  if (numCellsTol <= 0) {
111  for (int i = 0; i < gmb1->getNumberOfCells(); ++i) {
112  vtkSmartPointer<vtkPoints> points1 = gmb1->getCell(i)->GetPoints();
113  vtkSmartPointer<vtkPoints> points2 = gmb2->getCell(i)->GetPoints();
114 
115  if (points1->GetNumberOfPoints() != points2->GetNumberOfPoints()) {
116  std::cerr << "Meshes differ in cells" << std::endl;
117  return 1;
118  }
119 
120  for (int j = 0; j < points1->GetNumberOfPoints(); ++j) {
121  std::vector<double> coord1(3);
122  std::vector<double> coord2(3);
123  points1->GetPoint(j, coord1.data());
124  points2->GetPoint(j, coord2.data());
125  for (int k = 0; k < 3; ++k) {
126  if (!is_close(coord1[k], coord2[k], floor, relTol)) {
127  std::cerr << "Meshes differ in cells" << std::endl;
128  return 1;
129  }
130  }
131  }
132  }
133  }
134 
135  int numPointArr1 = gmb1->getNumberOfPointDataArrays();
136  int numPointArr2 = gmb2->getNumberOfPointDataArrays();
137 
138  if (numPointArr1 != numPointArr2) {
139  std::cerr << "Meshes have different numbers of point data" << std::endl;
140  std::cerr << "Mesh 1 has " << numPointArr1 << std::endl;
141  std::cerr << "Mesh 2 has " << numPointArr2 << std::endl;
142  return 1;
143  }
144 
145  for (int i = 0; i < numPointArr1; ++i) {
146  auto arr1 = gmb1->getPointDataArrayCopy(i);
147  auto arr2 = gmb2->getPointDataArrayCopy(arr1->GetName());
148  if (!arr2) {
149  std::cerr << "Mesh 2 does not have a point data array of name "
150  << arr1->GetName() << std::endl;
151  return 1;
152  }
153  /*
154  if (arr1->GetDataType() != arr2->GetDataType()) {
155  std::cerr << "For point data array " << arr1->GetName()
156  << ",\narrays do not have same data type:\n"
157  << arr1->GetDataTypeAsString() << " vs "
158  << arr2->GetDataTypeAsString() << std::endl;
159  return 1;
160  }
161  */
162  if (arr1->GetNumberOfComponents() != arr2->GetNumberOfComponents()) {
163  std::cerr << "For point data array " << arr1->GetName()
164  << ",\narrays do not have same number of components:\n"
165  << arr1->GetNumberOfComponents() << " vs "
166  << arr2->GetNumberOfComponents() << std::endl;
167  return 1;
168  }
169  if (numPointsTol <= 0) {
170  auto da1 = vtkDataArray::FastDownCast(arr1);
171  auto da2 = vtkDataArray::FastDownCast(arr2);
172  if (da1 && da2) {
173  int numComponents = da1->GetNumberOfComponents();
174  for (int j = 0; j < gmb1->getNumberOfPoints(); ++j) {
175  std::vector<double> comps1(numComponents);
176  std::vector<double> comps2(numComponents);
177  da1->GetTuple(j, comps1.data());
178  da2->GetTuple(j, comps2.data());
179  for (int k = 0; k < numComponents; ++k) {
180  if (!is_close(comps1[k], comps2[k], floor, relTol)) {
181  std::cerr << "For point data array " << da1->GetName()
182  << std::endl;
183  std::cerr << "Meshes differ in point data values at point " << j
184  << " component " << k << std::endl;
185  std::cerr << comps1[k] << " " << comps2[k] << std::endl;
186  return 1;
187  }
188  }
189  }
190  }
191  }
192  }
193 
194  int numCellArr1 = gmb1->getNumberOfCellDataArrays();
195  int numCellArr2 = gmb2->getNumberOfCellDataArrays();
196 
197  if (numCellArr1 != numCellArr2) {
198  std::cerr << "Meshes have different numbers of cell data" << std::endl;
199  std::cerr << "Mesh 1 has " << numCellArr1 << std::endl;
200  std::cerr << "Mesh 2 has " << numCellArr2 << std::endl;
201  return 1;
202  }
203 
204  for (int i = 0; i < numCellArr1; ++i) {
205  auto arr1 = gmb1->getCellDataArrayCopy(i);
206  if (arr1->GetName() == gmb1->getGeoEntArrayName()) {
207  continue;
208  }
209  auto arr2 = gmb2->getCellDataArrayCopy(arr1->GetName());
210  if (!arr2) {
211  std::cerr << "Mesh 2 does not have a cell data array of name "
212  << arr1->GetName() << std::endl;
213  return 1;
214  }
215  /*
216  if (arr1->GetDataType() != arr2->GetDataType()) {
217  std::cerr << "For cell data array " << arr1->GetName()
218  << ",\narrays do not have same data type:\n"
219  << arr1->GetDataTypeAsString() << " vs "
220  << arr2->GetDataTypeAsString() << std::endl;
221  return 1;
222  }
223  */
224  if (arr1->GetNumberOfComponents() != arr2->GetNumberOfComponents()) {
225  std::cerr << "For cell data array " << arr1->GetName()
226  << ",\narrays do not have same number of components:\n"
227  << arr1->GetNumberOfComponents() << " vs "
228  << arr2->GetNumberOfComponents() << std::endl;
229  return 1;
230  }
231  if (numCellsTol <= 0) {
232  auto da1 = vtkDataArray::FastDownCast(arr1);
233  auto da2 = vtkDataArray::FastDownCast(arr2);
234  if (da1 && da2) {
235  int numComponents = da1->GetNumberOfComponents();
236  for (int j = 0; j < gmb1->getNumberOfCells(); ++j) {
237  std::vector<double> comps1(numComponents);
238  std::vector<double> comps2(numComponents);
239  da1->GetTuple(j, comps1.data());
240  da2->GetTuple(j, comps2.data());
241  for (int k = 0; k < numComponents; ++k) {
242  if (!is_close(comps1[k], comps2[k], floor, relTol)) {
243  std::cerr << "For cell data array " << da1->GetName()
244  << std::endl;
245  std::cerr << "Meshes differ in cell data values at cell " << j
246  << " component " << k << std::endl;
247  std::cerr << comps1[k] << " " << comps2[k] << std::endl;
248  return 1;
249  }
250  }
251  }
252  }
253  }
254  }
255  std::cout << "Meshes are the same" << std::endl;
256  return 0;
257 }
bool is_close(double val1, double val2, double floor, double rel_tol)
Definition: diffMesh.C:43

◆ getOmega_hArrayFromVtkDataArray()

void NEM::MSH::getOmega_hArrayFromVtkDataArray ( Omega_h::Mesh *  oshMesh,
vtkSmartPointer< vtkDataArray > &  vtkArray,
Omega_h::Int  oshDim 
)

Definition at line 544 of file oshGeoMesh.C.

Referenced by NEM::MSH::oshGeoMesh::GM2osh().

546  {
547  switch (vtkArray->GetDataType()) {
548  case VTK_FLOAT: {
549  vtkSmartPointer<vtkFloatArray> vtkArray_tmp =
550  vtkFloatArray::FastDownCast(vtkArray);
551  Varray2Otag<vtkFloatArray, Omega_h::Real>(oshMesh, vtkArray_tmp, oshDim);
552  vtkArray = vtkArray_tmp;
553  break;
554  }
555  case VTK_DOUBLE: {
556  vtkSmartPointer<vtkDoubleArray> vtkArray_tmp =
557  vtkDoubleArray::FastDownCast(vtkArray);
558  Varray2Otag<vtkDoubleArray, Omega_h::Real>(oshMesh, vtkArray_tmp, oshDim);
559  vtkArray = vtkArray_tmp;
560  break;
561  }
562  case VTK_CHAR: {
563  vtkSmartPointer<vtkCharArray> vtkArray_tmp =
564  vtkCharArray::FastDownCast(vtkArray);
565  Varray2Otag2<vtkCharArray>(oshMesh, vtkArray_tmp, oshDim);
566  vtkArray = vtkArray_tmp;
567  break;
568  }
569  case VTK_SIGNED_CHAR: {
570  vtkSmartPointer<vtkSignedCharArray> vtkArray_tmp =
571  vtkSignedCharArray::FastDownCast(vtkArray);
572  Varray2Otag2<vtkSignedCharArray>(oshMesh, vtkArray_tmp, oshDim);
573  vtkArray = vtkArray_tmp;
574  break;
575  }
576  case VTK_SHORT: {
577  vtkSmartPointer<vtkShortArray> vtkArray_tmp =
578  vtkShortArray::FastDownCast(vtkArray);
579  Varray2Otag2<vtkShortArray>(oshMesh, vtkArray_tmp, oshDim);
580  vtkArray = vtkArray_tmp;
581  break;
582  }
583  case VTK_INT: {
584  vtkSmartPointer<vtkIntArray> vtkArray_tmp =
585  vtkIntArray::FastDownCast(vtkArray);
586  Varray2Otag2<vtkIntArray>(oshMesh, vtkArray_tmp, oshDim);
587  vtkArray = vtkArray_tmp;
588  break;
589  }
590  case VTK_LONG: {
591  vtkSmartPointer<vtkLongArray> vtkArray_tmp =
592  vtkLongArray::FastDownCast(vtkArray);
593  Varray2Otag2<vtkLongArray>(oshMesh, vtkArray_tmp, oshDim);
594  vtkArray = vtkArray_tmp;
595  break;
596  }
597  case VTK_ID_TYPE: {
598  vtkSmartPointer<vtkIdTypeArray> vtkArray_tmp =
599  vtkIdTypeArray::FastDownCast(vtkArray);
600  Varray2Otag2<vtkIdTypeArray>(oshMesh, vtkArray_tmp, oshDim);
601  vtkArray = vtkArray_tmp;
602  break;
603  }
604  case VTK_LONG_LONG: {
605  vtkSmartPointer<vtkLongLongArray> vtkArray_tmp =
606  vtkLongLongArray::FastDownCast(vtkArray);
607  Varray2Otag2<vtkLongLongArray>(oshMesh, vtkArray_tmp, oshDim);
608  vtkArray = vtkArray_tmp;
609  break;
610  }
611  /*
612  case VTK_VOID:
613  case VTK_BIT:
614  case VTK_UNSIGNED_CHAR:
615  case VTK_UNSIGNED_SHORT:
616  case VTK_UNSIGNED_INT:
617  case VTK_UNSIGNED_LONG:
618  case VTK_STRING:
619  case VTK_OPAQUE:
620  case VTK_UNSIGNED_LONG_LONG:
621  case VTK___INT64:
622  case VTK_UNSIGNED___INT64:
623  case VTK_VARIANT:
624  case VTK_OBJECT:
625  case VTK_UNICODE_STRING:
626  */
627  default: {
628  std::cerr << "VTK type " << vtkArray->GetDataTypeAsString()
629  << " is not supported as an Omega_h type." << std::endl;
630  exit(1);
631  }
632  }
633 }

◆ getOmega_hDimFromVTKType()

Omega_h::Int NEM::MSH::getOmega_hDimFromVTKType ( VTKCellType  vtkCellType)

Definition at line 237 of file oshGeoMesh.C.

Referenced by NEM::MSH::oshGeoMesh::GM2osh().

237  {
238  switch (vtkCellType) {
239  case VTK_TETRA:
240  case VTK_HEXAHEDRON: return 3;
241 
242  case VTK_TRIANGLE:
243  case VTK_QUAD: return 2;
244 
245  case VTK_LINE: return 1;
246 
247  case VTK_VERTEX: return 0;
248  case VTK_EMPTY_CELL:
249  default: {
250  std::cerr << "Omega_h supports only linear simplices and hypercubes. "
251  "VTKCellType "
252  << vtkCellType << " is not supported." << std::endl;
253  exit(1);
254  }
255  }
256 }

◆ getOmega_hFamilyFromVTKType()

Omega_h_Family NEM::MSH::getOmega_hFamilyFromVTKType ( VTKCellType  vtkCellType)

Definition at line 217 of file oshGeoMesh.C.

Referenced by NEM::MSH::oshGeoMesh::GM2osh().

217  {
218  switch (vtkCellType) {
219  case VTK_EMPTY_CELL:
220  case VTK_VERTEX:
221  case VTK_LINE:
222  case VTK_TRIANGLE:
223  case VTK_TETRA: return OMEGA_H_SIMPLEX;
224 
225  case VTK_QUAD:
226  case VTK_HEXAHEDRON: return OMEGA_H_HYPERCUBE;
227 
228  default: {
229  std::cerr << "Omega_h supports only linear simplices and hypercubes. "
230  "VTKCellType "
231  << vtkCellType << " is not supported." << std::endl;
232  exit(1);
233  }
234  }
235 }

◆ getVtkDataArrayFromOmega_hTag()

void NEM::MSH::getVtkDataArrayFromOmega_hTag ( const Omega_h::TagBase *  tagBase,
vtkSmartPointer< vtkDataArray > &  vtkArray 
)

Definition at line 272 of file oshGeoMesh.C.

References New().

Referenced by NEM::MSH::oshGeoMesh::osh2GM().

273  {
274  switch (tagBase->type()) {
275  case OMEGA_H_I8: {
276  vtkSmartPointer<vtkTypeInt8Array> vtkArray_tmp =
278  Otag2Varray<Omega_h::I8, vtkTypeInt8Array>(tagBase, vtkArray_tmp);
279  vtkArray = vtkArray_tmp;
280  break;
281  }
282  case OMEGA_H_I32: {
283  vtkSmartPointer<vtkTypeInt32Array> vtkArray_tmp =
285  Otag2Varray<Omega_h::I32, vtkTypeInt32Array>(tagBase, vtkArray_tmp);
286  vtkArray = vtkArray_tmp;
287  break;
288  }
289  case OMEGA_H_I64: {
290  vtkSmartPointer<vtkTypeInt64Array> vtkArray_tmp =
292  Otag2Varray<Omega_h::I64, vtkTypeInt64Array>(tagBase, vtkArray_tmp);
293  vtkArray = vtkArray_tmp;
294  break;
295  }
296  case OMEGA_H_F64: {
297  vtkSmartPointer<vtkTypeFloat64Array> vtkArray_tmp =
299  Otag2Varray<Omega_h::Real, vtkTypeFloat64Array>(tagBase, vtkArray_tmp);
300  vtkArray = vtkArray_tmp;
301  break;
302  }
303  }
304 }
geoMeshBase * New(MeshType meshType)
Create a new mesh object.

◆ getVTKTypeFromOmega_hFamilyDim()

VTKCellType NEM::MSH::getVTKTypeFromOmega_hFamilyDim ( Omega_h_Family  family,
Omega_h::Int  dim 
)

Definition at line 201 of file oshGeoMesh.C.

Referenced by NEM::MSH::oshGeoMesh::osh2GM().

202  {
203  return (
204  family == OMEGA_H_SIMPLEX
205  ? (dim == 3 ? VTK_TETRA
206  : (dim == 2 ? VTK_TRIANGLE
207  : (dim == 1 ? VTK_LINE
208  : (dim == 0 ? VTK_VERTEX
209  : VTK_EMPTY_CELL))))
210  : (dim == 3 ? VTK_HEXAHEDRON
211  : (dim == 2 ? VTK_QUAD
212  : (dim == 1 ? VTK_LINE
213  : (dim == 0 ? VTK_VERTEX
214  : VTK_EMPTY_CELL)))));
215 }

◆ is_close()

bool NEM::MSH::is_close ( double  val1,
double  val2,
double  floor,
double  rel_tol 
)

Definition at line 43 of file diffMesh.C.

Referenced by diffMesh().

43  {
44  if (std::fabs(val1) <= floor && std::fabs(val2) <= floor) {
45  return true;
46  } else {
47  return std::fabs((val1 - val2) / val1) < rel_tol;
48  }
49 }

◆ mesh_()

NEM::MSH::mesh_ ( gen_->  CreateMeshfalse)

Definition at line 123 of file smeshGeoMesh.C.

References NEM::MSH::geoMeshBase::report(), NEM::MSH::geoMeshBase::write(), and NEM::MSH::smeshGeoMesh::~smeshGeoMesh().

123 : geoMeshBase(), gen_{new SMESH_Gen{}}, mesh_(gen_->CreateMesh(false)) {}
Mesh mesh_
Definition: inpGeoMesh.C:280

◆ MeshTypeFromFilename()

MeshType NEM::MSH::MeshTypeFromFilename ( const std::string &  fileName)

Definition at line 53 of file geoMeshFactory.C.

References EXO_GEO_MESH, nemAux::find_ext(), FOAM_GEO_MESH, GMSH_GEO_MESH, INP_GEO_MESH, OSH_GEO_MESH, and VTK_GEO_MESH.

Referenced by New(), and Read().

53  {
54  std::string fileExt = nemAux::find_ext(fileName);
55 
56  if (fileExt == ".vtu" || fileExt == ".vtk") {
58  } else if (fileExt == ".msh" || fileExt == ".geo") {
60  } else if (fileExt == ".osh") {
62  } else if (fileExt == ".exo" || fileExt == ".e" || fileExt == ".gen" ||
63  fileExt == ".g") {
65  } else if (fileExt == ".inp") {
67  } else if (fileExt == ".foam") {
69  } else {
70  std::cerr << "File extension " << fileExt << " is not supported."
71  << std::endl;
72  exit(1);
73  }
74 }
OSH_GEO_MESH
Based on Omega_h::Mesh from Omega_h library.
INP_GEO_MESH
Abaqus/CalculiX mesh input format.
std::string find_ext(const std::string &fname)
VTK_GEO_MESH
Based on vtkUnstructuredGrid from VTK library.
FOAM_GEO_MESH
Based on fvMesh from OpenFOAM library.
GMSH_GEO_MESH
Mesh and geometry based on Gmsh&#39;s public API.
EXO_GEO_MESH
Interface to EXODUS II library.

◆ New() [1/2]

NEMOSYS_EXPORT geoMeshBase * NEM::MSH::New ( MeshType  type)
Parameters
typetype of mesh
Returns
pointer to base mesh class

Definition at line 116 of file geoMeshFactory.C.

References EXO_GEO_MESH, FOAM_GEO_MESH, GMSH_GEO_MESH, INP_GEO_MESH, NEM::MSH::gmshGeoMesh::New(), NEM::MSH::vtkGeoMesh::New(), NEM::MSH::smeshGeoMesh::New(), NEM::MSH::foamGeoMesh::New(), NEM::MSH::exoGeoMesh::New(), NEM::MSH::inpGeoMesh::New(), NEM::MSH::oshGeoMesh::New(), OSH_GEO_MESH, SMESH_GEO_MESH, and VTK_GEO_MESH.

Referenced by MeshManipulationFoam::addArtificialThicknessElements(), NEM::MSH::exoGeoMesh::addCellsToBlock(), MeshManipulationFoam::addCohesiveElements(), NEM::MSH::exoGeoMesh::addElemBlock(), NEM::MSH::mergeCells::AddNewCellsDataSet(), NEM::MSH::mergeCells::AddNewCellsUnstructuredGrid(), NEM::MSH::exoGeoMesh::addSideSet(), meshSrch::buildCellLocator(), meshBase::buildStaticCellLocator(), meshBase::buildStaticPointLocator(), Foam::AMRFoam::checkForMotion(), meshSrch::chkDuplElm(), COBALT::cobalt::cobalt(), OrderOfAccuracy::computeDiff(), NEM::ADP::GradSizeField::computeL2GradAtAllCells(), PatchRecovery::computeNodalError(), OrderOfAccuracy::computeRichardsonExtrapolation(), NEM::ADP::Z2ErrorSizeField::computeSizeField(), NEM::ADP::GradSizeField::computeSizeField(), NEM::ADP::ValSizeField::computeSizeField(), ConservativeVolumeTransfer::ConservativeVolumeTransfer(), GaussCubature::constructGaussMesh(), meshBase::convertHexToTetVTK(), meshBase::convertQuads(), ConservativeVolumeTransfer::convertSupermeshToUnstructuredGrid(), NEM::GEO::rocPack::createCohesiveElements(), NEM::GEO::rocPack::createVtkCell(), diffMesh(), NEM::DRV::SmartConversionDriver::execute(), NEM::DRV::CFMeshMeshGenDriver::execute(), NEM::DRV::NucMeshDriver::execute(), NEM::DRV::OmegahRefineDriver::execute(), NEM::DRV::HexPackMeshDriver::execute(), NEM::MSH::exoGeoMesh::exoReader2GM(), meshBase::exportExoToVtk(), meshBase::exportGmshToVtk(), meshBase::exportPntToVtk(), cgnsAnalyzer::exportToVTKMesh(), hdf5Reader::exportToVTKMesh(), meshBase::exportVolToVtk(), meshBase::extractSelectedCells(), vtkMesh::extractSurface(), FOAM::foamMesh::extractSurface(), ConservativeVolumeTransfer::extractSurface(), FETransfer::FETransfer(), meshSrch::FindCellsInPolyData(), meshSrch::FindCellsInSphere(), meshSrch::FindCellsInTriSrf(), meshSrch::FindCellsWithinBounds(), meshSrch::FindPntsOnEdge(), meshSrch::FindPntsOnTriSrf(), NEM::MSH::foamGeoMesh::foam2GM(), NEM::DRV::ConversionDriver::freeSurfaceSideSet(), selectEnclosedPoints::fromDataSet(), NEM::DRV::ConversionDriver::genExo(), NEM::MSH::geoMeshBase::geoMeshBase(), proteusHdf5::getBoundarySideSets(), FOAM::foamMesh::getCell(), vtkMesh::getCell(), meshSrch::getCellVec(), FOAM::foamMesh::getCellVec(), vtkMesh::getCellVec(), FETransfer::getClosestSourceCell(), NEM::MSH::exoGeoMesh::getExoReader(), cgnsAnalyzer::getSectionMesh(), NEM::MSH::exoGeoMesh::getSS(), MeshQuality::getStats(), getVtkDataArrayFromOmega_hTag(), NEM::MSH::foamGeoMesh::GM2foam(), NEM::MSH::oshGeoMesh::GM2osh(), NEM::MSH::gmshGeoMesh::gmsh2GM(), NEM::MSH::inpGeoMesh::inp2GM(), FOAM::foamMesh::inspectEdges(), vtkMesh::inspectEdges(), GaussCubature::integrateOverAllCells(), ConservativeVolumeTransfer::interpolateCellDataToPoints(), GaussCubature::interpolateToGaussPoints(), NEM::MSH::mergeCells::MapPointsToIdsUsingGlobalIds(), NEM::MSH::mergeCells::MapPointsToIdsUsingLocator(), vtkMesh::merge(), NEM::MSH::geoMeshBase::mergeGeoMesh(), NEM::DRV::MeshGenDriver::MeshGenDriver(), meshPartitioner::meshPartitioner(), MeshQuality::MeshQuality(), New(), NEM::MSH::oshGeoMesh::osh2GM(), meshBase::partition(), PATRAN::patran::patran(), MeshManipulationFoam::periodicMeshMapper(), PNTMesh::pntMesh::pntPopulate(), proteusHdf5::proteusHdf5(), gmshMesh::read(), ReadALegacyVTKFile(), ReadAnXMLOrSTLFile(), ReadDegenerateVTKFile(), readLegacyVTKCells(), readLegacyVTKData(), readLegacyVTKFieldData(), NEM::MSH::exoGeoMesh::reconstructGeo(), PatchRecovery::recoverNodalSolution(), NEM::MSH::dataSetRegionBoundaryFilter::RequestData(), NEM::MSH::gmshGeoMesh::resetNative(), NEM::MSH::exoGeoMesh::resetNative(), NEM::MSH::exoGeoMesh::scaleNodes(), vtkMesh::setCellDataArray(), vtkMesh::setPointDataArray(), NEM::MSH::exoGeoMesh::setSideSetObjId(), NEM::MSH::mergeCells::StartUGrid(), NEM::MSH::exoGeoMesh::stitch(), meshBase::stitchMB(), NEM::MSH::oshGeoMesh::takeGeoMesh(), NEM::MSH::exoGeoMesh::takeGeoMesh(), NEM::MSH::geoMeshBase::takeGeoMesh(), ConservativeVolumeTransfer::transfer(), FETransfer::transferCellData(), ConservativeSurfaceTransfer::transferPointData(), FETransfer::transferPointData(), NEM::MSH::vtkGeoMesh::vtk2GM(), vtkMesh::vtkMesh(), vtkStandardNewMacro(), COBALT::cobalt::write(), NEM::MSH::vtkGeoMesh::write(), gmshMesh::write(), FOAM::foamMesh::write(), PATRAN::patran::write2(), PATRAN::patran::write6(), PATRAN::patran::write8(), meshBase::writeCobalt(), and writeVTFile().

116  {
117  switch (meshType) {
120 #ifdef HAVE_GMSH
121  return gmshGeoMesh::New();
122 #else
123  std::cerr << "Please build NEMoSys with ENABLE_GMSH=ON to use this"
124  << " feature!" << std::endl;
125  throw;
126 #endif
127  }
132 #ifdef HAVE_CFMSH
133  return foamGeoMesh::New();
134 #else
135  std::cerr << "Please build NEMoSys with ENABLE_CFMSH=ON to use this"
136  << " feature!" << std::endl;
137  throw;
138 #endif
139  }
141 #ifdef HAVE_OCC
142  return smeshGeoMesh::New();
143 #else
144  std::cerr << "Requires ENABLE_OPENCASCADE=ON\n";
145  throw;
146 #endif
147  }
148  }
149  return nullptr;
150 }
OSH_GEO_MESH
Based on Omega_h::Mesh from Omega_h library.
SMESH_GEO_MESH
Based on SMESH_Mesh from Salome SMESH (see contrib/)
INP_GEO_MESH
Abaqus/CalculiX mesh input format.
VTK_GEO_MESH
Based on vtkUnstructuredGrid from VTK library.
geoMeshBase * New(const std::string &fileName)
Create a new mesh object.
FOAM_GEO_MESH
Based on fvMesh from OpenFOAM library.
GMSH_GEO_MESH
Mesh and geometry based on Gmsh&#39;s public API.
EXO_GEO_MESH
Interface to EXODUS II library.

◆ New() [2/2]

NEMOSYS_EXPORT geoMeshBase * NEM::MSH::New ( const std::string &  fileName)

Determines mesh type from file extension.

Parameters
fileNamename of file
Returns
pointer to base mesh class

Definition at line 152 of file geoMeshFactory.C.

References MeshTypeFromFilename(), and New().

152  {
153  return New(MeshTypeFromFilename(fileName));
154 }
MeshType MeshTypeFromFilename(const std::string &fileName)
geoMeshBase * New(const std::string &fileName)
Create a new mesh object.

◆ oshFace2vtkFace()

int NEM::MSH::oshFace2vtkFace ( int  oshFace,
Omega_h_Family  oshFamily,
Omega_h::Int  oshDim 
)

Definition at line 306 of file oshGeoMesh.C.

Referenced by NEM::MSH::oshGeoMesh::osh2GM().

307  {
308  switch (oshDim) {
309  case 2: return oshFace;
310  case 3:
311  switch (oshFamily) {
312  case OMEGA_H_SIMPLEX: return (oshFace + 3) % 4;
313  case OMEGA_H_HYPERCUBE:
314  switch (oshFace) {
315  case 0: return 4;
316  case 1: return 2;
317  case 2: return 1;
318  case 3: return 3;
319  case 4: return 0;
320  case 5: return 5;
321  default: return -1;
322  }
323  }
324  default: return -1;
325  }
326 }

◆ Otag2Varray()

template<typename OT , typename VT >
void NEM::MSH::Otag2Varray ( const Omega_h::TagBase *  tagBase,
VT *  vtkArray 
)

Definition at line 259 of file oshGeoMesh.C.

259  {
260  const Omega_h::Tag<OT> *tag = Omega_h::as<OT>(tagBase);
261  Omega_h::Int ncomps = tagBase->ncomps();
262  Omega_h::HostRead<OT> tagArray(tag->array());
263 
264  vtkArray->SetName(tagBase->name().c_str());
265  vtkArray->SetNumberOfComponents(ncomps);
266  vtkArray->SetNumberOfValues(tagArray.size());
267 
268  for (Omega_h::LO i = 0; i < tagArray.size(); ++i)
269  vtkArray->SetValue(i, tagArray.get(i));
270 }

◆ Read() [1/2]

NEMOSYS_EXPORT geoMeshBase * NEM::MSH::Read ( const std::string &  fileName)

Determines mesh type from file extension.

Parameters
fileNamename of file
Returns
pointer to base mesh class

Definition at line 76 of file geoMeshFactory.C.

References MeshTypeFromFilename().

Referenced by MeshManipulationFoam::addArtificialThicknessElements(), MeshManipulationFoam::addCohesiveElements(), cfmeshGen::createMeshFromSTL(), snappymeshGen::createMeshFromSTL(), blockMeshGen::createMeshFromSTL(), NEM::DRV::SmartConversionDriver::execute(), NEM::DRV::OmegahRefineDriver::execute(), NEM::DRV::HexPackMeshDriver::execute(), NEM::MSH::exoGeoMesh::getSS(), and MeshManipulationFoam::periodicMeshMapper().

76  {
77  return Read(fileName, MeshTypeFromFilename(fileName));
78 }
MeshType MeshTypeFromFilename(const std::string &fileName)
geoMeshBase * Read(const std::string &fileName, MeshType meshType)
Read a mesh from file.

◆ Read() [2/2]

NEMOSYS_EXPORT geoMeshBase * NEM::MSH::Read ( const std::string &  fileName,
MeshType  type 
)
Parameters
fileNamename of file
typetype of mesh
Returns
pointer to base mesh class

Definition at line 80 of file geoMeshFactory.C.

References EXO_GEO_MESH, FOAM_GEO_MESH, GMSH_GEO_MESH, INP_GEO_MESH, OSH_GEO_MESH, NEM::MSH::inpGeoMesh::Read(), SMESH_GEO_MESH, and VTK_GEO_MESH.

80  {
81  switch (meshType) {
82  case MeshType::VTK_GEO_MESH: return vtkGeoMesh::Read(fileName);
84 #ifdef HAVE_GMSH
85  return gmshGeoMesh::Read(fileName);
86 #else
87  std::cerr << "Please build NEMoSys with ENABLE_GMSH=ON to use this"
88  << " feature!" << std::endl;
89  throw;
90 #endif
91  }
92  case MeshType::OSH_GEO_MESH: return oshGeoMesh::Read(fileName);
93  case MeshType::EXO_GEO_MESH: return exoGeoMesh::Read(fileName);
94  case MeshType::INP_GEO_MESH: return inpGeoMesh::Read(fileName);
96 #ifdef HAVE_CFMSH
97  std::string fname = fileName;
98  fname.erase(fname.find_last_of('.'));
99  return foamGeoMesh::Read(fname);
100 #else
101  std::cerr << "Please build NEMoSys with ENABLE_CFMSH=ON to use this"
102  << " feature!" << std::endl;
103  throw;
104 #endif
105  }
107  std::cerr << "Unsupported.\n";
108  throw;
109  }
110  }
111  // Don't put inside switch statement so static analyzers can detect unhandled
112  // cases
113  return nullptr;
114 }
OSH_GEO_MESH
Based on Omega_h::Mesh from Omega_h library.
SMESH_GEO_MESH
Based on SMESH_Mesh from Salome SMESH (see contrib/)
INP_GEO_MESH
Abaqus/CalculiX mesh input format.
VTK_GEO_MESH
Based on vtkUnstructuredGrid from VTK library.
FOAM_GEO_MESH
Based on fvMesh from OpenFOAM library.
GMSH_GEO_MESH
Mesh and geometry based on Gmsh&#39;s public API.
geoMeshBase * Read(const std::string &fileName, MeshType meshType)
Read a mesh from file.
EXO_GEO_MESH
Interface to EXODUS II library.

◆ Varray2Otag()

template<typename VT , typename OT >
void NEM::MSH::Varray2Otag ( Omega_h::Mesh *  oshMesh,
VT *  vtkArray,
Omega_h::Int  oshDim 
)

Definition at line 518 of file oshGeoMesh.C.

Referenced by Varray2Otag2().

518  {
519  Omega_h::HostWrite<OT> h_oshArray;
520  h_oshArray = Omega_h::HostWrite<OT>(vtkArray->GetNumberOfValues(),
521  vtkArray->GetClassName());
522 
523  for (vtkIdType i = 0; i < vtkArray->GetNumberOfValues(); ++i)
524  h_oshArray.set(i, vtkArray->GetValue(i));
525 
526  Omega_h::Read<OT> oshArray(h_oshArray);
527 
528  oshMesh->add_tag(oshDim, vtkArray->GetName(),
529  vtkArray->GetNumberOfComponents(), oshArray);
530 }

◆ Varray2Otag2()

template<typename VT >
void NEM::MSH::Varray2Otag2 ( Omega_h::Mesh *  oshMesh,
VT *  vtkArray,
Omega_h::Int  oshDim 
)

Definition at line 533 of file oshGeoMesh.C.

References Varray2Otag().

533  {
534  Varray2Otag<
535  VT, typename std::conditional<
536  sizeof(typename VT::ValueType) <= 1, Omega_h::I8,
537  typename std::conditional<
538  sizeof(typename VT::ValueType) <= 4, Omega_h::I32,
539  typename std::enable_if<sizeof(typename VT::ValueType) <= 8,
540  Omega_h::I64>::type>::type>::type>(
541  oshMesh, vtkArray, oshDim);
542 }
void Varray2Otag(Omega_h::Mesh *oshMesh, VT *vtkArray, Omega_h::Int oshDim)
Definition: oshGeoMesh.C:518

◆ vtkStandardNewMacro() [1/9]

NEM::MSH::vtkStandardNewMacro ( vtkGeoMesh  )

Definition at line 54 of file vtkGeoMesh.C.

References nemAux::find_ext(), New(), and NEM::MSH::vtkGeoMesh::vtkGeoMesh().

57  {
58  vtkSmartPointer<vtkUnstructuredGrid> inUnstructuredGrid;
59 
60  std::string fileExt = nemAux::find_ext(fileName);
61  if (fileExt == ".vtk") { // Legacy type
62  vtkSmartPointer<vtkGenericDataObjectReader> reader =
64  reader->SetFileName(fileName.c_str());
65  reader->Update();
66 
67  inUnstructuredGrid = vtkUnstructuredGrid::SafeDownCast(reader->GetOutput());
68  } else if (fileExt.substr(0, 3) == ".vt") { // XML type
69  vtkSmartPointer<vtkXMLGenericDataObjectReader> reader =
71  reader->SetFileName(fileName.c_str());
72  reader->Update();
73 
74  inUnstructuredGrid = vtkUnstructuredGrid::SafeDownCast(reader->GetOutput());
75  /*
76  } else if (fileExt == ".stl") { // Allow STLs
77  vtkSmartPointer<vtkSTLReader> reader = vtkSmartPointer<vtkSTLReader>::New();
78  reader->SetFileName(fileName.c_str());
79  reader->Update();
80 
81  // TODO: This cast fails. STL reader outputs vtkPolyData!
82  inUnstructuredGrid = vtkUnstructuredGrid::SafeDownCast(reader->GetOutput());
83  */
84  } else { // Not a VTK type
85  std::cerr << "ERROR: File extension " << fileExt
86  << " is not supported by vtkGeoMesh." << std::endl;
87  exit(1);
88  }
89 
90  if (!inUnstructuredGrid) { // Not a mesh
91  std::cerr << "ERROR: Only VTK Unstructured Grid (vtu) is supported."
92  << std::endl;
93  exit(1);
94  }
95 
96  return new vtkGeoMesh(inUnstructuredGrid, geoArrayName);
97 }
geoMeshBase * New(MeshType meshType)
Create a new mesh object.
std::string find_ext(const std::string &fname)

◆ vtkStandardNewMacro() [2/9]

NEM::MSH::vtkStandardNewMacro ( mergeCells  )

Definition at line 59 of file mergeCells.C.

61  {
62  public:
63  std::map<vtkIdType, vtkIdType> IdTypeMap;
64 };

◆ vtkStandardNewMacro() [3/9]

NEM::MSH::vtkStandardNewMacro ( foamGeoMesh  )

Definition at line 63 of file foamGeoMesh.C.

References NEM::MSH::foamGeoMesh::foam2GM(), and NEM::MSH::foamGeoMesh::foamGeoMesh().

66  {
67  auto foamGM = new foamGeoMesh();
68  Foam::word regionName;
69 
70  if (fileName.empty()) {
71  regionName = Foam::fvMesh::defaultRegion;
72  } else {
73  regionName = fileName;
74  }
75 
76  bool writeDicts = false;
77  std::unique_ptr<getDicts> initFoam;
78  initFoam = std::unique_ptr<getDicts>(new getDicts());
79  foamGM->controlDict_ = initFoam->createControlDict(writeDicts);
80 
81  foamGM->runTime_ = std::unique_ptr<Foam::Time>(
82  new Foam::Time(foamGM->controlDict_.get(), ".", "."));
83 
84  foamGM->fmesh_.reset(new Foam::fvMesh(
85  Foam::IOobject(regionName, foamGM->runTime_->timeName(),
86  *foamGM->runTime_, Foam::IOobject::MUST_READ)));
87 
88  auto gmMesh = foam2GM(foamGM->fmesh_.get());
89 
90  foamGM->setGeoMesh(gmMesh);
91 
92  return foamGM;
93 }

◆ vtkStandardNewMacro() [4/9]

NEM::MSH::vtkStandardNewMacro ( oshGeoMesh  )

Definition at line 73 of file oshGeoMesh.C.

References NEM::MSH::OmegaHInterface::GetLibrary().

76  {
77  if (!lib) {
78  lib = OmegaHInterface::GetLibrary().get();
79  }
80  Omega_h::Mesh oshMesh = Omega_h::read_mesh_file(fileName, lib->world());
81 
82  return new oshGeoMesh(&oshMesh, lib);
83 }

◆ vtkStandardNewMacro() [5/9]

NEM::MSH::vtkStandardNewMacro ( gmshGeoMesh  )

Definition at line 98 of file gmshGeoMesh.C.

References NEM::MSH::gmshGeoMesh::gmshGeoMesh(), and NEM::MSH::GmshInterface::Initialize().

100  {
101  GmshInterface::Initialize();
102 
103  gmsh::open(fileName);
104 
105  std::string gmshMesh;
106  gmsh::model::getCurrent(gmshMesh);
107 
108  return new gmshGeoMesh(gmshMesh);
109 }

◆ vtkStandardNewMacro() [6/9]

NEM::MSH::vtkStandardNewMacro ( smeshGeoMesh  )

Definition at line 120 of file smeshGeoMesh.C.

123  : geoMeshBase(), gen_{new SMESH_Gen{}}, mesh_(gen_->CreateMesh(false)) {}
Mesh mesh_
Definition: inpGeoMesh.C:280

◆ vtkStandardNewMacro() [7/9]

NEM::MSH::vtkStandardNewMacro ( dataSetRegionBoundaryFilter  )

Definition at line 246 of file dataSetRegionBoundaryFilter.C.

249  {
250  info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkDataSet");
251  return 1;
252 }

◆ vtkStandardNewMacro() [8/9]

NEM::MSH::vtkStandardNewMacro ( exoGeoMesh  )

Definition at line 247 of file exoGeoMesh.C.

References NEM::MSH::exoGeoMesh::exoGeoMesh().

249  {
250  return new exoGeoMesh(fileName, timeStep);
251 }

◆ vtkStandardNewMacro() [9/9]

NEM::MSH::vtkStandardNewMacro ( inpGeoMesh  )

Definition at line 593 of file inpGeoMesh.C.

595  {
596  return new inpGeoMesh(fileName);
597 }