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.
meshBase.H File Reference

Go to the source code of this file.

Classes

class  meshBase
 A brief description of meshBase. More...
 
struct  sortNemId_tVec_compare
 sum comparison for vectors representing faces inserted into map More...
 

Typedefs

using nemId_t = std::size_t
 

Functions

NEMOSYS_EXPORT int diffMesh (meshBase *mesh1, meshBase *mesh2)
 compares two meshBase classes. More...
 
NEMOSYS_EXPORT void writePatchMap (const std::string &mapFile, const std::map< int, int > &patchMap)
 write patch map file for roc prep (trivial identity mapping) More...
 
NEMOSYS_EXPORT void writePatchMap (std::ofstream &outputStream, const std::map< int, int > &patchMap)
 write patch map file for roc prep (trivial identity mapping) More...
 

Typedef Documentation

◆ nemId_t

using nemId_t = std::size_t

Definition at line 51 of file meshBase.H.

Function Documentation

◆ diffMesh()

NEMOSYS_EXPORT int diffMesh ( meshBase mesh1,
meshBase mesh2 
)

used in testing

Parameters
mesh1<>
mesh2<>

Definition at line 1665 of file meshBase.C.

References meshBase::getCellVec(), meshBase::getDataSet(), meshBase::getNumberOfCells(), meshBase::getNumberOfPoints(), meshBase::getPoint(), and NEM::MSH::New().

1665  {
1666  // double tol = 1e-14;
1667  double tol = 3e-2;
1668 
1669  if (mesh1->getNumberOfPoints() != mesh2->getNumberOfPoints() ||
1670  mesh1->getNumberOfCells() != mesh2->getNumberOfCells()) {
1671  std::cerr << "Meshes don't have the same number of points or cells"
1672  << std::endl;
1673  return 1;
1674  }
1675 
1676  for (int i = 0; i < mesh1->getNumberOfPoints(); ++i) {
1677  std::vector<double> coord1 = mesh1->getPoint(i);
1678  std::vector<double> coord2 = mesh2->getPoint(i);
1679  for (int j = 0; j < 3; ++j) {
1680  if (std::fabs((coord1[j] - coord2[j]) / coord2[j]) > tol) {
1681  std::cerr << "Meshes differ in point coordinates" << std::endl;
1682  std::cerr << "Index " << i << " Component " << j << std::endl;
1683  std::cerr << "Coord 1 " << std::setprecision(15) << coord1[j]
1684  << " Coord 2 " << std::setprecision(15) << coord2[j]
1685  << std::endl;
1686  std::cerr << "Meshes differ in point coordinates" << std::endl;
1687  return 1;
1688  }
1689  }
1690  }
1691 
1692  for (int i = 0; i < mesh1->getNumberOfCells(); ++i) {
1693  std::vector<std::vector<double>> cell1 = mesh1->getCellVec(i);
1694  std::vector<std::vector<double>> cell2 = mesh2->getCellVec(i);
1695  if (cell1.size() != cell2.size()) {
1696  std::cerr << "Meshes differ in cells" << std::endl;
1697  return 1;
1698  }
1699  for (int j = 0; j < cell1.size(); ++j) {
1700  for (int k = 0; k < 3; ++k) {
1701  if (std::fabs((cell1[j][k] - cell2[j][k]) / cell2[j][k]) > tol) {
1702  std::cerr << "Meshes differ in cells" << std::endl;
1703  return 1;
1704  }
1705  }
1706  }
1707  }
1708 
1709  vtkSmartPointer<vtkPointData> pd1 = vtkSmartPointer<vtkPointData>::New();
1710  pd1 = mesh1->getDataSet()->GetPointData();
1711  vtkSmartPointer<vtkPointData> pd2 = vtkSmartPointer<vtkPointData>::New();
1712  pd2 = mesh2->getDataSet()->GetPointData();
1713  int numArr1 = pd1->GetNumberOfArrays();
1714  int numArr2 = pd2->GetNumberOfArrays();
1715 
1716  if (numArr1 != numArr2) {
1717  std::cerr << "Meshes have different numbers of point data" << std::endl;
1718  std::cerr << "Mesh 1 has " << numArr1 << std::endl;
1719  std::cerr << "Mesh 2 has " << numArr2 << std::endl;
1720  return 1;
1721  }
1722 
1723  for (int i = 0; i < numArr1; ++i) {
1724  vtkDataArray *da1 = pd1->GetArray(i);
1725  std::cerr << "checking array " << da1->GetName() << std::endl;
1726  vtkDataArray *da2 = pd2->GetArray(pd1->GetArrayName(i));
1727  double range[2];
1728  double abs_error;
1729  double rel_error;
1730  int numComponent = da1->GetNumberOfComponents();
1731  for (int j = 0; j < mesh1->getNumberOfPoints(); ++j) {
1732  double *comps1 = new double[numComponent];
1733  double *comps2 = new double[numComponent];
1734  da1->GetTuple(j, comps1);
1735  da2->GetTuple(j, comps2);
1736  for (int k = 0; k < numComponent; ++k) {
1737  da1->GetRange(range, k);
1738  abs_error = std::fabs(comps1[k] - comps2[k]);
1739  double max_val = std::max(std::abs(range[0]), std::abs(range[1]));
1740  rel_error = abs_error / std::max(1.0, max_val);
1741  if (rel_error > tol) {
1742  std::cerr << "For point data array " << da1->GetName() << std::endl;
1743  std::cerr << "Meshes differ in point data values at point " << j
1744  << " component " << k << std::endl;
1745  std::cerr << std::setprecision(15) << "Mesh 1 value : " << comps1[k]
1746  << std::endl;
1747  std::cerr << std::setprecision(15) << "Mesh 2 value : " << comps2[k]
1748  << std::endl;
1749  return 1;
1750  }
1751  }
1752  delete[] comps1;
1753  delete[] comps2;
1754  }
1755  }
1756  std::cerr << "Meshes are the same" << std::endl;
1757  return 0;
1758 }
geoMeshBase * New(MeshType meshType)
Create a new mesh object.
virtual std::vector< std::vector< double > > getCellVec(nemId_t id) const =0
get vector of coords of cell with id
vtkSmartPointer< vtkDataSet > getDataSet() const
get this meshes&#39; dataSet
Definition: meshBase.H:308
nemId_t getNumberOfPoints() const
return the number of points
Definition: meshBase.H:545
nemId_t getNumberOfCells() const
return the number of cells
Definition: meshBase.H:550
virtual std::vector< double > getPoint(nemId_t id) const =0
get point with id

