ViennaLS
Loading...
Searching...
No Matches
lsWriteVisualizationMesh.hpp
Go to the documentation of this file.
1#pragma once
2
3#ifdef VIENNALS_USE_VTK // this class needs vtk support
4
5#include <vtkAppendFilter.h>
6#include <vtkAppendPolyData.h>
7#include <vtkCellData.h>
8#include <vtkDataSetTriangleFilter.h>
9#include <vtkFloatArray.h>
10#include <vtkGeometryFilter.h>
11#include <vtkIncrementalOctreePointLocator.h>
12#include <vtkIntArray.h>
13#include <vtkPointData.h>
14#include <vtkProbeFilter.h>
15#include <vtkRectilinearGrid.h>
16#include <vtkSmartPointer.h>
17#include <vtkTableBasedClipDataSet.h>
18#include <vtkTriangleFilter.h>
19#include <vtkUnstructuredGrid.h>
20#include <vtkXMLPolyDataWriter.h>
21#include <vtkXMLUnstructuredGridWriter.h>
22
23#include <hrleDenseIterator.hpp>
24
25#include <lsDomain.hpp>
26#include <lsMaterialMap.hpp>
28#include <unordered_map>
29#include <utility>
30
31// #define LS_TO_VISUALIZATION_DEBUG
32#ifdef LS_TO_VISUALIZATION_DEBUG
33#include <vtkXMLRectilinearGridWriter.h>
34#endif
35
36namespace viennals {
37
38using namespace viennacore;
39
46template <class T, int D> class WriteVisualizationMesh {
47 typedef typename Domain<T, D>::DomainType hrleDomainType;
48 using LevelSetsType = std::vector<SmartPointer<Domain<T, D>>>;
49 LevelSetsType levelSets;
50 SmartPointer<MaterialMap> materialMap = nullptr;
51 std::string fileName;
52 bool extractVolumeMesh = true;
53 bool extractHullMesh = false;
54 bool bottomRemoved = false;
55 bool writeToFile = true;
56 double LSEpsilon = 1e-2;
57 std::unordered_map<std::string, std::vector<double>> metaData;
58
59 vtkSmartPointer<vtkUnstructuredGrid> volumeVTK;
60 vtkSmartPointer<vtkPolyData> hullVTK;
61
65 static void removeDuplicatePoints(vtkSmartPointer<vtkPolyData> &polyData,
66 const double tolerance) {
67
68 vtkSmartPointer<vtkPolyData> newPolyData =
69 vtkSmartPointer<vtkPolyData>::New();
70 vtkSmartPointer<vtkIncrementalOctreePointLocator> ptInserter =
71 vtkSmartPointer<vtkIncrementalOctreePointLocator>::New();
72 ptInserter->SetTolerance(tolerance);
73
74 vtkSmartPointer<vtkPoints> newPoints = vtkSmartPointer<vtkPoints>::New();
75
76 // get bounds
77 double gridBounds[6];
78 polyData->GetBounds(gridBounds);
79
80 // start point insertion from original points
81 std::vector<vtkIdType> newPointIds;
82 newPointIds.reserve(polyData->GetNumberOfPoints());
83 ptInserter->InitPointInsertion(newPoints, gridBounds);
84
85 // make new point list
86 for (vtkIdType pointId = 0; pointId < polyData->GetNumberOfPoints();
87 ++pointId) {
88 vtkIdType globalPtId = 0;
89 ptInserter->InsertUniquePoint(polyData->GetPoint(pointId), globalPtId);
90 newPointIds.push_back(globalPtId);
91 }
92
93 // now set the new points to the unstructured grid
94 newPolyData->SetPoints(newPoints);
95
96 // go through all cells and change point ids to match the new ids
97 vtkSmartPointer<vtkCellArray> oldCells = polyData->GetPolys();
98 vtkSmartPointer<vtkCellArray> newCells =
99 vtkSmartPointer<vtkCellArray>::New();
100
101 vtkSmartPointer<vtkIdList> cellPoints = vtkIdList::New();
102 oldCells->InitTraversal();
103 while (oldCells->GetNextCell(cellPoints)) {
104 for (vtkIdType pointId = 0; pointId < cellPoints->GetNumberOfIds();
105 ++pointId) {
106 cellPoints->SetId(pointId, newPointIds[cellPoints->GetId(pointId)]);
107 }
108 // insert same cell with new points
109 newCells->InsertNextCell(cellPoints);
110 }
111
112 newPolyData->SetPolys(newCells);
113
114 // conserve all cell data
115 // TODO transfer point data as well (do with "InsertTuples (vtkIdList
116 // *dstIds, vtkIdList *srcIds, vtkAbstractArray *source) override)" of
117 // vtkDataArray class)
118 newPolyData->GetCellData()->ShallowCopy(polyData->GetCellData());
119
120 // set ugrid to the newly created grid
121 polyData = newPolyData;
122 }
123
126 static void removeDuplicatePoints(vtkSmartPointer<vtkUnstructuredGrid> &ugrid,
127 const double tolerance) {
128
129 vtkSmartPointer<vtkUnstructuredGrid> newGrid =
130 vtkSmartPointer<vtkUnstructuredGrid>::New();
131 vtkSmartPointer<vtkIncrementalOctreePointLocator> ptInserter =
132 vtkSmartPointer<vtkIncrementalOctreePointLocator>::New();
133 ptInserter->SetTolerance(tolerance);
134
135 vtkSmartPointer<vtkPoints> newPoints = vtkSmartPointer<vtkPoints>::New();
136
137 // get bounds
138 double gridBounds[6];
139 ugrid->GetBounds(gridBounds);
140
141 // start point insertion from original points
142 std::vector<vtkIdType> newPointIds;
143 newPointIds.reserve(ugrid->GetNumberOfPoints());
144 ptInserter->InitPointInsertion(newPoints, gridBounds);
145
146 // make new point list
147 for (vtkIdType pointId = 0; pointId < ugrid->GetNumberOfPoints();
148 ++pointId) {
149 vtkIdType globalPtId = 0;
150 ptInserter->InsertUniquePoint(ugrid->GetPoint(pointId), globalPtId);
151 newPointIds.push_back(globalPtId);
152 }
153
154 // now set the new points to the unstructured grid
155 newGrid->SetPoints(newPoints);
156
157 // go through all cells and change point ids to match the new ids
158 for (vtkIdType cellId = 0; cellId < ugrid->GetNumberOfCells(); ++cellId) {
159 vtkSmartPointer<vtkIdList> cellPoints = vtkSmartPointer<vtkIdList>::New();
160 ugrid->GetCellPoints(cellId, cellPoints);
161 for (vtkIdType pointId = 0; pointId < cellPoints->GetNumberOfIds();
162 ++pointId) {
163 cellPoints->SetId(pointId, newPointIds[cellPoints->GetId(pointId)]);
164 }
165 // insert same cell with new points
166 newGrid->InsertNextCell(ugrid->GetCell(cellId)->GetCellType(),
167 cellPoints);
168 }
169
170 // conserve all cell data
171 // TODO transfer point data as well (do with "InsertTuples (vtkIdList
172 // *dstIds, vtkIdList *srcIds, vtkAbstractArray *source) override)" of
173 // vtkDataArray class)
174 newGrid->GetCellData()->ShallowCopy(ugrid->GetCellData());
175
176 // set ugrid to the newly created grid
177 ugrid = newGrid;
178 }
179
181 static void
182 removeDegenerateTetras(vtkSmartPointer<vtkUnstructuredGrid> &ugrid) {
183 vtkSmartPointer<vtkUnstructuredGrid> newGrid =
184 vtkSmartPointer<vtkUnstructuredGrid>::New();
185
186 // need to copy material numbers
187 vtkSmartPointer<vtkIntArray> materialNumberArray =
188 vtkSmartPointer<vtkIntArray>::New();
189 materialNumberArray->SetNumberOfComponents(1);
190 materialNumberArray->SetName("Material");
191
192 // see if material is defined
193 int arrayIndex;
194 vtkDataArray *matArray =
195 ugrid->GetCellData()->GetArray("Material", arrayIndex);
196 const int &materialArrayIndex = arrayIndex;
197
198 // go through all cells and delete those with duplicate entries
199 for (vtkIdType cellId = 0; cellId < ugrid->GetNumberOfCells(); ++cellId) {
200 vtkSmartPointer<vtkIdList> cellPoints = vtkSmartPointer<vtkIdList>::New();
201 ugrid->GetCellPoints(cellId, cellPoints);
202 bool isDuplicate = false;
203 for (vtkIdType pointId = 0; pointId < cellPoints->GetNumberOfIds();
204 ++pointId) {
205 for (vtkIdType nextId = pointId + 1;
206 nextId < cellPoints->GetNumberOfIds(); ++nextId) {
207 // if they are the same, remove the cell
208 if (cellPoints->GetId(pointId) == cellPoints->GetId(nextId))
209 isDuplicate = true;
210 }
211 }
212 if (!isDuplicate) {
213 // insert same cell if no degenerate points
214 newGrid->InsertNextCell(ugrid->GetCell(cellId)->GetCellType(),
215 cellPoints);
216 // if material was defined before, use it now
217 if (materialArrayIndex >= 0)
218 materialNumberArray->InsertNextValue(matArray->GetTuple1(cellId));
219 }
220 }
221
222 // just take the old points and point data
223 newGrid->SetPoints(ugrid->GetPoints());
224 newGrid->GetPointData()->ShallowCopy(ugrid->GetPointData());
225 // set material cell data
226 newGrid->GetCellData()->SetScalars(materialNumberArray);
227
228 ugrid = newGrid;
229 }
230
231 // This function takes a levelset and converts it to a vtkRectilinearGrid
232 // The full domain contains values, which are capped at numLayers * gridDelta
233 // gridExtraPoints layers of grid points are added to the domain according to
234 // boundary conditions
235 template <int gridExtraPoints = 0>
236 vtkSmartPointer<vtkRectilinearGrid>
237 LS2RectiLinearGrid(SmartPointer<Domain<T, D>> levelSet, const double LSOffset,
238 int infiniteMinimum = std::numeric_limits<int>::max(),
239 int infiniteMaximum = -std::numeric_limits<int>::max()) {
240
241 auto &grid = levelSet->getGrid();
242 auto &domain = levelSet->getDomain();
243 double gridDelta = grid.getGridDelta();
244 int numLayers = levelSet->getLevelSetWidth();
245
246 vtkSmartPointer<vtkFloatArray>
247 coords[3]; // needs to be 3 because vtk only knows 3D
248 int gridMin = 0, gridMax = 0;
249 // int openJumpDirection = -1; // direction above the open boundary
250 // direction
251
252 // fill grid with offset depending on orientation
253 for (unsigned i = 0; i < D; ++i) {
254 coords[i] = vtkSmartPointer<vtkFloatArray>::New();
255
256 if (grid.getBoundaryConditions(i) ==
257 Domain<T, D>::BoundaryType::INFINITE_BOUNDARY) {
258 // add one to gridMin and gridMax for numerical stability
259 gridMin = std::min(domain.getMinRunBreak(i), infiniteMinimum) -
260 1; // choose the smaller number so that for first levelset the
261 // overall minimum can be chosen
262 gridMax = std::max(domain.getMaxRunBreak(i), infiniteMaximum) + 1;
263
264 // openJumpDirection = i + 1;
265 } else {
266 gridMin = grid.getMinGridPoint(i) - gridExtraPoints;
267 gridMax = grid.getMaxGridPoint(i) + gridExtraPoints;
268 }
269
270 for (int x = gridMin; x <= gridMax; ++x) {
271 coords[i]->InsertNextValue(x * gridDelta);
272 }
273 }
274
275 // if we work in 2D, just add 1 grid point at origin
276 if (D == 2) {
277 coords[2] = vtkSmartPointer<vtkFloatArray>::New();
278 coords[2]->InsertNextValue(0);
279 }
280
281 vtkSmartPointer<vtkRectilinearGrid> rgrid =
282 vtkSmartPointer<vtkRectilinearGrid>::New();
283
284 rgrid->SetDimensions(coords[0]->GetNumberOfTuples(),
285 coords[1]->GetNumberOfTuples(),
286 coords[2]->GetNumberOfTuples());
287 rgrid->SetXCoordinates(coords[0]);
288 rgrid->SetYCoordinates(coords[1]);
289 rgrid->SetZCoordinates(coords[2]);
290
291 // bool fixBorderPoints = (gridExtraPoints != 0);
292
293 // need to save the current position one dimension above open boundary
294 // direction, so we can register a jump in the open boundary direction when
295 // it occurs, so we can fix the LS value as follows: if remove_bottom==true,
296 // we need to flip the sign otherwise the sign stays the same
297 // int currentOpenIndex =
298 // it.getIndices()[(openJumpDirection < D) ? openJumpDirection : 0];
299
300 // Make array to store signed distance function
301 auto const numGridPoints = rgrid->GetNumberOfPoints();
302 vtkSmartPointer<vtkFloatArray> signedDistances =
303 vtkSmartPointer<vtkFloatArray>::New();
304 signedDistances->SetNumberOfComponents(1);
305 signedDistances->SetNumberOfTuples(numGridPoints);
306 signedDistances->SetName("SignedDistances");
307
308#pragma omp parallel
309 {
310 // use dense iterator to got to every index location
311 viennahrle::ConstDenseIterator<typename Domain<T, D>::DomainType> it(
312 levelSet->getDomain());
313
314#pragma omp for
315 for (vtkIdType pointId = 0; pointId < numGridPoints; ++pointId) {
316 // iterate until all grid points have a signed distance value
317
318 double p[3];
319 rgrid->GetPoint(pointId, p);
320 // create index vector
321 viennahrle::Index<D> indices(grid.globalCoordinates2GlobalIndices(p));
322
323 // write the corresponding LSValue
324 T value;
325
326 // if indices are outside of domain map to local indices
327 if (grid.isOutsideOfDomain(indices)) {
328 indices = grid.globalIndices2LocalIndices(indices);
329 // fixBorderPoints = true;
330 // signedDistances->InsertNextValue(
331 // signedDistances->GetDataTypeValueMax());
332 }
333
334 it.goToIndices(indices);
335 if (it.getValue() == Domain<T, D>::POS_VALUE) {
336 value = numLayers;
337 } else if (it.getValue() == Domain<T, D>::NEG_VALUE) {
338 value = -numLayers;
339 } else {
340 value = it.getValue() + LSOffset;
341 }
342
343 // if (removeBottom) {
344 // // if we jump from one end of the domain to the other and are not
345 // // already in the new run, we need to fix the sign of the run
346 // if (currentOpenIndex != indices[openJumpDirection]) {
347 // currentOpenIndex = indices[openJumpDirection];
348 // if (indices >= it.getIndices()) {
349 // value = -value;
350 // }
351 // }
352 // }
353
354 signedDistances->SetValue(pointId, value * gridDelta);
355 // signedDistances->InsertNextValue(value * gridDelta);
356 }
357
358 // // advance iterator until it is at correct point
359 // while (compare(it_l.end_indices(), indices) < 0) {
360 // it_l.next();
361 // if (it_l.is_finished())
362 // break;
363 // }
364 // // now move iterator with pointId to get to next point
365 // switch (compare(it_l.end_indices(), indices)) {
366 // case 0:
367 // it_l.next();
368 // default:
369 // ++pointId;
370 // }
371 // }
372 }
373
374 // now need to go through again to fix border points, this is done by
375 // mapping existing points onto the points outside the domain according
376 // to the correct boundary conditions
377 // if (fixBorderPoints) {
378 // vtkIdType pointId = 0;
379 // while ((pointId < rgrid->GetNumberOfPoints())) {
380 // if (signedDistances->GetValue(pointId) ==
381 // signedDistances->GetDataTypeValueMax()) {
382 // double p[3];
383 // rgrid->GetPoint(pointId, p);
384
385 // // create index vector
386 // viennahrle::Index<D> indices(
387 // grid.globalCoordinates2GlobalIndices(p));
388
389 // // vector for mapped point inside domain
390 // viennahrle::Index<D> localIndices =
391 // grid.globalIndices2LocalIndices(indices);
392
393 // // now find ID of point we need to take value from
394 // int originalPointId = 0;
395 // for (int i = D - 1; i >= 0; --i) {
396 // originalPointId *=
397 // coords[i]->GetNumberOfTuples(); // extent in direction
398 // originalPointId += localIndices[i] - indices[i];
399 // }
400 // originalPointId += pointId;
401
402 // // now put value of mapped point in global point
403 // signedDistances->SetValue(pointId,
404 // signedDistances->GetValue(originalPointId));
405 // }
406 // ++pointId;
407 // }
408 // }
409
410 // Add the SignedDistances to the grid
411 rgrid->GetPointData()->SetScalars(signedDistances);
412
413 return rgrid;
414 }
415
416 void addMetaDataToVTK(vtkDataSet *data) const {
417 if (metaData.empty()) {
418 return;
419 }
420
421 // add metadata to field data
422 vtkSmartPointer<vtkFieldData> fieldData = data->GetFieldData();
423 for (const auto &meta : metaData) {
424 if (meta.second.empty())
425 continue; // skip empty metadata
426
427 vtkSmartPointer<vtkFloatArray> metaDataArray =
428 vtkSmartPointer<vtkFloatArray>::New();
429 metaDataArray->SetName(meta.first.c_str());
430 metaDataArray->SetNumberOfValues(meta.second.size());
431 for (size_t i = 0; i < meta.second.size(); ++i) {
432 metaDataArray->SetValue(i, meta.second[i]);
433 }
434 fieldData->AddArray(metaDataArray);
435 }
436 }
437
438public:
439 WriteVisualizationMesh() = default;
440
441 WriteVisualizationMesh(SmartPointer<Domain<T, D>> levelSet) {
442 levelSets.push_back(levelSet);
443 }
444
446 void insertNextLevelSet(SmartPointer<Domain<T, D>> levelSet) {
447 levelSets.push_back(levelSet);
448 }
449
450 void clearLevelSets() { levelSets.clear(); }
451
454 void setFileName(std::string passedFileName) {
455 fileName = std::move(passedFileName);
456 }
457
459 void setExtractHullMesh(bool passedExtractHullMesh) {
460 extractHullMesh = passedExtractHullMesh;
461 }
462
464 void setExtractVolumeMesh(bool passedExtractVolumeMesh) {
465 extractVolumeMesh = passedExtractVolumeMesh;
466 }
467
468 void setWriteToFile(bool passedWriteToFile) {
469 writeToFile = passedWriteToFile;
470 }
471
472 void setMaterialMap(SmartPointer<MaterialMap> passedMaterialMap) {
473 materialMap = passedMaterialMap;
474 }
475
476 void setWrappingLayerEpsilon(double epsilon) { LSEpsilon = epsilon; }
477
478 void setMetaData(const std::unordered_map<std::string, std::vector<double>>
479 &passedMetaData) {
480 metaData = passedMetaData;
481 }
482
483 void addMetaData(const std::string &key, double value) {
484 metaData[key] = std::vector<double>{value};
485 }
486
487 void addMetaData(const std::string &key, const std::vector<double> &values) {
488 metaData[key] = values;
489 }
490
491 void addMetaData(
492 const std::unordered_map<std::string, std::vector<double>> &newMetaData) {
493 for (const auto &pair : newMetaData) {
494 metaData[pair.first] = pair.second;
495 }
496 }
497
498 auto getVolumeMesh() const { return volumeVTK; }
499
500 auto getHullMesh() const { return hullVTK; }
501
502 void apply() {
503 // check if level sets have enough layers
504 for (unsigned i = 0; i < levelSets.size(); ++i) {
505 if (levelSets[i]->getLevelSetWidth() < 2) {
506 VIENNACORE_LOG_WARNING(
507 "WriteVisualizationMesh: Level Set " + std::to_string(i) +
508 " should have a width greater than 1! Conversion might fail!");
509 }
510 }
511
512 const double gridDelta = levelSets[0]->getGrid().getGridDelta();
513
514 // store volume for each material
515 std::vector<vtkSmartPointer<vtkUnstructuredGrid>> materialMeshes;
516 std::vector<unsigned> materialIds;
517
518 int totalMinimum = std::numeric_limits<int>::max();
519 int totalMaximum = -std::numeric_limits<int>::max();
520 for (auto &it : levelSets) {
521 if (it->getNumberOfPoints() == 0) {
522 continue;
523 }
524 auto &grid = it->getGrid();
525 auto &domain = it->getDomain();
526 for (unsigned i = 0; i < D; ++i) {
527 if (grid.getBoundaryConditions(i) ==
528 Domain<T, D>::BoundaryType::INFINITE_BOUNDARY) {
529 totalMinimum = std::min(totalMinimum, domain.getMinRunBreak(i));
530 totalMaximum = std::max(totalMaximum, domain.getMaxRunBreak(i));
531 }
532 }
533 }
534
535 // create volume mesh for largest LS
536 // Use vtkClipDataSet to slice the grid
537 vtkSmartPointer<vtkTableBasedClipDataSet> clipper =
538 vtkSmartPointer<vtkTableBasedClipDataSet>::New();
539 auto topGrid = vtkSmartPointer<vtkRectilinearGrid>::New();
540 // if (bottomRemoved) {
541 topGrid =
542 LS2RectiLinearGrid(levelSets.back(), 0, totalMinimum, totalMaximum);
543 // } else {
544 // topGrid = LS2RectiLinearGrid<false>(levelSets.back(), 0,
545 // totalMinimum,
546 // totalMaximum);
547 // }
548#ifdef LS_TO_VISUALIZATION_DEBUG
549 {
550 auto gwriter = vtkSmartPointer<vtkXMLRectilinearGridWriter>::New();
551 gwriter->SetFileName("./grid_0.vtr");
552 gwriter->SetInputData(topGrid);
553 gwriter->Write();
554 std::cout << "Wrote grid 0" << std::endl;
555 }
556#endif
557 clipper->SetInputData(topGrid);
558 clipper->InsideOutOn();
559 clipper->SetValue(0.0);
560 clipper->GenerateClippedOutputOn(); // TODO remove
561 clipper->Update();
562
563#ifdef LS_TO_VISUALIZATION_DEBUG
564 {
565 auto gwriter = vtkSmartPointer<vtkXMLUnstructuredGridWriter>::New();
566 gwriter->SetFileName("./clipped.vtu");
567 gwriter->SetInputData(clipper->GetClippedOutput());
568 gwriter->Write();
569 std::cout << "Wrote clipped" << std::endl;
570 }
571#endif
572
573 const bool useMaterialMap = materialMap != nullptr;
574 materialMeshes.emplace_back(clipper->GetOutput());
575 materialIds.push_back(useMaterialMap ? materialMap->getMaterialId(0) : 0);
576
577#ifdef LS_TO_VISUALIZATION_DEBUG
578 {
579 auto gwriter = vtkSmartPointer<vtkXMLUnstructuredGridWriter>::New();
580 gwriter->SetFileName("./probed_0.vtu");
581 gwriter->SetInputData(materialMeshes.front());
582 gwriter->Write();
583 }
584#endif
585
586 unsigned counter = 1;
587
588 // now cut large volume mesh with all the smaller ones
589 for (typename LevelSetsType::const_reverse_iterator it =
590 ++levelSets.rbegin();
591 it != levelSets.rend(); ++it) {
592 if (it->get()->getNumberOfPoints() == 0)
593 continue; // ignore empty levelSets
594
595 // create grid of next LS with slight offset and project into current
596 // mesh
597 vtkSmartPointer<vtkRectilinearGrid> rgrid =
598 vtkSmartPointer<vtkRectilinearGrid>::New();
599 // if (bottomRemoved) {
600 rgrid = LS2RectiLinearGrid<1>(*it, -LSEpsilon * counter, totalMinimum,
601 totalMaximum);
602 // } else {
603 // rgrid = LS2RectiLinearGrid<false, 1>(*it, -LSEpsilon * counter,
604 // totalMinimum, totalMaximum);
605 // }
606
607#ifdef LS_TO_VISUALIZATION_DEBUG
608 {
609 vtkSmartPointer<vtkXMLRectilinearGridWriter> gwriter =
610 vtkSmartPointer<vtkXMLRectilinearGridWriter>::New();
611 gwriter->SetFileName(
612 ("./grid_" + std::to_string(counter) + ".vtr").c_str());
613 gwriter->SetInputData(rgrid);
614 gwriter->Write();
615 std::cout << "Wrote grid " << counter << std::endl;
616 }
617#endif
618
619 // now transfer implicit values to mesh points
620 vtkSmartPointer<vtkProbeFilter> probeFilter =
621 vtkSmartPointer<vtkProbeFilter>::New();
622 probeFilter->SetInputData(materialMeshes.back()); // last element
623 probeFilter->SetSourceData(rgrid);
624 probeFilter->Update();
625
626#ifdef LS_TO_VISUALIZATION_DEBUG
627 {
628 vtkSmartPointer<vtkXMLUnstructuredGridWriter> gwriter =
629 vtkSmartPointer<vtkXMLUnstructuredGridWriter>::New();
630 gwriter->SetFileName(
631 ("./probed_" + std::to_string(counter) + ".vtu").c_str());
632 gwriter->SetInputData(probeFilter->GetOutput());
633 gwriter->Write();
634 std::cout << "Wrote unstructured grid " << counter << std::endl;
635 }
636#endif
637
638 // now clip the mesh and save the clipped as the 1st layer and use the
639 // inverse for the next layer clipping Use vtkTabelBasedClipDataSet to
640 // slice the grid
641 vtkSmartPointer<vtkTableBasedClipDataSet> insideClipper =
642 vtkSmartPointer<vtkTableBasedClipDataSet>::New();
643 insideClipper->SetInputConnection(probeFilter->GetOutputPort());
644 insideClipper->GenerateClippedOutputOn();
645 insideClipper->Update();
646
647 materialMeshes.back() = insideClipper->GetOutput();
648 materialMeshes.emplace_back(insideClipper->GetClippedOutput());
649 auto material = counter;
650 if (useMaterialMap)
651 material = materialMap->getMaterialId(counter);
652 materialIds.push_back(material);
653
654 ++counter;
655 }
656
657 vtkSmartPointer<vtkAppendFilter> appendFilter =
658 vtkSmartPointer<vtkAppendFilter>::New();
659
660 vtkSmartPointer<vtkAppendPolyData> hullAppendFilter =
661 vtkSmartPointer<vtkAppendPolyData>::New();
662
663 for (unsigned i = 0; i < materialMeshes.size(); ++i) {
664
665 // write material number in mesh
666 vtkSmartPointer<vtkIntArray> materialNumberArray =
667 vtkSmartPointer<vtkIntArray>::New();
668 materialNumberArray->SetNumberOfComponents(1);
669 materialNumberArray->SetName("Material");
670 for (unsigned j = 0;
671 j <
672 materialMeshes[materialMeshes.size() - 1 - i]->GetNumberOfCells();
673 ++j) {
674 materialNumberArray->InsertNextValue(materialIds[i]);
675 }
676 materialMeshes[materialMeshes.size() - 1 - i]->GetCellData()->SetScalars(
677 materialNumberArray);
678
679 // delete all point data, so it is not in ouput
680 // TODO this includes signed distance information which could be
681 // conserved for debugging also includes wheter a cell was vaild for
682 // cutting by the grid
683 vtkSmartPointer<vtkPointData> pointData =
684 materialMeshes[materialMeshes.size() - 1 - i]->GetPointData();
685 const int numberOfArrays = pointData->GetNumberOfArrays();
686 for (int j = 0; j < numberOfArrays; ++j) {
687 pointData->RemoveArray(0); // remove first array until none are left
688 }
689
690 // if hull mesh should be exported, create hull for each layer and put
691 // them together
692 if (extractHullMesh) {
693 vtkSmartPointer<vtkGeometryFilter> geoFilter =
694 vtkSmartPointer<vtkGeometryFilter>::New();
695 geoFilter->SetInputData(materialMeshes[materialMeshes.size() - 1 - i]);
696 geoFilter->Update();
697 hullAppendFilter->AddInputData(geoFilter->GetOutput());
698 }
699
700 appendFilter->AddInputData(materialMeshes[materialMeshes.size() - 1 - i]);
701 }
702
703 // do not need tetrahedral volume mesh if we do not print volume
704 volumeVTK = vtkSmartPointer<vtkUnstructuredGrid>::New();
705 hullVTK = vtkSmartPointer<vtkPolyData>::New();
706 if (extractVolumeMesh) {
707 appendFilter->Update();
708
709 // remove degenerate points and remove cells which collapse to zero
710 // volume then
711 volumeVTK = appendFilter->GetOutput();
712#ifdef LS_TO_VISUALIZATION_DEBUG
713 {
714 std::cout << "Before duplicate removal: " << std::endl;
715 std::cout << "Points: " << volumeVTK->GetNumberOfPoints() << std::endl;
716 std::cout << "Cells: " << volumeVTK->GetNumberOfCells() << std::endl;
717 vtkSmartPointer<vtkXMLUnstructuredGridWriter> gwriter =
718 vtkSmartPointer<vtkXMLUnstructuredGridWriter>::New();
719 gwriter->SetFileName("before_removal.vtu");
720 gwriter->SetInputData(appendFilter->GetOutput());
721 gwriter->Update();
722 }
723#endif
724
725 // use 1/1000th of grid spacing for contraction of two similar points,
726 // so that tetrahedralisation works correctly
727 removeDuplicatePoints(volumeVTK, 1e-3 * gridDelta);
728
729#ifdef LS_TO_VISUALIZATION_DEBUG
730 {
731 std::cout << "After duplicate removal: " << std::endl;
732 std::cout << "Points: " << volumeVTK->GetNumberOfPoints() << std::endl;
733 std::cout << "Cells: " << volumeVTK->GetNumberOfCells() << std::endl;
734 vtkSmartPointer<vtkXMLUnstructuredGridWriter> gwriter =
735 vtkSmartPointer<vtkXMLUnstructuredGridWriter>::New();
736 gwriter->SetFileName("after_removal.vtu");
737 gwriter->SetInputData(volumeVTK);
738 gwriter->Update();
739 }
740#endif
741
742 // change all 3D cells into tetras and all 2D cells to triangles
743 vtkSmartPointer<vtkDataSetTriangleFilter> triangleFilter =
744 vtkSmartPointer<vtkDataSetTriangleFilter>::New();
745 triangleFilter->SetInputData(volumeVTK);
746 triangleFilter->Update();
747 volumeVTK = triangleFilter->GetOutput();
748
749 // now that only tetras are left, remove the ones with degenerate points
750 removeDegenerateTetras(volumeVTK);
751
752 // add meta data to volume mesh
753 addMetaDataToVTK(volumeVTK);
754
755 if (writeToFile) {
756 auto writer = vtkSmartPointer<vtkXMLUnstructuredGridWriter>::New();
757 writer->SetFileName((fileName + "_volume.vtu").c_str());
758 writer->SetInputData(volumeVTK);
759 writer->Write();
760 }
761 }
762
763 // Now make hull mesh if necessary
764 if (extractHullMesh) {
765 hullAppendFilter->Update();
766 hullVTK = hullAppendFilter->GetOutput();
767 // use 1/1000th of grid spacing for contraction of two similar points
768 removeDuplicatePoints(hullVTK, 1e-3 * gridDelta);
769
770 vtkSmartPointer<vtkTriangleFilter> hullTriangleFilter =
771 vtkSmartPointer<vtkTriangleFilter>::New();
772 hullTriangleFilter->SetInputData(hullVTK);
773 hullTriangleFilter->Update();
774
775 hullVTK = hullTriangleFilter->GetOutput();
776
777 // add meta data to hull mesh
778 addMetaDataToVTK(hullVTK);
779
780 if (writeToFile) {
781 auto writer = vtkSmartPointer<vtkXMLPolyDataWriter>::New();
782 writer->SetFileName((fileName + "_hull.vtp").c_str());
783 writer->SetInputData(hullVTK);
784 writer->Write();
785 } else {
786 hullTriangleFilter->Update();
787 }
788 }
789 }
790};
791
792// add all template specialisations for this class
793PRECOMPILE_PRECISION_DIMENSION(WriteVisualizationMesh)
794
795} // namespace viennals
796
797#endif // VIENNALS_USE_VTK
constexpr int D
Definition Epitaxy.cpp:11
double T
Definition Epitaxy.cpp:12
#define PRECOMPILE_PRECISION_DIMENSION(className)
Definition lsPreCompileMacros.hpp:24
float gridDelta
Definition AirGapDeposition.py:21
writer
Definition AirGapDeposition.py:89
int counter
Definition Deposition.py:75
Definition lsAdvect.hpp:37