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.
COBALT::cobalt Class Reference

Detailed Description

Definition at line 37 of file cobalt.H.

Public Member Functions

 cobalt (std::shared_ptr< meshBase > fullMesh, const std::string &inFnameVtk, const std::string &outFnameCgr, const std::string &outFnameCgi)
 
void write () const
 
void writePatchMap (const std::string &mapFile, const std::map< nemId_t, nemId_t > &patchMap) const
 
void writePatchMap (std::ofstream &outputStream, const std::map< nemId_t, nemId_t > &patchMap) const
 

Private Attributes

std::shared_ptr< meshBasevolMeshBase
 
std::shared_ptr< meshBasesurfMeshBase
 
std::string inFnameVtk
 
std::string outFnameCgr
 
std::string outFnameCgi
 

Constructor & Destructor Documentation

◆ cobalt()

COBALT::cobalt::cobalt ( std::shared_ptr< meshBase fullMesh,
const std::string &  inFnameVtk,
const std::string &  outFnameCgr,
const std::string &  outFnameCgi 
)

Definition at line 48 of file cobalt.C.

References meshBase::CreateShared(), inFnameVtk, NEM::MSH::New(), outFnameCgi, outFnameCgr, surfMeshBase, and volMeshBase.

52 {
53  inFnameVtk = _inFnameVtk; // vtk input file
54  outFnameCgr = _outFnameCgr; // cobalt output grid file
55  outFnameCgi = _outFnameCgi; // cobalt output patch mapping file
56 
57  // declare array to store element type
58  vtkSmartPointer<vtkDataArray> elementTypeArray
60  elementTypeArray->SetNumberOfComponents(1);
61  elementTypeArray->SetName("elemType");
62 
63  // get element types
64  double tmp[1]; // storing inserted tuples
65  for (int iCell = 0;
66  iCell < fullMesh->getDataSet()->GetNumberOfCells(); iCell++)
67  {
68  if (fullMesh->getDataSet()->GetCellType(iCell) != VTK_TETRA
69  && fullMesh->getDataSet()->GetCellType(iCell) != VTK_TRIANGLE)
70  {
71  std::cout << "Error: only triangle and tetrahedral elements supported."
72  << std::endl;
73  exit(1);
74  }
75  tmp[0] = fullMesh->getDataSet()->GetCellType(iCell);
76  elementTypeArray->InsertNextTuple(tmp);
77  }
78  fullMesh->getDataSet()->GetCellData()->AddArray(elementTypeArray);
79 
80  // threshold to extract only volume cells
81  vtkSmartPointer<vtkThreshold> volThreshold = vtkSmartPointer<vtkThreshold>::New();
82  volThreshold->SetInputData(fullMesh->getDataSet());
83  volThreshold->SetInputArrayToProcess(0, 0, 0,
84  vtkDataObject::FIELD_ASSOCIATION_CELLS,
85  "elemType");
86  volThreshold->ThresholdBetween(VTK_TETRA, VTK_TETRA);
87  volThreshold->Update();
88  // store output in unstructured grid
89  vtkSmartPointer<vtkUnstructuredGrid> volUG = volThreshold->GetOutput();
90  // create meshBase object
91  volMeshBase = meshBase::CreateShared(volUG, "extractedVolume.vtu");
92  //volUG_mb->write(); // write to file
93 
94  // threshold to extract only surface cells
95  vtkSmartPointer<vtkThreshold> surfThreshold = vtkSmartPointer<vtkThreshold>::New();
96  surfThreshold->SetInputData(fullMesh->getDataSet());
97  surfThreshold->SetInputArrayToProcess(0, 0, 0,
98  vtkDataObject::FIELD_ASSOCIATION_CELLS,
99  "elemType");
100  surfThreshold->ThresholdBetween(VTK_TRIANGLE, VTK_TRIANGLE);
101  surfThreshold->Update();
102  // store output in unstructured grid
103  vtkSmartPointer<vtkUnstructuredGrid> surfUG = surfThreshold->GetOutput();
104  // geometry filter to polyData object
105  vtkSmartPointer<vtkGeometryFilter> geomFilter = vtkGeometryFilter::New();
106  geomFilter->SetInputData(surfUG);
107  geomFilter->Update();
108  vtkSmartPointer<vtkPolyData> surfPD = geomFilter->GetOutput();
109  // create meshBase object
110  surfMeshBase = meshBase::CreateShared(surfPD, "extractedSurface.vtp");
111  //surfPD_mb->write(); // write to file
112 }
std::shared_ptr< meshBase > volMeshBase
Definition: cobalt.H:53
std::string outFnameCgr
Definition: cobalt.H:56
std::string inFnameVtk
Definition: cobalt.H:55
geoMeshBase * New(MeshType meshType)
Create a new mesh object.
static std::shared_ptr< meshBase > CreateShared(const std::string &fname)
Create shared ptr from fname.
Definition: meshBase.C:171
std::string outFnameCgi
Definition: cobalt.H:57
std::shared_ptr< meshBase > surfMeshBase
Definition: cobalt.H:54

Member Function Documentation

◆ write()

void COBALT::cobalt::write ( ) const