◆ writePatchMap() [1/2]

NEMOSYS_EXPORT void writePatchMap ( const std::string &  mapFile,
const std::map< int, int > &  patchMap 
)
Parameters
mapFile<>
patchMap<>

Definition at line 1376 of file meshBase.C.

References writePatchMap().

Referenced by meshBase::writeCobalt(), and writePatchMap().

1377  {
1378  std::ofstream outputStream(mapFile);
1379  if (!outputStream.good()) {
1380  std::cout << "Error opening file " << mapFile << std::endl;
1381  exit(1);
1382  }
1383  writePatchMap(outputStream, patchMap);
1384 }
void writePatchMap(const std::string &mapFile, const std::map< int, int > &patchMap)
write patch map file for roc prep (trivial identity mapping)
Definition: meshBase.C:1376

◆ writePatchMap() [2/2]

NEMOSYS_EXPORT void writePatchMap ( std::ofstream &  outputStream,
const std::map< int, int > &  patchMap 
)
Parameters
outputStream<>
patchMap<>

Definition at line 1388 of file meshBase.C.

1389  {
1390  outputStream << patchMap.size() << std::endl;
1391  outputStream << patchMap.size() << std::endl;
1392  auto it = patchMap.begin();
1393  int normPatchNo = 1;
1394  while (it != patchMap.end()) {
1395  for (int i = 0; i < 2; ++i) {
1396  outputStream << std::setw(2) << std::left << it->first << " ";
1397  }
1398  outputStream << std::setw(2) << std::left << normPatchNo << " ";
1399  outputStream << std::endl;
1400  ++it;
1401  normPatchNo++;
1402  }
1403 }