TTK
Loading...
Searching...
No Matches
ReebSpace.h
Go to the documentation of this file.
1
28
29#pragma once
30
31// base code includes
32#include <FiberSurface.h>
33#include <Geometry.h>
34#include <JacobiSet.h>
35#include <Triangulation.h>
36
37#include <map>
38#include <set>
39
40namespace ttk {
41
42 class ReebSpace : virtual public Debug {
43 public:
45 domainVolume, // 0
46 rangeArea, // 1
47 hyperVolume // 2
48 };
49
50 struct Sheet0 {
52 // 0: extrema-sheet, 1: saddle-sheet
53 char type_{};
54 std::vector<SimplexId> sheet1List_{};
55 std::vector<SimplexId> sheet3List_{};
56 };
57
58 struct Sheet1 {
60 std::vector<SimplexId> edgeList_{};
61 // 0: extrema-sheet, 1: saddle-sheet (can be inferred by the number
62 // of 2-sheets attached to it?)
63 std::vector<SimplexId> sheet0List_{};
64 // NB: the corresponding 2-sheet should have the same global
65 // identifier.
66 std::vector<SimplexId> sheet3List_{};
67 };
68
69 struct Sheet2 {
70 bool pruned_{};
72 // for each point, x, y, z coordinates
73 // this is how coincident point will be merged afterwards
74 // this list is meant to be temporary. after the merge, we'll use a
75 // global list
76 std::vector<std::vector<FiberSurface::Vertex>> vertexList_{};
77 std::vector<std::vector<FiberSurface::Triangle>> triangleList_{};
78 std::vector<SimplexId> sheet3List_{};
79 };
80
81 struct Sheet3 {
83 bool pruned_{};
85 std::vector<SimplexId> vertexList_{};
86 std::vector<SimplexId> tetList_{};
87 std::vector<SimplexId> sheet0List_{};
88 std::vector<SimplexId> sheet1List_{};
89 std::vector<SimplexId> sheet2List_{};
90 std::vector<SimplexId> sheet3List_{};
91 std::vector<SimplexId> preMergedSheets_{};
92 };
93
94 ReebSpace();
95
96 inline bool empty() const {
97 return (currentData_.vertex2sheet0_.size() == 0);
98 }
99
100 template <class dataTypeU, class dataTypeV, typename triangulationType>
101 inline int execute(const dataTypeU *const uField,
102 const dataTypeV *const vField,
103 const triangulationType &triangulation);
104
105 inline const Sheet0 *get0sheet(const SimplexId &sheetId) const {
106
107#ifndef TTK_ENABLE_KAMIKAZE
108 if((sheetId < 0)
109 || (sheetId >= (SimplexId)currentData_.sheet0List_.size()))
110 return nullptr;
111#endif
112
113 return &(currentData_.sheet0List_[sheetId]);
114 }
115
116 inline const Sheet1 *get1sheet(const SimplexId &sheetId) const {
117#ifndef TTK_ENABLE_KAMIKAZE
118 if((sheetId < 0)
119 || (sheetId >= (SimplexId)currentData_.sheet1List_.size()))
120 return nullptr;
121#endif
122
123 return &(currentData_.sheet1List_[sheetId]);
124 }
125
126 // warning, these are the originals
127 inline const Sheet2 *get2sheet(const SimplexId &sheetId) const {
128#ifndef TTK_ENABLE_KAMIKAZE
129 if((sheetId < 0)
130 || (sheetId >= (SimplexId)originalData_.sheet2List_.size()))
131 return nullptr;
132#endif
133
134 return &(originalData_.sheet2List_[sheetId]);
135 }
136
137 inline const Sheet3 *get3sheet(const SimplexId &sheetId) const {
138#ifndef TTK_ENABLE_KAMIKAZE
139 if((sheetId < 0)
140 || (sheetId >= (SimplexId)currentData_.sheet3List_.size()))
141 return nullptr;
142#endif
143
144 return &(currentData_.sheet3List_[sheetId]);
145 }
146
147 inline const std::vector<SimplexId> *get0sheetSegmentation() const {
148 return &currentData_.vertex2sheet0_;
149 }
150
151 const std::vector<SimplexId> *get1sheetSegmentation() const {
152 return &currentData_.edge2sheet1_;
153 }
154
155 const std::vector<SimplexId> *get3sheetVertexSegmentation() const {
156 return &currentData_.vertex2sheet3_;
157 }
158
159 const std::vector<SimplexId> *get3sheetTetSegmentation() const {
160 return &currentData_.tet2sheet3_;
161 }
162
163 const std::vector<char> *getEdgeTypes() const {
164 return &currentData_.edgeTypes_;
165 }
166
167 inline const std::vector<FiberSurface::Vertex> *
169 return &(fiberSurfaceVertexList_);
170 }
171
172 inline int getJacobi2Edge(const SimplexId &jacobiEdgeId) const {
173 if((jacobiEdgeId < 0)
174 || (jacobiEdgeId >= (SimplexId)jacobi2edges_.size()))
175 return -1;
176 return jacobi2edges_[jacobiEdgeId];
177 }
178
180 return currentData_.sheet2List_.size();
181 }
182
183 // inline std::vector<long long int>* getSheetTriangulationCells(){
184 // return &sheet3cells_;
185 // }
186 //
187 // inline std::vector<float>* getSheetTriangulationPoints(){
188 // return &sheet3points_;
189 // }
190
191 template <class dataTypeU, class dataTypeV>
192 inline int
193 perturb(const dataTypeU *const uField,
194 const dataTypeV *const vField,
195 const dataTypeU &uEpsilon = Geometry::powIntTen(-DBL_DIG),
196 const dataTypeV &vEpsilon = Geometry::powIntTen(-DBL_DIG));
197
198 inline void setExpand3Sheets(const bool &onOff) {
199 expand3sheets_ = onOff;
200 }
201
202 inline bool setRangeDrivenOctree(const bool &onOff) {
203 if(onOff != withRangeDrivenOctree_) {
205 return true;
206 }
207
209 return false;
210 }
211
220 inline void setSosOffsetsU(const SimplexId *const sosOffsetsU) {
221 sosOffsetsU_ = sosOffsetsU;
222 }
223
232 inline void setSosOffsetsV(const SimplexId *const sosOffsetsV) {
233 sosOffsetsV_ = sosOffsetsV;
234 }
235
236 inline void setTetNumber(const SimplexId &tetNumber) {
237 tetNumber_ = tetNumber;
238 }
239
243 inline void setVertexNumber(const SimplexId &vertexNumber) {
244 vertexNumber_ = vertexNumber;
245 }
246
247 // WARNING: if you plan to use the range driven octree, make sure
248 // that you provided pointers to the u and v fields.
249 inline int
251 if(triangulation) {
252 triangulation->preconditionVertexStars();
253 triangulation->preconditionEdges();
254 triangulation->preconditionVertexEdges();
255
256 // trigger the jacobiSet pre-processing on the triangulation.
258
259 // trigger the fiberSurface pre-processing on the triangulation.
261
262 vertexNumber_ = triangulation->getNumberOfVertices();
263 edgeNumber_ = triangulation->getNumberOfEdges();
264 tetNumber_ = triangulation->getNumberOfCells();
265 }
266
267 return 0;
268 }
269
270 template <class dataTypeU, class dataTypeV, typename triangulationType>
271 inline int simplify(const dataTypeU *const uField,
272 const dataTypeV *const vField,
273 const triangulationType &triangulation,
274 const double &simplificationThreshold,
275 const SimplificationCriterion &criterion
277
278 private:
279 // output segmentation
280 struct ReebSpaceData {
281 SimplificationCriterion simplificationCriterion_{};
282 double simplificationThreshold_{};
283
284 std::vector<SimplexId> edge2sheet1_{};
285 std::vector<char> edgeTypes_{};
286 std::vector<SimplexId> tet2sheet3_{};
287 std::vector<SimplexId> vertex2sheet0_{};
288 std::vector<SimplexId> vertex2sheet3_{};
289
290 // structure
291 std::vector<Sheet0> sheet0List_{};
292 std::vector<Sheet1> sheet1List_{};
293 std::vector<Sheet2> sheet2List_{};
294 std::vector<Sheet3> sheet3List_{};
295
296 // std::vector<float> sheet3points_;
297 // std::vector<long long int> sheet3cells_;
298 };
299
300 protected:
301 template <typename triangulationType>
303 const std::vector<std::pair<SimplexId, char>> &jacobiSet,
304 std::vector<std::pair<SimplexId, SimplexId>> &jacobiSetClassification,
305 const triangulationType &triangulation);
306
307 template <typename triangulationType>
308 int compute1sheets(
309 const std::vector<std::pair<SimplexId, char>> &jacobiSet,
310 std::vector<std::pair<SimplexId, SimplexId>> &jacobiSetClassification,
311 const triangulationType &triangulation);
312
313 template <class dataTypeU, class dataTypeV, typename triangulationType>
314 inline int compute2sheets(
315 const std::vector<std::pair<SimplexId, SimplexId>> &jacobiEdges,
316 const dataTypeU *const uField,
317 const dataTypeV *const vField,
318 const triangulationType &triangulation);
319
320 template <class dataTypeU, class dataTypeV, typename triangulationType>
321 inline int compute2sheetChambers(const dataTypeU *const uField,
322 const dataTypeV *const vField,
323 const triangulationType &triangulation);
324
325 template <typename triangulationType>
326 int compute3sheet(
327 const SimplexId &vertexId,
328 const std::vector<std::vector<std::array<SimplexId, 3>>> &tetTriangles,
329 const triangulationType &triangulation);
330
331 template <typename triangulationType>
332 int compute3sheets(
333 std::vector<std::vector<std::array<SimplexId, 3>>> &tetTriangles,
334 const triangulationType &triangulation);
335
336 template <class dataTypeU, class dataTypeV, typename triangulationType>
337 inline int
338 computeGeometricalMeasures(Sheet3 &sheet,
339 const dataTypeU *const uField,
340 const dataTypeV *const vField,
341 const triangulationType &triangulation);
342
343 int connect3sheetTo0sheet(ReebSpaceData &data,
344 const SimplexId &sheet3Id,
345 const SimplexId &sheet0Id);
346
347 int connect3sheetTo1sheet(ReebSpaceData &data,
348 const SimplexId &sheet3Id,
349 const SimplexId &sheet1Id);
350
351 int connect3sheetTo2sheet(ReebSpaceData &data,
352 const SimplexId &sheet3Id,
353 const SimplexId &sheet2Id);
354
355 int connect3sheetTo3sheet(ReebSpaceData &data,
356 const SimplexId &sheet3Id,
357 const SimplexId &otherSheet3Id);
358
359 template <typename triangulationType>
360 int connectSheets(const triangulationType &triangulation);
361
362 int disconnect1sheetFrom0sheet(ReebSpaceData &data,
363 const SimplexId &sheet1Id,
364 const SimplexId &sheet0Id,
365 const SimplexId &biggerId);
366
367 int disconnect3sheetFrom0sheet(ReebSpaceData &data,
368 const SimplexId &sheet3Id,
369 const SimplexId &sheet0Id);
370
371 template <typename triangulationType>
372 int disconnect3sheetFrom1sheet(ReebSpaceData &data,
373 const SimplexId &sheet3Id,
374 const SimplexId &sheet1Id,
375 const SimplexId &biggerId,
376 const triangulationType &triangulation);
377
378 int disconnect3sheetFrom2sheet(ReebSpaceData &data,
379 const SimplexId &sheet3Id,
380 const SimplexId &sheet2Id);
381
382 int disconnect3sheetFrom3sheet(ReebSpaceData &data,
383 const SimplexId &sheet3Id,
384 const SimplexId &other3SheetId);
385
386 template <typename triangulationType>
387 int flush(const triangulationType &triangulation);
388
389 template <typename triangulationType>
390 int mergeSheets(const SimplexId &smallerId,
391 const SimplexId &biggerId,
392 const triangulationType &triangulation);
393
394 int preMergeSheets(const SimplexId &sheetId0, const SimplexId &sheetId1);
395
397
398 int printConnectivity(const ReebSpaceData &data) const;
399
400 template <typename triangulationType>
401 int simplifySheets(const double &simplificationThreshold,
402 const SimplificationCriterion &simplificationCriterion,
403 const triangulationType &triangulation);
404
405 template <typename triangulationType>
406 int simplifySheet(const SimplexId &sheetId,
407 const SimplificationCriterion &simplificationCriterion,
408 const triangulationType &triangulation);
409
410 // int triangulateTetrahedron(const int &tetId,
411 // const std::vector<std::vector<int> > &triangles,
412 // std::vector<long long int> &outputTets);
413 //
414 // int triangulateThreeSheets();
415
418
420
423 ReebSpaceData originalData_{}, currentData_{};
424
425 // information that does not get simplified
426 std::vector<std::pair<SimplexId, char>> jacobiSetEdges_{};
427 std::vector<SimplexId> jacobi2edges_{};
428
431 std::vector<FiberSurface::Vertex> fiberSurfaceVertexList_{};
432 };
433} // namespace ttk
434
435// if the package is not a template, comment the following line
436// #include <ReebSpace.cpp>
437
438template <class dataTypeU, class dataTypeV, typename triangulationType>
439inline int ttk::ReebSpace::execute(const dataTypeU *const uField,
440 const dataTypeV *const vField,
441 const triangulationType &triangulation) {
442
443 flush(triangulation);
444 fiberSurface_.setInputField(uField, vField);
445
446#ifdef TTK_ENABLE_FIBER_SURFACE_WITH_RANGE_OCTREE
447 fiberSurface_.flushOctree();
449 fiberSurface_.buildOctree<dataTypeU, dataTypeV>(&triangulation);
450 }
451#endif
452
453 Timer t;
454
455 // 1) compute the jacobi set
458 jacobiSet_.execute(jacobiSetEdges_, uField, vField, triangulation);
459
460 // 2) compute the list saddle 1-sheets
461 // + list of saddle 0-sheets
462 std::vector<std::pair<SimplexId, SimplexId>> jacobiSetClassification;
463 compute1sheetsOnly(jacobiSetEdges_, jacobiSetClassification, triangulation);
464 // at this stage, jacobiSetClassification contains the list of saddle edges
465 // along with their 1-sheet Id.
466
467 compute2sheets(jacobiSetClassification, uField, vField, triangulation);
468 // compute2sheetChambers<dataTypeU, dataTypeV>();
469
470 std::vector<std::vector<std::array<SimplexId, 3>>> tetTriangles;
471 compute3sheets(tetTriangles, triangulation);
472
473 this->printMsg(
474 "Data-set processed", 1.0, t.getElapsedTime(), this->threadNumber_);
475
476 // post-process for further interaction
477 if((totalArea_ == -1) || (totalVolume_ == -1) || (totalHyperVolume_ == -1)) {
478
479 Timer tm;
480
481#ifdef TTK_ENABLE_OPENMP
482#pragma omp parallel for num_threads(threadNumber_)
483#endif
484 for(size_t i = 0; i < originalData_.sheet3List_.size(); i++) {
486 originalData_.sheet3List_[i], uField, vField, triangulation);
487 }
488
489 for(size_t i = 0; i < originalData_.sheet3List_.size(); i++) {
490 totalArea_ += originalData_.sheet3List_[i].rangeArea_;
491 totalVolume_ += originalData_.sheet3List_[i].domainVolume_;
492 totalHyperVolume_ += originalData_.sheet3List_[i].hyperVolume_;
493 }
494
495 this->printMsg("Computed geometrical measures", 1.0, tm.getElapsedTime(),
496 this->threadNumber_);
497 }
498
499 fiberSurface_.finalize<dataTypeU, dataTypeV>();
500
502
503 return 0;
504}
505
506template <class dataTypeU, class dataTypeV, typename triangulationType>
508 const std::vector<std::pair<SimplexId, SimplexId>> &jacobiEdges,
509 const dataTypeU *const uField,
510 const dataTypeV *const vField,
511 const triangulationType &triangulation) {
512
513 Timer t;
514
515 // at this point, they have exactly the same size
516 originalData_.sheet2List_.resize(originalData_.sheet1List_.size());
517 for(size_t i = 0; i < originalData_.sheet2List_.size(); i++) {
518 originalData_.sheet2List_[i].sheet1Id_ = i;
519 originalData_.sheet2List_[i].pruned_ = false;
520 originalData_.sheet2List_[i].vertexList_.resize(
521 originalData_.sheet1List_[originalData_.sheet2List_[i].sheet1Id_]
522 .edgeList_.size());
523 originalData_.sheet2List_[i].triangleList_.resize(
524 originalData_.sheet1List_[originalData_.sheet2List_[i].sheet1Id_]
525 .edgeList_.size());
526 for(size_t j = 0; j < originalData_.sheet2List_[i].vertexList_.size();
527 j++) {
528 originalData_.sheet2List_[i].vertexList_[j].clear();
529 originalData_.sheet2List_[i].triangleList_[j].clear();
530 }
531 }
532
533 fiberSurface_.setGlobalVertexList(&fiberSurfaceVertexList_);
534 fiberSurface_.setPolygonEdgeNumber(jacobiEdges.size());
535
536 std::vector<SimplexId> edge2polygonEdgeId(edgeNumber_, -1);
537 jacobi2edges_.resize(jacobiEdges.size());
538
539 for(size_t i = 0; i < jacobiEdges.size(); i++) {
540
541 SimplexId const edgeId = jacobiEdges[i].first;
542
543 edge2polygonEdgeId[edgeId] = i;
544 jacobi2edges_[i] = edgeId;
545 }
546
547#ifdef TTK_ENABLE_OPENMP
548#pragma omp parallel for num_threads(threadNumber_)
549#endif
550 for(size_t i = 0; i < originalData_.sheet2List_.size(); i++) {
551
552 for(size_t j = 0;
553 j < originalData_.sheet1List_[originalData_.sheet2List_[i].sheet1Id_]
554 .edgeList_.size();
555 j++) {
556
557 SimplexId const edgeId
558 = originalData_.sheet1List_[originalData_.sheet2List_[i].sheet1Id_]
559 .edgeList_[j];
560
561 fiberSurface_.setTriangleList(
562 edge2polygonEdgeId[edgeId],
563 &(originalData_.sheet2List_[i].triangleList_[j]));
564 fiberSurface_.setVertexList(
565 edge2polygonEdgeId[edgeId],
566 &(originalData_.sheet2List_[i].vertexList_[j]));
567 }
568 }
569
570#ifdef TTK_ENABLE_OPENMP
571#pragma omp parallel for num_threads(threadNumber_)
572#endif
573 for(size_t i = 0; i < jacobiEdges.size(); i++) {
574
575 SimplexId const edgeId = jacobiEdges[i].first;
576
577 std::pair<double, double> rangePoint0, rangePoint1;
578
579 SimplexId vertexId0 = -1, vertexId1 = -1;
580 triangulation.getEdgeVertex(edgeId, 0, vertexId0);
581 triangulation.getEdgeVertex(edgeId, 1, vertexId1);
582
583 rangePoint0.first = uField[vertexId0];
584 rangePoint0.second = vField[vertexId0];
585
586 rangePoint1.first = uField[vertexId1];
587 rangePoint1.second = vField[vertexId1];
588
589 if(originalData_.edgeTypes_[edgeId] == 1) {
590 std::vector<SimplexId> edgeSeeds(
591 triangulation.getEdgeStarNumber(edgeId), -1);
592 for(size_t j = 0; j < edgeSeeds.size(); j++) {
593 triangulation.getEdgeStar(edgeId, j, edgeSeeds[j]);
594 }
595
596 fiberSurface_.computeContour<dataTypeU, dataTypeV>(
597 rangePoint0, rangePoint1, edgeSeeds, &triangulation,
598 edge2polygonEdgeId[edgeId]);
599 } else {
600#ifdef TTK_ENABLE_FIBER_SURFACE_WITH_RANGE_OCTREE
601 if(withRangeDrivenOctree_) {
602 fiberSurface_.computeSurfaceWithOctree<dataTypeU, dataTypeV>(
603 rangePoint0, rangePoint1, &triangulation, edge2polygonEdgeId[edgeId]);
604 } else {
605 fiberSurface_.computeSurface<dataTypeU, dataTypeV>(
606 rangePoint0, rangePoint1, &triangulation, edge2polygonEdgeId[edgeId]);
607 }
608#else
609 fiberSurface_.computeSurface<dataTypeU, dataTypeV>(
610 rangePoint0, rangePoint1, edge2polygonEdgeId[edgeId]);
611#endif
612 }
613 }
614
615 // #ifdef TTK_ENABLE_OPENMP
616 // #pragma omp parallel for num_threads(threadNumber_)
617 // #endif
618 // for(int i = 0; i < (int) originalData_.sheet2List_.size(); i++){
619 //
620 // int sheet1Id = originalData_.sheet2List_[i].sheet1Id_;
621 //
622 // std::vector<bool> inList(tetNumber_, false);
623 // std::vector<int> seedList;
624 // std::vector<std::pair<std::pair<double, double>, std::pair<double,
625 // double> > > edgeList;
626 // std::vector<int> jacobiEdgeIdList;
627 //
628 // for(int j = 0;
629 // j < (int) originalData_.sheet1List_[sheet1Id].edgeList_.size(); j++){
630 // int edgeId = originalData_.sheet1List_[sheet1Id].edgeList_[j];
631 //
632 // if(originalData_.edgeTypes_[edgeId] == 1){
633 // for(int k = 0; k < (int) (*edgeStars_)[edgeId].size(); k++){
634 // int tetId = (*edgeStars_)[edgeId][k];
635 // if(!inList[tetId]){
636 // seedList.push_back(tetId);
637 // }
638 // }
639 // }
640 //
641 // if((originalData_.edgeTypes_[edgeId] != 1)
642 // &&(originalData_.edgeTypes_[edgeId] != -1)){
643 // edgeList.resize(edgeList.size() + 1);
644 //
645 // edgeList.back().first.first =
646 // ((dataTypeU *) uField_)[(*edgeList_)[edgeId].first];
647 // edgeList.back().first.second =
648 // ((dataTypeV *) vField_)[(*edgeList_)[edgeId].first];
649 //
650 // edgeList.back().second.first =
651 // ((dataTypeU *) uField_)[(*edgeList_)[edgeId].second];
652 // edgeList.back().second.second =
653 // ((dataTypeV *) vField_)[(*edgeList_)[edgeId].second];
654 //
655 // jacobiEdgeIdList.push_back(edge2polygonEdgeId[edgeId]);
656 // }
657 // }
658 //
659 // if(seedList.size()){
660 // fiberSurface_.computeContour<dataTypeU, dataTypeV>(
661 // edgeList, seedList,
662 // &jacobiEdgeIdList);
663 // }
664 // }
665
666 this->printMsg(
667 "Computed fiber surfaces", 1.0, t.getElapsedTime(), this->threadNumber_);
668
669 return 0;
670}
671
672template <class dataTypeU, class dataTypeV, typename triangulationType>
674 const dataTypeU *const uField,
675 const dataTypeV *const vField,
676 const triangulationType &triangulation) {
677
678 this->printWrn("Computing chambers' pre-images.");
679 this->printWrn("This will take a LONG time.");
680
681 Timer t;
682
683 // at this point, they have exactly the same size
684 originalData_.sheet2List_.resize(threadNumber_);
685 for(size_t i = 0; i < originalData_.sheet2List_.size(); i++) {
686 originalData_.sheet2List_[i].sheet1Id_ = i;
687 originalData_.sheet2List_[i].pruned_ = false;
688 originalData_.sheet2List_[i].vertexList_.resize(1);
689 originalData_.sheet2List_[i].triangleList_.resize(1);
690 }
691
692 fiberSurface_.setGlobalVertexList(&fiberSurfaceVertexList_);
693 fiberSurface_.setTetNumber(tetNumber_);
694 fiberSurface_.setPolygonEdgeNumber(threadNumber_);
695
696 for(size_t i = 0; i < originalData_.sheet2List_.size(); i++) {
697
698 fiberSurface_.setTriangleList(
699 i, &(originalData_.sheet2List_[i].triangleList_[0]));
700 fiberSurface_.setVertexList(
701 i, &(originalData_.sheet2List_[i].vertexList_[0]));
702 }
703
704#ifdef TTK_ENABLE_OPENMP
705#pragma omp parallel for num_threads(threadNumber_)
706#endif
707 for(SimplexId i = 0; i < (SimplexId)edgeNumber_; i++) {
708
709#ifdef TTK_ENABLE_OPENMP
710 ThreadId threadId = omp_get_thread_num();
711#else
712 ThreadId threadId = 0;
713#endif
714
715 SimplexId edgeId = i;
716
717 std::pair<double, double> rangePoint0, rangePoint1;
718
719 SimplexId vertexId0 = -1, vertexId1 = -1;
720 triangulation.getEdgeVertex(edgeId, 0, vertexId0);
721 triangulation.getEdgeVertex(edgeId, 1, vertexId1);
722
723 rangePoint0.first = uField[vertexId0];
724 rangePoint0.second = vField[vertexId0];
725
726 rangePoint1.first = uField[vertexId1];
727 rangePoint1.second = vField[vertexId1];
728
729 fiberSurface_.computeSurface<dataTypeU, dataTypeV>(
730 rangePoint0, rangePoint1, threadId);
731
732 // clear the memory now.. otherwise we'll swap and never end
733 originalData_.sheet2List_[threadId].triangleList_[0].clear();
734 originalData_.sheet2List_[threadId].vertexList_[0].clear();
735 }
736
737 this->printMsg("Computed chambers pre-image boundaries", 1.0,
738 t.getElapsedTime(), this->threadNumber_);
739
740 return 0;
741}
742
743template <class dataTypeU, class dataTypeV, typename triangulationType>
745 Sheet3 &sheet,
746 const dataTypeU *const uField,
747 const dataTypeV *const vField,
748 const triangulationType &triangulation) {
749
750 sheet.domainVolume_ = 0;
751 sheet.rangeArea_ = 0;
752 sheet.hyperVolume_ = 0;
753
754 for(size_t i = 0; i < sheet.tetList_.size(); i++) {
755
756 SimplexId const tetId = sheet.tetList_[i];
757 std::array<std::pair<double, double>, 3> domainBox{};
758 std::array<std::pair<double, double>, 2> rangeBox;
759 std::array<std::array<float, 3>, 4> domainPoints{};
760 std::array<std::array<float, 2>, 4> rangePoints{};
761
762 for(int j = 0; j < 4; j++) {
763
764 SimplexId vertexId = -1;
765 triangulation.getCellVertex(tetId, j, vertexId);
766
767 triangulation.getVertexPoint(
768 vertexId, domainPoints[j][0], domainPoints[j][1], domainPoints[j][2]);
769
770 rangePoints[j][0] = uField[vertexId];
771 rangePoints[j][1] = vField[vertexId];
772 }
773
774 Geometry::getBoundingBox(domainPoints, domainBox);
775 Geometry::getBoundingBox(rangePoints, rangeBox);
776
777 sheet.domainVolume_ += (domainBox[0].second - domainBox[0].first)
778 * (domainBox[1].second - domainBox[1].first)
779 * (domainBox[2].second - domainBox[2].first);
780
781 sheet.rangeArea_ += (rangeBox[0].second - rangeBox[0].first)
782 * (rangeBox[1].second - rangeBox[1].first);
783 }
784
785 if(sheet.domainVolume_) {
786 sheet.hyperVolume_ = sheet.rangeArea_ / sheet.domainVolume_;
787 } else {
788 sheet.hyperVolume_ = 0;
789 }
790
791 return 0;
792}
793
794template <class dataTypeU, class dataTypeV>
795inline int ttk::ReebSpace::perturb(const dataTypeU *const uField,
796 const dataTypeV *const vField,
797 const dataTypeU &uEpsilon,
798 const dataTypeV &vEpsilon) {
799
800 jacobiSet_.setVertexNumber(vertexNumber_);
801 jacobiSet_.perturb(uField, vField, uEpsilon, vEpsilon);
802
803 return 0;
804}
805
806template <class dataTypeU, class dataTypeV, typename triangulationType>
807inline int ttk::ReebSpace::simplify(const dataTypeU *const uField,
808 const dataTypeV *const vField,
809 const triangulationType &triangulation,
810 const double &simplificationThreshold,
811 const SimplificationCriterion &criterion) {
812
813 if((totalArea_ == -1) || (totalVolume_ == -1) || (totalHyperVolume_ == -1)) {
814
815 Timer t;
816
817#ifdef TTK_ENABLE_OPENMP
818#pragma omp parallel for num_threads(threadNumber_)
819#endif
820 for(size_t i = 0; i < originalData_.sheet3List_.size(); i++) {
821 computeGeometricalMeasures(
822 originalData_.sheet3List_[i], uField, vField, triangulation);
823 }
824
825 for(size_t i = 0; i < originalData_.sheet3List_.size(); i++) {
826 totalArea_ += originalData_.sheet3List_[i].rangeArea_;
827 totalVolume_ += originalData_.sheet3List_[i].domainVolume_;
828 totalHyperVolume_ += originalData_.sheet3List_[i].hyperVolume_;
829 }
830
831 this->printMsg("Computed geometrical measures", 1.0, t.getElapsedTime(),
832 this->threadNumber_);
833 }
834
835 if(!hasConnectedSheets_) {
836 connectSheets(triangulation);
837 prepareSimplification();
838 }
839
840 {
841 std::stringstream msg;
842 msg << "Simplifying (";
843 switch(criterion) {
844 case SimplificationCriterion::domainVolume:
845 msg << "'Domain Volume'";
846 break;
847 case SimplificationCriterion::rangeArea:
848 msg << "'Range Area'";
849 break;
850 case SimplificationCriterion::hyperVolume:
851 msg << "'HyperVolume'";
852 break;
853 }
854 msg << ", thr: " << simplificationThreshold << ").";
855 this->printMsg(msg.str());
856 }
857
858 if(!((criterion == currentData_.simplificationCriterion_)
859 && (simplificationThreshold > currentData_.simplificationThreshold_))) {
860 prepareSimplification();
861 }
862
863 simplifySheets(simplificationThreshold, criterion, triangulation);
864
865 return 0;
866}
867
868template <typename triangulationType>
870 const std::vector<std::pair<SimplexId, char>> &jacobiSet,
871 std::vector<std::pair<SimplexId, SimplexId>> &jacobiClassification,
872 const triangulationType &triangulation) {
873
874 // bfs on the saddle jacobi edges only to identify saddle 1-sheets as well as
875 // saddle 0-sheets
876
877 Timer t;
878
879 // alloc
880 jacobiClassification.reserve(jacobiSet.size());
881
882 // markup the saddle jacobi and visisted edges
883 for(size_t i = 0; i < jacobiSet.size(); i++) {
884 originalData_.edgeTypes_[jacobiSet[i].first] = jacobiSet[i].second;
885 }
886
887 std::vector<bool> visitedEdges(triangulation.getNumberOfEdges(), false);
888
889 for(size_t i = 0; i < jacobiSet.size(); i++) {
890
891 if(visitedEdges[jacobiSet[i].first] == false) {
892
893 SimplexId const sheet1Id = originalData_.sheet1List_.size();
894 originalData_.sheet1List_.resize(originalData_.sheet1List_.size() + 1);
895 originalData_.sheet1List_.back().hasSaddleEdges_ = false;
896 originalData_.sheet1List_.back().pruned_ = false;
897
898 std::queue<SimplexId> edgeQueue;
899 edgeQueue.push(jacobiSet[i].first);
900
901 do {
902
903 SimplexId const edgeId = edgeQueue.front();
904 edgeQueue.pop();
905
906 if(!visitedEdges[edgeId]) {
907
908 jacobiClassification.emplace_back(edgeId, sheet1Id);
909
910 originalData_.sheet1List_.back().edgeList_.push_back(edgeId);
911 originalData_.edge2sheet1_[edgeId] = sheet1Id;
912
913 if(originalData_.edgeTypes_[edgeId] == 1) {
914 originalData_.sheet1List_.back().hasSaddleEdges_ = true;
915 }
916
917 SimplexId vertexId0 = -1;
918 triangulation.getEdgeVertex(edgeId, 0, vertexId0);
919 SimplexId vertexId1 = -1;
920 triangulation.getEdgeVertex(edgeId, 1, vertexId1);
921
922 SimplexId neighborNumber = 0;
923
924 SimplexId vertexEdgeNumber
925 = triangulation.getVertexEdgeNumber(vertexId0);
926
927 for(SimplexId j = 0; j < vertexEdgeNumber; j++) {
928
929 SimplexId vertexEdgeId = -1;
930 triangulation.getVertexEdge(vertexId0, j, vertexEdgeId);
931
932 if(originalData_.edgeTypes_[vertexEdgeId] != -1) {
933 if(vertexEdgeId != edgeId) {
934
935 neighborNumber++;
936
937 if(!visitedEdges[vertexEdgeId]) {
938 edgeQueue.push(vertexEdgeId);
939 }
940 }
941 }
942 }
943 if((neighborNumber > 2)
944 && (originalData_.vertex2sheet0_[vertexId0] == -1)) {
945 // vertexId0 is a 0-sheet
946 // mark up the segmentation
947 originalData_.vertex2sheet0_[vertexId0]
948 = originalData_.sheet0List_.size();
949
950 originalData_.sheet0List_.resize(originalData_.sheet0List_.size()
951 + 1);
952 originalData_.sheet0List_.back().vertexId_ = vertexId0;
953 originalData_.sheet0List_.back().type_ = 1;
954 originalData_.sheet0List_.back().pruned_ = false;
955 }
956
957 neighborNumber = 0;
958 vertexEdgeNumber = triangulation.getVertexEdgeNumber(vertexId1);
959 for(SimplexId j = 0; j < vertexEdgeNumber; j++) {
960
961 SimplexId vertexEdgeId = -1;
962 triangulation.getVertexEdge(vertexId1, j, vertexEdgeId);
963
964 if(originalData_.edgeTypes_[vertexEdgeId] != -1) {
965 if(vertexEdgeId != edgeId) {
966
967 neighborNumber++;
968
969 if(!visitedEdges[vertexEdgeId]) {
970 edgeQueue.push(vertexEdgeId);
971 }
972 }
973 }
974 }
975
976 if((neighborNumber > 2)
977 && (originalData_.vertex2sheet0_[vertexId1] == -1)) {
978 // vertexId1 is a 0-sheet
979 // mark up the segmentation
980
981 originalData_.vertex2sheet0_[vertexId1]
982 = originalData_.sheet0List_.size();
983
984 originalData_.sheet0List_.resize(originalData_.sheet0List_.size()
985 + 1);
986 originalData_.sheet0List_.back().vertexId_ = vertexId1;
987 originalData_.sheet0List_.back().type_ = 1;
988 originalData_.sheet0List_.back().pruned_ = false;
989 }
990 }
991
992 visitedEdges[edgeId] = true;
993
994 } while(edgeQueue.size());
995 }
996 }
997
998 this->printMsg(std::vector<std::vector<std::string>>{
999 {"#1-sheets", std::to_string(originalData_.sheet1List_.size())},
1000 {"#0-sheets", std::to_string(originalData_.sheet0List_.size())}});
1001 this->printMsg(
1002 "Extracted 0- and 1-sheets", 1.0, t.getElapsedTime(), this->threadNumber_);
1003
1004 return 0;
1005}
1006
1007template <typename triangulationType>
1009 const std::vector<std::pair<SimplexId, char>> &jacobiSet,
1010 std::vector<std::pair<SimplexId, SimplexId>> &jacobiClassification,
1011 const triangulationType &triangulation) {
1012
1013 // bfs on the saddle jacobi edges only to identify saddle 1-sheets as well as
1014 // saddle 0-sheets
1015
1016 Timer t;
1017
1018 // alloc
1019 jacobiClassification.reserve(jacobiSet.size());
1020
1021 // markup the saddle jacobi and visisted edges
1022 for(size_t i = 0; i < jacobiSet.size(); i++) {
1023 originalData_.edgeTypes_[jacobiSet[i].first] = jacobiSet[i].second;
1024 }
1025
1026 std::vector<bool> visitedEdges(triangulation.getNumberOfEdges(), false);
1027 for(size_t i = 0; i < jacobiSet.size(); i++) {
1028
1029 if(/*(saddleEdge[jacobiSet[i].first] == 1)&&*/
1030 (visitedEdges[jacobiSet[i].first] == false)) {
1031 // saddle, non-visited edge
1032
1033 // we have a seed here
1034 SimplexId const sheet1Id = originalData_.sheet1List_.size();
1035 originalData_.sheet1List_.resize(originalData_.sheet1List_.size() + 1);
1036 originalData_.sheet1List_.back().hasSaddleEdges_ = false;
1037 originalData_.sheet1List_.back().pruned_ = false;
1038
1039 std::queue<SimplexId> edgeQueue;
1040 edgeQueue.push(jacobiSet[i].first);
1041
1042 do {
1043
1044 SimplexId edgeId = edgeQueue.front();
1045 edgeQueue.pop();
1046
1047 if(!visitedEdges[edgeId]) {
1048
1049 jacobiClassification.emplace_back(edgeId, sheet1Id);
1050
1051 if(originalData_.edgeTypes_[edgeId] == 1) {
1052 // saddle edge
1053 originalData_.sheet1List_.back().hasSaddleEdges_ = true;
1054 }
1055 originalData_.sheet1List_.back().edgeList_.push_back(edgeId);
1056 originalData_.edge2sheet1_[edgeId] = sheet1Id;
1057
1058 // next grab its neighbors
1059 // if saddleEdge and not visited continue (only if one)
1060 SimplexId vertexId0 = -1;
1061 triangulation.getEdgeVertex(edgeId, 0, vertexId0);
1062 SimplexId vertexId1 = -1;
1063 triangulation.getEdgeVertex(edgeId, 1, vertexId1);
1064
1065 if(originalData_.vertex2sheet0_[vertexId0]
1066 >= static_cast<SimplexId>(originalData_.sheet0List_.size())) {
1067 // WEIRD BUG after multiple runs
1068 originalData_.vertex2sheet0_[vertexId0] = -1;
1069 }
1070
1071 // make sure these are not already visited 0-sheets
1072 if(originalData_.vertex2sheet0_[vertexId0] == -1) {
1073
1074 // inspect neighbor edges
1075 SimplexId neighborNumber = 0;
1076 SimplexId nextEdgeId = -1;
1077 SimplexId vertexEdgeNumber
1078 = triangulation.getVertexEdgeNumber(vertexId0);
1079 for(SimplexId j = 0; j < vertexEdgeNumber; j++) {
1080
1081 SimplexId vertexEdgeId = -1;
1082 triangulation.getVertexEdge(vertexId0, j, vertexEdgeId);
1083
1084 if(originalData_.edgeTypes_[vertexEdgeId] != -1) {
1085 neighborNumber++;
1086 if(vertexEdgeId != edgeId
1087 /*&&(saddleEdge[(*vertexEdgeList_)[vertexId0][j]] == 1)*/) {
1088 nextEdgeId = vertexEdgeId;
1089 }
1090 }
1091 }
1092
1093 if((neighborNumber == 2) && (nextEdgeId != -1)) {
1094 // this is a regular vertex and we can add the next edge to the
1095 // queue if it's not been visited before.
1096 if(!visitedEdges[nextEdgeId])
1097 edgeQueue.push(nextEdgeId);
1098 } else {
1099 // we hit a 0-sheet that hasn't been visited before.
1100 // let's create it
1101
1102 // mark up the segmentation
1103 originalData_.vertex2sheet0_[vertexId0]
1104 = originalData_.sheet0List_.size();
1105
1106 originalData_.sheet0List_.resize(originalData_.sheet0List_.size()
1107 + 1);
1108 originalData_.sheet0List_.back().vertexId_ = vertexId0;
1109 originalData_.sheet0List_.back().type_ = 1;
1110 originalData_.sheet0List_.back().pruned_ = false;
1111 }
1112 }
1113
1114 if(originalData_.vertex2sheet0_[vertexId0] != -1) {
1115 // we hit a 0-sheet
1116 // attach it to the one sheet and reciprocally
1117
1118 // attach the 0-sheet to the 1-sheet
1119 originalData_.sheet1List_[sheet1Id].sheet0List_.push_back(
1120 originalData_.vertex2sheet0_[vertexId0]);
1121
1122 // attach the 1-sheet to the 0-sheet
1123 originalData_.sheet0List_[originalData_.vertex2sheet0_[vertexId0]]
1124 .sheet1List_.push_back(sheet1Id);
1125 }
1126
1127 if(originalData_.vertex2sheet0_[vertexId1]
1128 >= static_cast<SimplexId>(originalData_.sheet0List_.size())) {
1129 // WEIRD BUG after multiple runs
1130 originalData_.vertex2sheet0_[vertexId1] = -1;
1131 }
1132 // now do the same for the other vertex
1133 // make sure these are not already visited 0-sheets
1134 if(originalData_.vertex2sheet0_[vertexId1] == -1) {
1135
1136 // inspect neighbor edges
1137 SimplexId neighborNumber = 0;
1138 SimplexId nextEdgeId = -1;
1139
1140 SimplexId vertexEdgeNumber
1141 = triangulation.getVertexEdgeNumber(vertexId1);
1142
1143 for(SimplexId j = 0; j < vertexEdgeNumber; j++) {
1144
1145 SimplexId vertexEdgeId = -1;
1146 triangulation.getVertexEdge(vertexId1, j, vertexEdgeId);
1147
1148 if(originalData_.edgeTypes_[vertexEdgeId] != -1) {
1149 neighborNumber++;
1150 if(vertexEdgeId != edgeId
1151 /*&&(saddleEdge[(*vertexEdgeList_)[vertexId1][j]] == 1)*/) {
1152 nextEdgeId = vertexEdgeId;
1153 }
1154 }
1155 }
1156
1157 if((neighborNumber == 2) && (nextEdgeId != -1)) {
1158 // this is a regular vertex and we can add the next edge to the
1159 // queue if it's not been visited before.
1160 if(!visitedEdges[nextEdgeId])
1161 edgeQueue.push(nextEdgeId);
1162 } else {
1163 // we hit a 0-sheet that hasn't been visited before.
1164 // let's create it
1165
1166 // mark up the segmentation
1167 originalData_.vertex2sheet0_[vertexId1]
1168 = originalData_.sheet0List_.size();
1169
1170 originalData_.sheet0List_.resize(originalData_.sheet0List_.size()
1171 + 1);
1172 originalData_.sheet0List_.back().vertexId_ = vertexId1;
1173 originalData_.sheet0List_.back().type_ = 1;
1174 originalData_.sheet0List_.back().pruned_ = false;
1175 }
1176 }
1177
1178 if(originalData_.vertex2sheet0_[vertexId1] != -1) {
1179 // we hit a 0-sheet
1180 // attach it to the one sheet and reciprocally
1181
1182 // attach the 0-sheet to the 1-sheet
1183 originalData_.sheet1List_[sheet1Id].sheet0List_.push_back(
1184 originalData_.vertex2sheet0_[vertexId1]);
1185
1186 // attach the 1-sheet to the 0-sheet
1187 originalData_.sheet0List_[originalData_.vertex2sheet0_[vertexId1]]
1188 .sheet1List_.push_back(sheet1Id);
1189 }
1190
1191 visitedEdges[edgeId] = true;
1192 }
1193
1194 } while(edgeQueue.size());
1195 }
1196 }
1197
1198 this->printMsg(std::vector<std::vector<std::string>>{
1199 {"#1-sheets", std::to_string(originalData_.sheet1List_.size())},
1200 {"#0-sheets", std::to_string(originalData_.sheet0List_.size())}});
1201 this->printMsg(
1202 "Extracted 0- and 1-sheets", 1.0, t.getElapsedTime(), this->threadNumber_);
1203
1204 return 0;
1205}
1206
1207template <typename triangulationType>
1209 const SimplexId &vertexId,
1210 const std::vector<std::vector<std::array<SimplexId, 3>>> &tetTriangles,
1211 const triangulationType &triangulation) {
1212
1213 SimplexId const sheetId = originalData_.sheet3List_.size();
1214 originalData_.sheet3List_.resize(originalData_.sheet3List_.size() + 1);
1215 originalData_.sheet3List_.back().pruned_ = false;
1216 originalData_.sheet3List_.back().preMerger_ = -1;
1217 originalData_.sheet3List_.back().Id_ = sheetId;
1218
1219 std::queue<SimplexId> vertexQueue;
1220 vertexQueue.push(vertexId);
1221
1222 do {
1223
1224 SimplexId const localVertexId = vertexQueue.front();
1225 vertexQueue.pop();
1226
1227 if(originalData_.vertex2sheet3_[localVertexId] == -1) {
1228 // not visited yet
1229
1230 originalData_.sheet3List_.back().vertexList_.push_back(localVertexId);
1231 originalData_.vertex2sheet3_[localVertexId] = sheetId;
1232
1233 SimplexId const vertexStarNumber
1234 = triangulation.getVertexStarNumber(localVertexId);
1235
1236 for(SimplexId i = 0; i < vertexStarNumber; i++) {
1237 SimplexId tetId = -1;
1238 triangulation.getVertexStar(localVertexId, i, tetId);
1239
1240 if(tetTriangles[tetId].empty()) {
1241 // if(originalData_.tet2sheet3_[tetId] == -1){
1242 // originalData_.tet2sheet3_[tetId] = sheetId;
1243 // originalData_.sheet3List_.back().tetList_.push_back(tetId);
1244 // }
1245 for(int j = 0; j < 4; j++) {
1246
1247 SimplexId tetVertexId = -1;
1248 triangulation.getCellVertex(tetId, j, tetVertexId);
1249
1250 if(originalData_.vertex2sheet3_[tetVertexId] == -1) {
1251 vertexQueue.push(tetVertexId);
1252 }
1253 }
1254 } else {
1255 // fiber surface in there
1256 for(int j = 0; j < 4; j++) {
1257 SimplexId otherVertexId = -1;
1258 triangulation.getCellVertex(tetId, j, otherVertexId);
1259 if((otherVertexId != localVertexId)
1260 && (originalData_.vertex2sheet3_[otherVertexId] == -1)) {
1261 // we need to see if the edge <localVertexId, otherVertexId> is
1262 // cut by a fiber surface triangle or not.
1263
1264 bool isCut = false;
1265 for(size_t k = 0; k < tetTriangles[tetId].size(); k++) {
1266 SimplexId l = 0, m = 0, n = 0;
1267 l = tetTriangles[tetId][k][0];
1268 m = tetTriangles[tetId][k][1];
1269 n = tetTriangles[tetId][k][2];
1270
1271 for(int p = 0; p < 3; p++) {
1272 std::pair<SimplexId, SimplexId> meshEdge;
1273
1274 if(fiberSurfaceVertexList_.size()) {
1275 // the fiber surfaces have been merged
1276 meshEdge
1277 = fiberSurfaceVertexList_[originalData_.sheet2List_[l]
1278 .triangleList_[m][n]
1279 .vertexIds_[p]]
1280 .meshEdge_;
1281 } else {
1282 // the fiber surfaces have not been merged
1283 meshEdge = originalData_.sheet2List_[l]
1284 .vertexList_[m][originalData_.sheet2List_[l]
1285 .triangleList_[m][n]
1286 .vertexIds_[p]]
1287 .meshEdge_;
1288 }
1289
1290 if(((meshEdge.first == localVertexId)
1291 && (meshEdge.second == otherVertexId))
1292 || ((meshEdge.second == localVertexId)
1293 && (meshEdge.first == otherVertexId))) {
1294 isCut = true;
1295 break;
1296 }
1297 }
1298
1299 if(isCut)
1300 break;
1301 }
1302
1303 if(!isCut) {
1304 // add the vertex to the queue
1305 vertexQueue.push(otherVertexId);
1306 }
1307 }
1308 }
1309 }
1310 }
1311 }
1312
1313 } while(vertexQueue.size());
1314
1315 return 0;
1316}
1317
1318template <typename triangulationType>
1320 std::vector<std::vector<std::array<SimplexId, 3>>> &tetTriangles,
1321 const triangulationType &triangulation) {
1322
1323 Timer t;
1324
1325 tetTriangles.resize(tetNumber_);
1326
1327 for(size_t i = 0; i < originalData_.sheet2List_.size(); i++) {
1328 for(size_t j = 0; j < originalData_.sheet2List_[i].triangleList_.size();
1329 j++) {
1330 for(size_t k = 0;
1331 k < originalData_.sheet2List_[i].triangleList_[j].size(); k++) {
1332
1333 SimplexId const tetId
1334 = originalData_.sheet2List_[i].triangleList_[j][k].tetId_;
1335
1336 tetTriangles[tetId].emplace_back();
1337 tetTriangles[tetId].back()[0] = i;
1338 tetTriangles[tetId].back()[1] = j;
1339 tetTriangles[tetId].back()[2] = k;
1340 }
1341 }
1342 }
1343
1344 // mark all the jacobi edge vertices
1345 for(size_t i = 0; i < originalData_.sheet1List_.size(); i++) {
1346 for(size_t j = 0; j < originalData_.sheet1List_[i].edgeList_.size(); j++) {
1347
1348 SimplexId const edgeId = originalData_.sheet1List_[i].edgeList_[j];
1349
1350 SimplexId vertexId0 = -1, vertexId1 = -1;
1351 triangulation.getEdgeVertex(edgeId, 0, vertexId0);
1352 triangulation.getEdgeVertex(edgeId, 1, vertexId1);
1353
1354 originalData_.vertex2sheet3_[vertexId0] = -2 - i;
1355 originalData_.vertex2sheet3_[vertexId1] = -2 - i;
1356 }
1357 }
1358
1359 for(SimplexId i = 0; i < vertexNumber_; i++) {
1360 if(originalData_.vertex2sheet3_[i] == -1) {
1361 compute3sheet(i, tetTriangles, triangulation);
1362 }
1363 }
1364
1365 // for 3-sheet expansion
1366 std::vector<std::vector<std::pair<SimplexId, bool>>> neighborList(
1367 originalData_.sheet3List_.size());
1368 // end of 3-sheet expansion
1369
1370 // add the tets in parallel
1371#ifdef TTK_ENABLE_OPENMP
1372#pragma omp parallel for num_threads(threadNumber_)
1373#endif
1374 for(size_t i = 0; i < originalData_.sheet3List_.size(); i++) {
1375 for(size_t j = 0; j < originalData_.sheet3List_[i].vertexList_.size();
1376 j++) {
1377
1378 SimplexId const vertexId = originalData_.sheet3List_[i].vertexList_[j];
1379 SimplexId const sheetId = originalData_.vertex2sheet3_[vertexId];
1380
1381 SimplexId const vertexStarNumber
1382 = triangulation.getVertexStarNumber(vertexId);
1383
1384 for(SimplexId k = 0; k < vertexStarNumber; k++) {
1385 SimplexId tetId = -1;
1386 triangulation.getVertexStar(vertexId, k, tetId);
1387 if(tetTriangles[tetId].empty()) {
1388 if(originalData_.tet2sheet3_[tetId] == -1) {
1389 originalData_.tet2sheet3_[tetId] = i;
1390 originalData_.sheet3List_[i].tetList_.push_back(tetId);
1391 }
1392 } else {
1393
1394 // expending here 3-sheets.
1395 if(expand3sheets_) {
1396 for(int l = 0; l < 4; l++) {
1397 SimplexId otherVertexId = -1;
1398
1399 triangulation.getCellVertex(tetId, l, otherVertexId);
1400
1401 if(vertexId != otherVertexId) {
1402
1403 SimplexId const otherSheetId
1404 = originalData_.vertex2sheet3_[otherVertexId];
1405
1406 if((sheetId != otherSheetId) && (otherSheetId >= 0)) {
1407
1408 bool inThere = false;
1409 for(size_t m = 0; m < neighborList[sheetId].size(); m++) {
1410 if(neighborList[sheetId][m].first == otherSheetId) {
1411 inThere = true;
1412 break;
1413 }
1414 }
1415
1416 if(!inThere) {
1417 neighborList[sheetId].emplace_back(otherSheetId, true);
1418 }
1419
1420 for(size_t m = 0; m < tetTriangles[tetId].size(); m++) {
1421
1422 // see if this guy is a saddle
1423 SimplexId const x = tetTriangles[tetId][m][0];
1424 SimplexId const y = tetTriangles[tetId][m][1];
1425 SimplexId const z = tetTriangles[tetId][m][2];
1426
1428 = &(originalData_.sheet2List_[x].triangleList_[y][z]);
1429
1430 bool cuttingTriangle = false;
1431 for(int n = 0; n < 3; n++) {
1433 = &(originalData_.sheet2List_[x]
1434 .vertexList_[y][tr->vertexIds_[n]]);
1435
1436 if(((v->meshEdge_.first == vertexId)
1437 && (v->meshEdge_.second == otherVertexId))
1438 || ((v->meshEdge_.second == vertexId)
1439 && (v->meshEdge_.first == otherVertexId))) {
1440
1441 cuttingTriangle = true;
1442 break;
1443 }
1444 }
1445
1446 if(cuttingTriangle) {
1447
1448 SimplexId const polygonId = tr->polygonEdgeId_;
1449 SimplexId const edgeId = jacobi2edges_[polygonId];
1450 if(originalData_.edgeTypes_[edgeId] == 1) {
1451 // this is a saddle Jacobi edge
1452
1453 inThere = false;
1454 for(size_t n = 0; n < neighborList[sheetId].size();
1455 n++) {
1456
1457 if(neighborList[sheetId][n].first == otherSheetId) {
1458 if(neighborList[sheetId][n].second == true) {
1459 neighborList[sheetId][n].second = false;
1460 }
1461 inThere = true;
1462 break;
1463 }
1464 }
1465
1466 if(!inThere) {
1467 neighborList[sheetId].emplace_back(
1468 otherSheetId, false);
1469 }
1470 }
1471 }
1472 }
1473 }
1474 }
1475 }
1476 }
1477 // end of the 3-sheet expansion
1478 }
1479 }
1480 }
1481 }
1482
1483 SimplexId totalSheetNumber = 0;
1484
1485 if(expand3sheets_) {
1486
1487 totalSheetNumber = originalData_.sheet3List_.size();
1488
1489 // expending sheets
1490 for(size_t i = 0; i < originalData_.sheet3List_.size(); i++) {
1491 if(originalData_.sheet3List_[i].pruned_ == false) {
1492
1493 for(size_t j = 0; j < neighborList[i].size(); j++) {
1494
1495 if(neighborList[i][j].second) {
1496
1497 bool isForbidden = false;
1498 SimplexId neighborId = neighborList[i][j].first;
1499
1500 // get the sheet where this guys has been merged
1501 while(originalData_.sheet3List_[neighborId].preMerger_ != -1) {
1502
1503 neighborId = originalData_.sheet3List_[neighborId].preMerger_;
1504 }
1505
1506 // make sure that no forbidden neighbor has been merged in the
1507 // candidate neighbor
1508 for(size_t k = 0;
1509 k
1510 < originalData_.sheet3List_[neighborId].preMergedSheets_.size();
1511 k++) {
1512
1513 SimplexId const subNeighborId
1514 = originalData_.sheet3List_[neighborId].preMergedSheets_[k];
1515
1516 for(size_t l = 0; l < neighborList[i].size(); l++) {
1517 if((neighborList[i][l].first == subNeighborId)
1518 && (!neighborList[i][l].second)) {
1519 isForbidden = true;
1520 break;
1521 }
1522 }
1523 if(isForbidden)
1524 break;
1525 }
1526
1527 // make sure that neighborId is not a candidate for a merge with a
1528 // sheet that is forbidden for i
1529 if(!isForbidden) {
1530 for(size_t k = 0; k < neighborList[i].size(); k++) {
1531 if(!neighborList[i][k].second) {
1532 SimplexId const forbiddenNeighbor = neighborList[i][k].first;
1533
1534 // make sure forbiddenNeighbor is not a valid merger for
1535 // neighborId
1536 for(size_t l = 0; l < neighborList[neighborId].size(); l++) {
1537 if((forbiddenNeighbor == neighborList[neighborId][l].first)
1538 && (neighborList[neighborId][l].second)) {
1539 isForbidden = true;
1540 break;
1541 }
1542 }
1543 if(isForbidden)
1544 break;
1545 }
1546 }
1547 }
1548
1549 if((neighborId != static_cast<SimplexId>(i)) && (!isForbidden)
1550 && (originalData_.sheet3List_[neighborId].pruned_ == false)
1551 && (originalData_.sheet3List_[neighborId].vertexList_.size()
1552 > originalData_.sheet3List_[i].vertexList_.size())) {
1553
1554 // these can be merged
1555 preMergeSheets(i, neighborId);
1556 totalSheetNumber--;
1557 break;
1558 }
1559 }
1560 }
1561 }
1562 }
1563 } else {
1564 totalSheetNumber = originalData_.sheet3List_.size();
1565 }
1566 // end of the 3-sheet expansion
1567
1568 this->printMsg("Computed " + std::to_string(totalSheetNumber) + " 3-sheets",
1569 1.0, t.getElapsedTime(), this->threadNumber_);
1570
1571 return 0;
1572}
1573
1574template <typename triangulationType>
1575int ttk::ReebSpace::connectSheets(const triangulationType &triangulation) {
1576
1577 Timer tm;
1578
1579 for(size_t i = 0; i < originalData_.sheet2List_.size(); i++) {
1580 for(size_t j = 0; j < originalData_.sheet2List_[i].triangleList_.size();
1581 j++) {
1582
1583 for(size_t k = 0;
1584 k < originalData_.sheet2List_[i].triangleList_[j].size(); k++) {
1585
1586 SimplexId const tetId
1587 = originalData_.sheet2List_[i].triangleList_[j][k].tetId_;
1588
1589 for(int l = 0; l < 4; l++) {
1590 SimplexId vertexId = -1;
1591 triangulation.getCellVertex(tetId, l, vertexId);
1592
1593 SimplexId const sheet3Id = originalData_.vertex2sheet3_[vertexId];
1594
1595 if(sheet3Id >= 0) {
1596 connect3sheetTo2sheet(originalData_, sheet3Id, i);
1597 }
1598 }
1599 }
1600 }
1601 }
1602
1603 // connect 3-sheets together
1604 for(SimplexId i = 0; i < vertexNumber_; i++) {
1605 if(originalData_.vertex2sheet3_[i] >= 0) {
1606
1607 SimplexId const vertexEdgeNumber = triangulation.getVertexEdgeNumber(i);
1608
1609 for(SimplexId j = 0; j < vertexEdgeNumber; j++) {
1610 SimplexId edgeId = -1;
1611 triangulation.getVertexEdge(i, j, edgeId);
1612 SimplexId otherVertexId = -1;
1613 triangulation.getEdgeVertex(edgeId, 0, otherVertexId);
1614 if(otherVertexId == i) {
1615 triangulation.getEdgeVertex(edgeId, 1, otherVertexId);
1616 }
1617
1618 if(originalData_.vertex2sheet3_[otherVertexId] >= 0) {
1619 if(originalData_.vertex2sheet3_[otherVertexId]
1620 != originalData_.vertex2sheet3_[i]) {
1621
1622 connect3sheetTo3sheet(originalData_,
1623 originalData_.vertex2sheet3_[i],
1624 originalData_.vertex2sheet3_[otherVertexId]);
1625 }
1626 }
1627
1628 if(originalData_.vertex2sheet0_[otherVertexId] != -1) {
1629 connect3sheetTo0sheet(originalData_, originalData_.vertex2sheet3_[i],
1630 originalData_.vertex2sheet0_[otherVertexId]);
1631 }
1632
1633 if(originalData_.vertex2sheet3_[otherVertexId] < -1) {
1634 SimplexId const sheet1Id
1635 = -2 - originalData_.vertex2sheet3_[otherVertexId];
1636 connect3sheetTo1sheet(
1637 originalData_, originalData_.vertex2sheet3_[i], sheet1Id);
1638 }
1639 }
1640 }
1641 }
1642
1643 this->printMsg(
1644 "Sheet connectivity established.", 1.0, tm.getElapsedTime(), 1);
1645
1646 printConnectivity(originalData_);
1647
1648 hasConnectedSheets_ = true;
1649
1650 return 0;
1651}
1652
1653template <typename triangulationType>
1655 ReebSpaceData &data,
1656 const SimplexId &sheet3Id,
1657 const SimplexId &sheet1Id,
1658 const SimplexId &biggerId,
1659 const triangulationType &triangulation) {
1660
1661 std::vector<SimplexId> newList;
1662
1663 newList.reserve(data.sheet1List_[sheet1Id].sheet3List_.size());
1664 for(size_t i = 0; i < data.sheet1List_[sheet1Id].sheet3List_.size(); i++) {
1665 SimplexId const other3SheetId = data.sheet1List_[sheet1Id].sheet3List_[i];
1666 if((other3SheetId != sheet3Id) && (!data.sheet3List_[other3SheetId].pruned_)
1667 && (data.sheet3List_[other3SheetId].tetList_.size())) {
1668 newList.push_back(data.sheet1List_[sheet1Id].sheet3List_[i]);
1669 }
1670 }
1671
1672 data.sheet1List_[sheet1Id].sheet3List_ = newList;
1673
1674 if((data.sheet1List_[sheet1Id].hasSaddleEdges_)
1675 && (data.sheet1List_[sheet1Id].sheet3List_.size() == 1)) {
1676 // this guy is no longer separating any body.
1677 data.sheet1List_[sheet1Id].pruned_ = true;
1678 data.sheet2List_[sheet1Id].pruned_ = true;
1679
1680 // update the segmentation
1681 for(size_t i = 0; i < data.sheet1List_[sheet1Id].edgeList_.size(); i++) {
1682
1683 SimplexId vertexId = -1;
1684 triangulation.getEdgeVertex(
1685 data.sheet1List_[sheet1Id].edgeList_[i], 0, vertexId);
1686 data.vertex2sheet3_[vertexId] = biggerId;
1687
1688 vertexId = -1;
1689 triangulation.getEdgeVertex(
1690 data.sheet1List_[sheet1Id].edgeList_[i], 1, vertexId);
1691 data.vertex2sheet3_[vertexId] = biggerId;
1692 }
1693
1694 for(size_t i = 0; i < data.sheet1List_[sheet1Id].sheet0List_.size(); i++) {
1695
1696 disconnect1sheetFrom0sheet(
1697 data, sheet1Id, data.sheet1List_[sheet1Id].sheet0List_[i], biggerId);
1698 }
1699 }
1700
1701 return 0;
1702}
1703
1704template <typename triangulationType>
1705int ttk::ReebSpace::flush(const triangulationType &triangulation) {
1706
1707 totalArea_ = -1;
1708 totalVolume_ = -1;
1709 totalHyperVolume_ = -1;
1710 hasConnectedSheets_ = false;
1711
1712 // store the segmentation for later purpose
1713 originalData_.vertex2sheet0_.resize(vertexNumber_);
1714 for(SimplexId i = 0; i < vertexNumber_; i++)
1715 originalData_.vertex2sheet0_[i] = -1;
1716
1717 originalData_.vertex2sheet3_.resize(vertexNumber_);
1718 for(SimplexId i = 0; i < vertexNumber_; i++)
1719 originalData_.vertex2sheet3_[i] = -1;
1720
1721 originalData_.edge2sheet1_.resize(triangulation.getNumberOfEdges(), -1);
1722 originalData_.edgeTypes_.resize(triangulation.getNumberOfEdges(), -1);
1723
1724 originalData_.tet2sheet3_.resize(tetNumber_, -1);
1725
1726 jacobi2edges_.clear();
1727 jacobiSetEdges_.clear();
1728
1729 originalData_.sheet0List_.clear();
1730 originalData_.sheet1List_.clear();
1731 originalData_.sheet2List_.clear();
1732 originalData_.sheet3List_.clear();
1733
1734 fiberSurfaceVertexList_.clear();
1735
1736 return 0;
1737}
1738
1739template <typename triangulationType>
1741 const SimplexId &biggerId,
1742 const triangulationType &triangulation) {
1743
1744 // 1. add the vertices and tets of smaller to bigger
1745 for(size_t i = 0; i < currentData_.sheet3List_[smallerId].vertexList_.size();
1746 i++) {
1747 SimplexId const vertexId
1748 = currentData_.sheet3List_[smallerId].vertexList_[i];
1749 currentData_.sheet3List_[biggerId].vertexList_.push_back(vertexId);
1750 currentData_.vertex2sheet3_[vertexId] = biggerId;
1751 }
1752 for(size_t i = 0; i < currentData_.sheet3List_[smallerId].tetList_.size();
1753 i++) {
1754 SimplexId const tetId = currentData_.sheet3List_[smallerId].tetList_[i];
1755 currentData_.sheet3List_[biggerId].tetList_.push_back(tetId);
1756 currentData_.tet2sheet3_[tetId] = biggerId;
1757 }
1758
1759 // 2. update bigger's score and re-insert it in the candidate list
1760 currentData_.sheet3List_[biggerId].domainVolume_
1761 += currentData_.sheet3List_[smallerId].domainVolume_;
1762 currentData_.sheet3List_[biggerId].rangeArea_
1763 += currentData_.sheet3List_[smallerId].rangeArea_;
1764 currentData_.sheet3List_[biggerId].hyperVolume_
1765 += currentData_.sheet3List_[smallerId].hyperVolume_;
1766
1767 // 3. add smaller's connections to bigger's
1768 for(size_t i = 0; i < currentData_.sheet3List_[smallerId].sheet3List_.size();
1769 i++) {
1770
1771 SimplexId const otherSheetId
1772 = currentData_.sheet3List_[smallerId].sheet3List_[i];
1773
1774 if(otherSheetId != biggerId) {
1775 connect3sheetTo3sheet(currentData_, biggerId, otherSheetId);
1776 }
1777 }
1778
1779 for(size_t i = 0; i < currentData_.sheet3List_[smallerId].sheet2List_.size();
1780 i++) {
1781
1782 SimplexId const otherSheetId
1783 = currentData_.sheet3List_[smallerId].sheet2List_[i];
1784
1785 if(otherSheetId != biggerId) {
1786 connect3sheetTo2sheet(currentData_, biggerId, otherSheetId);
1787 }
1788 }
1789
1790 for(size_t i = 0; i < currentData_.sheet3List_[smallerId].sheet1List_.size();
1791 i++) {
1792
1793 SimplexId const otherSheetId
1794 = currentData_.sheet3List_[smallerId].sheet1List_[i];
1795
1796 if(otherSheetId != biggerId) {
1797 connect3sheetTo1sheet(currentData_, biggerId, otherSheetId);
1798 }
1799 }
1800
1801 for(size_t i = 0; i < currentData_.sheet3List_[smallerId].sheet0List_.size();
1802 i++) {
1803
1804 SimplexId const otherSheetId
1805 = currentData_.sheet3List_[smallerId].sheet0List_[i];
1806
1807 if(otherSheetId != biggerId) {
1808 connect3sheetTo0sheet(currentData_, biggerId, otherSheetId);
1809 }
1810 }
1811
1812 currentData_.sheet3List_[smallerId].pruned_ = true;
1813
1814 // 4. disconnect smaller from all of its connections.
1815 for(size_t i = 0; i < currentData_.sheet3List_[smallerId].sheet3List_.size();
1816 i++) {
1817 disconnect3sheetFrom3sheet(
1818 currentData_, smallerId,
1819 currentData_.sheet3List_[smallerId].sheet3List_[i]);
1820 }
1821
1822 for(size_t i = 0; i < currentData_.sheet3List_[smallerId].sheet2List_.size();
1823 i++) {
1824 disconnect3sheetFrom2sheet(
1825 currentData_, smallerId,
1826 currentData_.sheet3List_[smallerId].sheet2List_[i]);
1827 }
1828
1829 for(size_t i = 0; i < currentData_.sheet3List_[smallerId].sheet1List_.size();
1830 i++) {
1831 disconnect3sheetFrom1sheet(
1832 currentData_, smallerId,
1833 currentData_.sheet3List_[smallerId].sheet1List_[i], biggerId,
1834 triangulation);
1835 }
1836
1837 for(size_t i = 0; i < currentData_.sheet3List_[smallerId].sheet0List_.size();
1838 i++) {
1839 disconnect3sheetFrom0sheet(
1840 currentData_, smallerId,
1841 currentData_.sheet3List_[smallerId].sheet0List_[i]);
1842 }
1843
1844 return 0;
1845}
1846
1847template <typename triangulationType>
1849 const SimplificationCriterion &criterion,
1850 const triangulationType &triangulation) {
1851
1852 SimplexId candidateId = -1;
1853 double maximumScore = -1;
1854
1855 // see the adjacent 3-sheets
1856 for(size_t i = 0; i < currentData_.sheet3List_[sheetId].sheet3List_.size();
1857 i++) {
1858 SimplexId const otherSheetId
1859 = currentData_.sheet3List_[sheetId].sheet3List_[i];
1860 if((!currentData_.sheet3List_[otherSheetId].pruned_)
1861 && (sheetId != otherSheetId)) {
1862
1863 double otherScore = 0;
1864
1865 switch(criterion) {
1866 case SimplificationCriterion::domainVolume:
1867 otherScore = currentData_.sheet3List_[otherSheetId].domainVolume_;
1868 break;
1869 case SimplificationCriterion::rangeArea:
1870 otherScore = currentData_.sheet3List_[otherSheetId].rangeArea_;
1871 break;
1872 case SimplificationCriterion::hyperVolume:
1873 otherScore = currentData_.sheet3List_[otherSheetId].hyperVolume_;
1874 break;
1875 }
1876
1877 if((maximumScore < 0) || (otherScore > maximumScore)) {
1878 candidateId = otherSheetId;
1879 maximumScore = otherScore;
1880 }
1881 }
1882 }
1883
1884 // see adjacent 3-sheets along 1-sheets
1885 // for(int i = 0;
1886 // i < (int) currentData_.sheet3List_[sheetId].sheet1List_.size(); i++){
1887 //
1888 // int sheet1Id = currentData_.sheet3List_[sheetId].sheet1List_[i];
1889 //
1890 // for(int j = 0;
1891 // j < (int) currentData_.sheet1List_[sheet1Id].sheet3List_.size();
1892 // j++){
1893 //
1894 // int otherSheetId = currentData_.sheet1List_[sheet1Id].sheet3List_[j];
1895 // if((otherSheetId != sheetId)
1896 // &&(!currentData_.sheet3List_[otherSheetId].pruned_)){
1897 //
1898 // double otherScore = 0;
1899 //
1900 // switch(criterion){
1901 // case domainVolume:
1902 // otherScore =
1903 // currentData_.sheet3List_[otherSheetId].domainVolume_;
1904 // break;
1905 // case rangeArea:
1906 // otherScore =
1907 // currentData_.sheet3List_[otherSheetId].rangeArea_;
1908 // break;
1909 // case hyperVolume:
1910 // otherScore =
1911 // currentData_.sheet3List_[otherSheetId].hyperVolume_;
1912 // break;
1913 // }
1914 //
1915 // if((maximumScore < 0)||(otherScore > maximumScore)){
1916 // candidateId = otherSheetId;
1917 // maximumScore = otherScore;
1918 // }
1919 // }
1920 // }
1921 // }
1922
1923 // see adjacent 3-sheets on 0-sheets
1924 // for(int i = 0;
1925 // i < (int) currentData_.sheet3List_[sheetId].sheet0List_.size(); i++){
1926 //
1927 // int sheet0Id = currentData_.sheet3List_[sheetId].sheet0List_[i];
1928 //
1929 // for(int j = 0;
1930 // j < (int) currentData_.sheet0List_[sheet0Id].sheet3List_.size();
1931 // j++){
1932 //
1933 // int otherSheetId = currentData_.sheet0List_[sheet0Id].sheet3List_[j];
1934 // if((otherSheetId != sheetId)
1935 // &&(!currentData_.sheet3List_[otherSheetId].pruned_)){
1936 //
1937 // double otherScore = 0;
1938 //
1939 // switch(criterion){
1940 // case domainVolume:
1941 // otherScore =
1942 // currentData_.sheet3List_[otherSheetId].domainVolume_;
1943 // break;
1944 // case rangeArea:
1945 // otherScore =
1946 // currentData_.sheet3List_[otherSheetId].rangeArea_;
1947 // break;
1948 // case hyperVolume:
1949 // otherScore =
1950 // currentData_.sheet3List_[otherSheetId].hyperVolume_;
1951 // break;
1952 // }
1953 //
1954 // if((maximumScore < 0)||(otherScore > maximumScore)){
1955 // candidateId = otherSheetId;
1956 // maximumScore = otherScore;
1957 // }
1958 // }
1959 // }
1960 // }
1961
1962 // mark as pruned if so
1963 if(candidateId != -1) {
1964 mergeSheets(sheetId, candidateId, triangulation);
1965 }
1966
1967 currentData_.sheet3List_[sheetId].pruned_ = true;
1968
1969 return 0;
1970}
1971
1972template <typename triangulationType>
1974 const double &simplificationThreshold,
1975 const SimplificationCriterion &simplificationCriterion,
1976 const triangulationType &triangulation) {
1977
1978 Timer t;
1979
1980 if(!currentData_.sheet3List_.size())
1981 return -1;
1982
1983 SimplexId simplifiedSheets = 0;
1984 double lastThreshold = -1;
1985
1986 for(size_t it = 0; it < originalData_.sheet3List_.size(); it++) {
1987
1988 // do while, avoiding infinite loop
1989
1990 double minValue = -1;
1991 SimplexId minId = -1;
1992
1993 for(size_t i = 0; i < currentData_.sheet3List_.size(); i++) {
1994 if(!currentData_.sheet3List_[i].pruned_) {
1995
1996 double value = 0;
1997 switch(simplificationCriterion) {
1998
1999 case SimplificationCriterion::domainVolume:
2000 value = currentData_.sheet3List_[i].domainVolume_ / totalVolume_;
2001 break;
2002
2003 case SimplificationCriterion::rangeArea:
2004 value = currentData_.sheet3List_[i].rangeArea_ / totalArea_;
2005 break;
2006
2007 case SimplificationCriterion::hyperVolume:
2008 value
2009 = currentData_.sheet3List_[i].hyperVolume_ / totalHyperVolume_;
2010 break;
2011 }
2012
2013 if((minId == -1) || (value < minValue)) {
2014 minValue = value;
2015 minId = i;
2016 }
2017 }
2018 }
2019
2020 if((minId != -1) && (minValue < simplificationThreshold)) {
2021 simplifySheet(minId, simplificationCriterion, triangulation);
2022 simplifiedSheets++;
2023 lastThreshold = minValue;
2024 } else {
2025 break;
2026 }
2027 }
2028
2029 currentData_.simplificationThreshold_ = simplificationThreshold;
2030 currentData_.simplificationCriterion_ = simplificationCriterion;
2031
2032 SimplexId simplificationId = 0;
2033 for(size_t i = 0; i < currentData_.sheet3List_.size(); i++) {
2034 if(!currentData_.sheet3List_[i].pruned_) {
2035 currentData_.sheet3List_[i].simplificationId_ = simplificationId;
2036 simplificationId++;
2037 }
2038 }
2039 for(size_t i = 0; i < currentData_.sheet3List_.size(); i++) {
2040 if(currentData_.sheet3List_[i].pruned_) {
2041 // find where it merged
2042 if(currentData_.sheet3List_[i].vertexList_.size()) {
2043 SimplexId const vertexId = currentData_.sheet3List_[i].vertexList_[0];
2044 SimplexId const sheetId = currentData_.vertex2sheet3_[vertexId];
2045 if(sheetId != static_cast<SimplexId>(i)) {
2046 currentData_.sheet3List_[i].simplificationId_
2047 = currentData_.sheet3List_[sheetId].simplificationId_;
2048 }
2049 }
2050 }
2051 }
2052
2053 // orphans
2054 for(size_t i = 0; i < currentData_.sheet1List_.size(); i++) {
2055 if((!currentData_.sheet1List_[i].pruned_)
2056 && (currentData_.sheet1List_[i].hasSaddleEdges_)) {
2057 SimplexId nonSimplified = 0;
2058 for(size_t j = 0; j < currentData_.sheet1List_[i].sheet3List_.size();
2059 j++) {
2060 SimplexId const sheet3Id = currentData_.sheet1List_[i].sheet3List_[j];
2061
2062 if(!currentData_.sheet3List_[sheet3Id].pruned_) {
2063 nonSimplified++;
2064 }
2065 }
2066
2067 if(nonSimplified < 2) {
2068 currentData_.sheet1List_[i].pruned_ = true;
2069 }
2070 }
2071 }
2072
2073 // TODO: update segmentation for 1-sheets and 0-sheets?...
2074
2075 printConnectivity(currentData_);
2076
2077 this->printMsg(std::vector<std::vector<std::string>>{
2078 {"#3-sheets simplified", std::to_string(simplifiedSheets)},
2079 {"Last 3-sheet threshold", std::to_string(lastThreshold)},
2080 {"#3-sheets left", std::to_string(simplificationId)}});
2081
2082 this->printMsg(
2083 "3-sheets simplified", 1.0, t.getElapsedTime(), this->threadNumber_);
2084
2085 return 0;
2086}
AbstractTriangulation is an interface class that defines an interface for efficient traversal methods...
virtual SimplexId getNumberOfCells() const
virtual SimplexId getNumberOfVertices() const
virtual SimplexId getNumberOfEdges() const
Minimalist debugging class.
Definition Debug.h:88
TTK processing package that computes fiber surfaces.
int setInputField(const void *uField, const void *vField)
int finalize(const bool &mergeDuplicatedVertices=false, const bool &removeSmallEdges=false, const bool &edgeFlips=false, const bool &intersectionRemesh=false)
void preconditionTriangulation(AbstractTriangulation *triangulation)
TTK processing package for the computation of the Jacobi set of bivariate volumetric data.
Definition JacobiSet.h:39
void setSosOffsetsV(const SimplexId *const sosOffsets)
Definition JacobiSet.h:103
void preconditionTriangulation(AbstractTriangulation *const triangulation)
Definition JacobiSet.h:118
void setSosOffsetsU(const SimplexId *const sosOffsets)
Definition JacobiSet.h:91
int execute(std::vector< std::pair< SimplexId, char > > &jacobiSet, const dataTypeU *const uField, const dataTypeV *const vField, const triangulationType &triangulation, std::vector< char > *isPareto=nullptr)
TTK processing package that efficiently computes the Reeb space of bivariate volumetric data.
Definition ReebSpace.h:42
void setSosOffsetsV(const SimplexId *const sosOffsetsV)
Definition ReebSpace.h:232
int disconnect3sheetFrom2sheet(ReebSpaceData &data, const SimplexId &sheet3Id, const SimplexId &sheet2Id)
std::vector< SimplexId > jacobi2edges_
Definition ReebSpace.h:427
SimplexId vertexNumber_
Definition ReebSpace.h:416
const std::vector< SimplexId > * get0sheetSegmentation() const
Definition ReebSpace.h:147
std::vector< std::pair< SimplexId, char > > jacobiSetEdges_
Definition ReebSpace.h:426
bool setRangeDrivenOctree(const bool &onOff)
Definition ReebSpace.h:202
bool expand3sheets_
Definition ReebSpace.h:421
const std::vector< char > * getEdgeTypes() const
Definition ReebSpace.h:163
int compute3sheets(std::vector< std::vector< std::array< SimplexId, 3 > > > &tetTriangles, const triangulationType &triangulation)
Definition ReebSpace.h:1319
const SimplexId * sosOffsetsU_
Definition ReebSpace.h:419
int connectSheets(const triangulationType &triangulation)
Definition ReebSpace.h:1575
int getJacobi2Edge(const SimplexId &jacobiEdgeId) const
Definition ReebSpace.h:172
ReebSpaceData originalData_
Definition ReebSpace.h:423
double totalHyperVolume_
Definition ReebSpace.h:417
ReebSpaceData currentData_
Definition ReebSpace.h:423
int preconditionTriangulation(AbstractTriangulation *const triangulation)
Definition ReebSpace.h:250
FiberSurface fiberSurface_
Definition ReebSpace.h:429
bool empty() const
Definition ReebSpace.h:96
const Sheet3 * get3sheet(const SimplexId &sheetId) const
Definition ReebSpace.h:137
int connect3sheetTo3sheet(ReebSpaceData &data, const SimplexId &sheet3Id, const SimplexId &otherSheet3Id)
Definition ReebSpace.cpp:82
SimplexId edgeNumber_
Definition ReebSpace.h:416
void setExpand3Sheets(const bool &onOff)
Definition ReebSpace.h:198
std::vector< FiberSurface::Vertex > fiberSurfaceVertexList_
Definition ReebSpace.h:431
int simplifySheets(const double &simplificationThreshold, const SimplificationCriterion &simplificationCriterion, const triangulationType &triangulation)
Definition ReebSpace.h:1973
int connect3sheetTo0sheet(ReebSpaceData &data, const SimplexId &sheet3Id, const SimplexId &sheet0Id)
Definition ReebSpace.cpp:7
int printConnectivity(const ReebSpaceData &data) const
int connect3sheetTo2sheet(ReebSpaceData &data, const SimplexId &sheet3Id, const SimplexId &sheet2Id)
Definition ReebSpace.cpp:57
int mergeSheets(const SimplexId &smallerId, const SimplexId &biggerId, const triangulationType &triangulation)
Definition ReebSpace.h:1740
int simplify(const dataTypeU *const uField, const dataTypeV *const vField, const triangulationType &triangulation, const double &simplificationThreshold, const SimplificationCriterion &criterion=SimplificationCriterion::rangeArea)
Definition ReebSpace.h:807
int compute2sheetChambers(const dataTypeU *const uField, const dataTypeV *const vField, const triangulationType &triangulation)
Definition ReebSpace.h:673
int prepareSimplification()
int connect3sheetTo1sheet(ReebSpaceData &data, const SimplexId &sheet3Id, const SimplexId &sheet1Id)
Definition ReebSpace.cpp:32
const std::vector< SimplexId > * get3sheetTetSegmentation() const
Definition ReebSpace.h:159
const std::vector< SimplexId > * get3sheetVertexSegmentation() const
Definition ReebSpace.h:155
const std::vector< SimplexId > * get1sheetSegmentation() const
Definition ReebSpace.h:151
bool withRangeDrivenOctree_
Definition ReebSpace.h:422
const Sheet0 * get0sheet(const SimplexId &sheetId) const
Definition ReebSpace.h:105
int preMergeSheets(const SimplexId &sheetId0, const SimplexId &sheetId1)
double totalArea_
Definition ReebSpace.h:417
const Sheet2 * get2sheet(const SimplexId &sheetId) const
Definition ReebSpace.h:127
bool hasConnectedSheets_
Definition ReebSpace.h:421
int disconnect1sheetFrom0sheet(ReebSpaceData &data, const SimplexId &sheet1Id, const SimplexId &sheet0Id, const SimplexId &biggerId)
int execute(const dataTypeU *const uField, const dataTypeV *const vField, const triangulationType &triangulation)
Definition ReebSpace.h:439
SimplexId tetNumber_
Definition ReebSpace.h:416
void setTetNumber(const SimplexId &tetNumber)
Definition ReebSpace.h:236
int compute2sheets(const std::vector< std::pair< SimplexId, SimplexId > > &jacobiEdges, const dataTypeU *const uField, const dataTypeV *const vField, const triangulationType &triangulation)
Definition ReebSpace.h:507
void setSosOffsetsU(const SimplexId *const sosOffsetsU)
Definition ReebSpace.h:220
JacobiSet jacobiSet_
Definition ReebSpace.h:430
const Sheet1 * get1sheet(const SimplexId &sheetId) const
Definition ReebSpace.h:116
const std::vector< FiberSurface::Vertex > * getFiberSurfaceVertices() const
Definition ReebSpace.h:168
double totalVolume_
Definition ReebSpace.h:417
int disconnect3sheetFrom3sheet(ReebSpaceData &data, const SimplexId &sheet3Id, const SimplexId &other3SheetId)
const SimplexId * sosOffsetsV_
Definition ReebSpace.h:419
int compute1sheets(const std::vector< std::pair< SimplexId, char > > &jacobiSet, std::vector< std::pair< SimplexId, SimplexId > > &jacobiSetClassification, const triangulationType &triangulation)
Definition ReebSpace.h:1008
int computeGeometricalMeasures(Sheet3 &sheet, const dataTypeU *const uField, const dataTypeV *const vField, const triangulationType &triangulation)
Definition ReebSpace.h:744
int flush(const triangulationType &triangulation)
Definition ReebSpace.h:1705
int compute1sheetsOnly(const std::vector< std::pair< SimplexId, char > > &jacobiSet, std::vector< std::pair< SimplexId, SimplexId > > &jacobiSetClassification, const triangulationType &triangulation)
Definition ReebSpace.h:869
int disconnect3sheetFrom1sheet(ReebSpaceData &data, const SimplexId &sheet3Id, const SimplexId &sheet1Id, const SimplexId &biggerId, const triangulationType &triangulation)
Definition ReebSpace.h:1654
int compute3sheet(const SimplexId &vertexId, const std::vector< std::vector< std::array< SimplexId, 3 > > > &tetTriangles, const triangulationType &triangulation)
Definition ReebSpace.h:1208
void setVertexNumber(const SimplexId &vertexNumber)
Definition ReebSpace.h:243
int simplifySheet(const SimplexId &sheetId, const SimplificationCriterion &simplificationCriterion, const triangulationType &triangulation)
Definition ReebSpace.h:1848
int disconnect3sheetFrom0sheet(ReebSpaceData &data, const SimplexId &sheet3Id, const SimplexId &sheet0Id)
SimplexId getNumberOf2sheets() const
Definition ReebSpace.h:179
int perturb(const dataTypeU *const uField, const dataTypeV *const vField, const dataTypeU &uEpsilon=Geometry::powIntTen(-DBL_DIG), const dataTypeV &vEpsilon=Geometry::powIntTen(-DBL_DIG))
Definition ReebSpace.h:795
double getElapsedTime()
Definition Timer.h:15
int getBoundingBox(const Container &points, std::array< std::pair< T, T >, dim > &bBox)
Definition Geometry.h:273
T powIntTen(const int n)
Compute the nth power of ten.
Definition Geometry.h:405
The Topology ToolKit.
int ThreadId
Identifier type for threads (i.e. with OpenMP).
Definition DataTypes.h:26
int SimplexId
Identifier type for simplices of any dimension.
Definition DataTypes.h:22
std::array< SimplexId, 3 > vertexIds_
std::pair< SimplexId, SimplexId > meshEdge_
std::vector< SimplexId > sheet1List_
Definition ReebSpace.h:54
std::vector< SimplexId > sheet3List_
Definition ReebSpace.h:55
std::vector< SimplexId > sheet3List_
Definition ReebSpace.h:66
std::vector< SimplexId > edgeList_
Definition ReebSpace.h:60
std::vector< SimplexId > sheet0List_
Definition ReebSpace.h:63
std::vector< std::vector< FiberSurface::Vertex > > vertexList_
Definition ReebSpace.h:76
std::vector< std::vector< FiberSurface::Triangle > > triangleList_
Definition ReebSpace.h:77
std::vector< SimplexId > sheet3List_
Definition ReebSpace.h:78
std::vector< SimplexId > sheet2List_
Definition ReebSpace.h:89
std::vector< SimplexId > sheet3List_
Definition ReebSpace.h:90
std::vector< SimplexId > tetList_
Definition ReebSpace.h:86
std::vector< SimplexId > sheet0List_
Definition ReebSpace.h:87
SimplexId simplificationId_
Definition ReebSpace.h:82
std::vector< SimplexId > vertexList_
Definition ReebSpace.h:85
std::vector< SimplexId > preMergedSheets_
Definition ReebSpace.h:91
std::vector< SimplexId > sheet1List_
Definition ReebSpace.h:88
printMsg(debug::output::BOLD+" | | | | | . \\ | | (__| | / __/| |_| / __/|__ _|"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)