Definition at line 115 of file cobalt.C.

References nemAux::l2_Norm(), NEM::MSH::New(), outFnameCgi, outFnameCgr, surfMeshBase, volMeshBase, and writePatchMap().

Referenced by NEM::DRV::VtkToCobaltConversionDriver::execute().

116 {
117  std::ofstream outputStream(outFnameCgr);
118  if (!outputStream.good())
119  {
120  std::cout << "Cannot open file " << outFnameCgr << std::endl;
121  exit(1);
122  }
123 
124  if (!surfMeshBase)
125  {
126  std::cout << "surface mesh is empty!" << std::endl;
127  exit(1);
128  }
129  auto surfPatchNoArr =
130  surfMeshBase->getDataSet()->GetCellData()->GetArray("patchNo");
131  if (!surfPatchNoArr) {
132  surfPatchNoArr =
133  surfMeshBase->getDataSet()->GetCellData()->GetArray("PhysGrpId");
134  if (!surfPatchNoArr) {
135  std::cout << "surface mesh must have cell array with name \"patchNo\" or "
136  "\"PhysGrpId\" denoting patch number."
137  << std::endl;
138  exit(1);
139  }
140  }
141 
142  vtkSmartPointer<vtkIdList> facePtIds;
143  vtkSmartPointer<vtkIdList> sharedCellPtIds = vtkSmartPointer<vtkIdList>::New();
144  vtkSmartPointer<vtkGenericCell> genCell1 = vtkSmartPointer<vtkGenericCell>::New();
145  vtkSmartPointer<vtkGenericCell> genCell2 = vtkSmartPointer<vtkGenericCell>::New();
146  // Signed to represent patch numbers as negative integers
147  std::map<std::vector<nemId_t>,
148  std::pair<nemId_t, long long int>,
149  sortNemId_tVec_compare> faceMap;
150  // building cell locator for looking up patch number in remeshed surface mesh
151  vtkSmartPointer<vtkStaticCellLocator> surfCellLocator = surfMeshBase->buildStaticCellLocator();
152  // maximum number of vertices per face (to be found in proceeding loop)
153  vtkIdType nVerticesPerFaceMax = 0;
154  // maximum number of faces per cell (to be found in proceeding loop)
155  int nFacesPerCellMax = 0;
156 
157  for (vtkIdType i = 0; i < volMeshBase->getNumberOfCells(); ++i)
158  {
159  // get cell i
160  volMeshBase->getDataSet()->GetCell(i, genCell1);
161  // get faces, find cells sharing it. if no cell shares it,
162  // use the locator of the surfMeshBase to find the patch number
163  int numFaces = genCell1->GetNumberOfFaces();
164  nFacesPerCellMax = (nFacesPerCellMax < numFaces
165  ? numFaces
166  : nFacesPerCellMax);
167  for (int j = 0; j < numFaces; ++j)
168  {
169  vtkCell *face = genCell1->GetFace(j);
170  // bool shared = false;
171  vtkIdType numVerts = face->GetNumberOfPoints();
172  nVerticesPerFaceMax = (nVerticesPerFaceMax < numVerts
173  ? numVerts
174  : nVerticesPerFaceMax);
175  facePtIds = face->GetPointIds();
176  volMeshBase->getDataSet()->GetCellNeighbors(i, facePtIds,
177  sharedCellPtIds);
178  std::vector<nemId_t> facePntIds(numVerts);
179  for (vtkIdType k = 0; k < numVerts; ++k)
180  {
181  facePntIds[k] = face->GetPointId(k) + 1;
182  }
183  if (sharedCellPtIds->GetNumberOfIds())
184  {
185  faceMap.emplace(facePntIds, std::pair<nemId_t, long long int>(
186  i + 1, sharedCellPtIds->GetId(0) + 1));
187  }
188  else
189  {
190  double bounds[6];
191  face->GetBounds(bounds);
192  auto cellIds = vtkSmartPointer<vtkIdList>::New();
193  surfCellLocator->FindCellsWithinBounds(bounds, cellIds);
194  vtkIdType closestCellId;
195  if (cellIds->GetNumberOfIds() == 0) {
196  std::cerr << "No surface cell found." << std::endl;
197  continue;
198  } else {
199  closestCellId = cellIds->GetId(0);
200  if (cellIds->GetNumberOfIds() > 1) {
201  double p1[3], p2[3], p3[3];
202  face->GetPoints()->GetPoint(0, p1);
203  face->GetPoints()->GetPoint(1, p2);
204  face->GetPoints()->GetPoint(2, p3);
205  std::vector<double> faceCenter(3);
206  for (vtkIdType k = 0; k < numVerts; ++k) {
207  faceCenter[k] = (p1[k] + p2[k] + p3[k]) / 3.0;
208  }
209  double minDist = nemAux::l2_Norm(
210  faceCenter - surfMeshBase->getCellCenter(cellIds->GetId(0)));
211  for (vtkIdType k = 1; k < cellIds->GetNumberOfIds(); ++k) {
212  auto candidateCellCenter =
213  surfMeshBase->getCellCenter(cellIds->GetId(k));
214  auto dist = nemAux::l2_Norm(faceCenter - candidateCellCenter);
215  if (dist < minDist) {
216  minDist = dist;
217  closestCellId = cellIds->GetId(k);
218  }
219  }
220  }
221  }
222  double patchNo = surfPatchNoArr->GetComponent(closestCellId, 0);
223  faceMap.emplace(facePntIds,
224  std::pair<nemId_t, long long int>(i + 1, -1 * patchNo));
225  }
226  }
227  }
228 
229  std::map<nemId_t, nemId_t> patchMap;
230  for (nemId_t i = 0; i < surfMeshBase->getNumberOfCells(); ++i)
231  {
232  double patchNo[1];
233  surfPatchNoArr->GetTuple(i, patchNo);
234  patchMap.insert(std::pair<nemId_t, nemId_t>(patchNo[0], i));
235  }
236 
237  // write patch mapping file
238  writePatchMap(outFnameCgi, patchMap);
239  // write cobalt file
240  outputStream << 3 << " " << 1 << " " << patchMap.size() << std::endl;
241  outputStream << volMeshBase->getNumberOfPoints() << " " << faceMap.size()
242  << " " << volMeshBase->getNumberOfCells() << " "
243  << nVerticesPerFaceMax << " " << nFacesPerCellMax << std::endl;
244  for (nemId_t i = 0; i < volMeshBase->getNumberOfPoints(); ++i)
245  {
246  std::vector<double> pnt(volMeshBase->getPoint(i));
247  outputStream << std::setw(21) << std::fixed << std::setprecision(15)
248  << pnt[0] << " " << pnt[1] << " " << pnt[2] << std::endl;
249  }
250 
251  auto it = faceMap.begin();
252  while (it != faceMap.end())
253  {
254  outputStream << it->first.size() << " ";
255  auto faceIdIter = it->first.begin();
256  while (faceIdIter != it->first.end())
257  {
258  outputStream << *faceIdIter << " ";
259  ++faceIdIter;
260  }
261  outputStream << it->second.first << " " << it->second.second << std::endl;
262  ++it;
263  }
264 }
std::shared_ptr< meshBase > volMeshBase
Definition: cobalt.H:53
std::string outFnameCgr
Definition: cobalt.H:56
std::size_t nemId_t
Definition: meshBase.H:51
void writePatchMap(const std::string &mapFile, const std::map< nemId_t, nemId_t > &patchMap) const
Definition: cobalt.C:267
geoMeshBase * New(MeshType meshType)
Create a new mesh object.
sum comparison for vectors representing faces inserted into map
Definition: meshBase.H:765
std::string outFnameCgi
Definition: cobalt.H:57
std::shared_ptr< meshBase > surfMeshBase
Definition: cobalt.H:54
T l2_Norm(const std::vector< T > &x)

