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