◆ writePatchMap() [1/2]

void COBALT::cobalt::writePatchMap ( const std::string &  mapFile,
const std::map< nemId_t, nemId_t > &  patchMap 
) const

Definition at line 267 of file cobalt.C.

Referenced by write().

269 {
270  std::ofstream outputStream(mapFile);
271  if (!outputStream.good())
272  {
273  std::cout << "Error opening file " << mapFile << std::endl;
274  exit(1);
275  }
276  writePatchMap(outputStream, patchMap);
277 }
void writePatchMap(const std::string &mapFile, const std::map< nemId_t, nemId_t > &patchMap) const
Definition: cobalt.C:267

◆ writePatchMap() [2/2]

void COBALT::cobalt::writePatchMap ( std::ofstream &  outputStream,
const std::map< nemId_t, nemId_t > &  patchMap 
) const

Definition at line 280 of file cobalt.C.

282 {
283  outputStream << patchMap.size() << std::endl;
284  outputStream << patchMap.size() << std::endl;
285  auto it = patchMap.begin();
286  size_t normPatchNo = 1;
287  while (it != patchMap.end())
288  {
289  for (int i = 0; i < 2; ++i)
290  {
291  outputStream << std::setw(2) << std::left << it->first << " ";
292  }
293  outputStream << std::setw(2) << std::left << normPatchNo << " ";
294  outputStream << std::endl;
295  ++it;
296  normPatchNo++;
297  }
298 }

Member Data Documentation

◆ inFnameVtk

std::string COBALT::cobalt::inFnameVtk
private

Definition at line 55 of file cobalt.H.

Referenced by cobalt().

◆ outFnameCgi

std::string COBALT::cobalt::outFnameCgi
private

Definition at line 57 of file cobalt.H.

Referenced by cobalt(), and write().

◆ outFnameCgr

std::string COBALT::cobalt::outFnameCgr
private

Definition at line 56 of file cobalt.H.

Referenced by cobalt(), and write().

◆ surfMeshBase

std::shared_ptr<meshBase> COBALT::cobalt::surfMeshBase
private

Definition at line 54 of file cobalt.H.

Referenced by cobalt(), and write().

◆ volMeshBase

std::shared_ptr<meshBase> COBALT::cobalt::volMeshBase
private

Definition at line 53 of file cobalt.H.

Referenced by cobalt(), and write().


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