TTK
Loading...
Searching...
No Matches
ttkContourForests.cpp
Go to the documentation of this file.
1#include <vtkAppendPolyData.h>
2#include <vtkCellArray.h>
3#include <vtkCellData.h>
4#include <vtkDataSet.h>
5#include <vtkDoubleArray.h>
6#include <vtkGenericCell.h>
7#include <vtkInformation.h>
8#include <vtkIntArray.h>
9#include <vtkLine.h>
10#include <vtkLineSource.h>
11#include <vtkNew.h>
12#include <vtkObjectFactory.h>
13#include <vtkPointData.h>
14
15#include "ttkContourForests.h"
16
17#include <ttkMacros.h>
18#include <ttkUtils.h>
19
20using namespace std;
21using namespace ttk;
22using namespace cf;
23
25
27 // VTK Interface //
28 this->SetNumberOfInputPorts(1);
29 this->SetNumberOfOutputPorts(3);
30
31 vtkWarningMacro(
32 "Contour Forests is deprecated, please use FTM Tree instead.");
33}
34
36 samples_.clear();
37 barycenters_.clear();
38
39 skeletonNodes_->Delete();
40 skeletonNodes_ = vtkPolyData::New();
41 skeletonArcs_->Delete();
42 skeletonArcs_ = vtkPolyData::New();
43}
44
46 segmentation_->Delete();
47}
48
50 tree_ = nullptr;
51}
52
54 vtkInformation *info) {
55 if(port == 0) {
56 info->Set(vtkDataObject::DATA_TYPE_NAME(), "vtkDataSet");
57 return 1;
58 }
59 return 0;
60}
61
63 vtkInformation *info) {
64 if(port == 0 || port == 1) {
65 info->Set(vtkDataObject::DATA_TYPE_NAME(), "vtkPolyData");
66 return 1;
67 } else if(port == 2) {
69 return 1;
70 }
71 return 0;
72}
73
75 toComputeSkeleton_ = true;
76 toComputeContourTree_ = true;
77 ttkAlgorithm::Modified();
78}
79
81 toUpdateVertexSoSoffsets_ = true;
82 toComputeContourTree_ = true;
83 toUpdateTree_ = true;
84 toComputeSkeleton_ = true;
85 toComputeSegmentation_ = true;
86
87 ForceInputOffsetScalarField = onOff;
88 Modified();
89}
90
92 if(treeType >= 0 && treeType <= 2) {
93 toUpdateTree_ = true;
94 if(treeType_ == TreeType::Contour
95 || static_cast<TreeType>(treeType) == TreeType::Contour) {
96 toComputeContourTree_ = true;
97 }
98
99 toComputeSkeleton_ = true;
100 toComputeSegmentation_ = true;
101
102 treeType_ = static_cast<TreeType>(treeType);
103
104 Modified();
105 }
106}
107
109 toComputeSkeleton_ = true;
110
111 showMin_ = state;
112 Modified();
113}
114
116 toComputeSkeleton_ = true;
117
118 showMax_ = state;
119 Modified();
120}
121
123 toComputeSkeleton_ = true;
124
125 showSaddle1_ = state;
126 Modified();
127}
128
130 toComputeSkeleton_ = true;
131
132 showSaddle2_ = state;
133 Modified();
134}
135
137 toComputeSkeleton_ = true;
138
139 showArc_ = state;
140 Modified();
141}
142
143void ttkContourForests::SetArcResolution(int arcResolution) {
144 if(arcResolution >= 0) {
145 toComputeSkeleton_ = true;
146
147 arcResolution_ = arcResolution;
148 Modified();
149 }
150}
151
153 partitionNum_ = partitionNum;
154
155 toComputeContourTree_ = true;
156 toComputeSkeleton_ = true;
157 toUpdateTree_ = true;
158 Modified();
159}
160
162 lessPartition_ = l;
163 Modified();
164}
165
166void ttkContourForests::SetSkeletonSmoothing(double skeletonSmoothing) {
167 if(skeletonSmoothing >= 0) {
168 toComputeSkeleton_ = true;
169
170 skeletonSmoothing_ = skeletonSmoothing;
171 Modified();
172 }
173}
174
176 simplificationType_ = type;
177 Modified();
178}
179
181 double simplificationThreshold) {
182 if(simplificationThreshold >= 0.0 && simplificationThreshold <= 1.0) {
183 simplificationThresholdBuffer_ = simplificationThreshold;
184 toComputeContourTree_ = true;
185 toUpdateTree_ = true;
186 toComputeSkeleton_ = true;
187 toComputeSegmentation_ = true;
188
189 Modified();
190 }
191}
192
193bool ttkContourForests::isCoincident(float p1[], double p2[]) {
194 double sPrev[3];
195 double sNext[3];
196 for(unsigned int k = 0; k < 3; k++) {
197 sPrev[k] = p2[k] - p1[k];
198 sNext[k] = sPrev[k];
199 }
200
201 return (vtkMath::Normalize(sNext) == 0.0);
202}
203
204bool ttkContourForests::isCoincident(double p1[], double p2[]) {
205 double sPrev[3];
206 double sNext[3];
207 for(unsigned int k = 0; k < 3; k++) {
208 sPrev[k] = p2[k] - p1[k];
209 sNext[k] = sPrev[k];
210 }
211
212 return (vtkMath::Normalize(sNext) == 0.0);
213}
214
216 vtkNew<vtkAppendPolyData> app{};
217
218 vtkDoubleArray *scalars{};
219 ttkSimplexIdTypeArray *identifierScalars{};
220 vtkIntArray *typeScalars{};
221 ttkSimplexIdTypeArray *sizeScalars{};
222 vtkDoubleArray *spanScalars{};
223 int type = static_cast<int>(TreeComponent::Arc);
224
225 float point1[3];
226 vector<double> point2(3);
227 // get skeleton scalars
228 vector<vector<double>> skeletonScalars{};
229 getSkeletonScalars(vertexScalars_, skeletonScalars);
230
231 double inputScalar;
232 SuperArc *a;
233 SimplexId regionSize;
234 double regionSpan;
235 SimplexId currentZone = 0;
236 SimplexId regionId;
237
238 for(SimplexId i = 0; i < (SimplexId)tree_->getNumberOfSuperArcs(); ++i) {
239 a = tree_->getSuperArc(i);
240
241 if(a->isVisible()) {
242 SimplexId upNodeId = tree_->getSuperArc(i)->getUpNodeId();
243 SimplexId upVertex = tree_->getNode(upNodeId)->getVertexId();
244 float coordUp[3];
245 triangulation_->getVertexPoint(
246 upVertex, coordUp[0], coordUp[1], coordUp[2]);
247
248 SimplexId downNodeId = tree_->getSuperArc(i)->getDownNodeId();
249 SimplexId downVertex = tree_->getNode(downNodeId)->getVertexId();
250 float coordDown[3];
251 triangulation_->getVertexPoint(
252 downVertex, coordDown[0], coordDown[1], coordDown[2]);
253
254 regionSize = tree_->getNumberOfVisibleRegularNode(i);
255 regionSpan = Geometry::distance(coordUp, coordDown, 3);
256 regionId = currentZone++;
257
258 // Line //
259 if(barycenters_[static_cast<int>(treeType_)][i].size()) {
260 // init: min
261 SimplexId downNodeVId;
262 if(treeType_ == TreeType::Split)
263 downNodeVId = tree_->getNode(a->getUpNodeId())->getVertexId();
264 else
265 downNodeVId = tree_->getNode(a->getDownNodeId())->getVertexId();
266
267 triangulation_->getVertexPoint(
268 downNodeVId, point1[0], point1[1], point1[2]);
271 line->SetPoint1(point1);
272
273 const auto nbBarycenter
274 = barycenters_[static_cast<int>(treeType_)][i].size();
275 for(unsigned int j = 0; j < nbBarycenter; ++j) {
276 point2 = barycenters_[static_cast<int>(treeType_)][i][j];
277 line->SetPoint2(point2.data());
278
279 if(!isCoincident(point1, point2.data())) {
280 line->Update();
281 vtkPolyData *lineData = line->GetOutput();
282
283 // Point data //
284 {
285 inputScalar = skeletonScalars[i][j];
286
287 scalars = vtkDoubleArray::New();
288 scalars->SetName(vtkInputScalars_->GetName());
289 for(unsigned int k = 0; k < 2; ++k)
290 scalars->InsertTuple1(k, inputScalar);
291 lineData->GetPointData()->AddArray(scalars);
292 scalars->Delete();
293 }
294
295 // Cell data //
296 // Identifier
297 identifierScalars = ttkSimplexIdTypeArray::New();
298 identifierScalars->SetName("SegmentationId");
299 for(unsigned int k = 0; k < 2; ++k)
300 identifierScalars->InsertTuple1(k, regionId);
301 lineData->GetCellData()->AddArray(identifierScalars);
302 identifierScalars->Delete();
303 // Type
304 typeScalars = vtkIntArray::New();
305 typeScalars->SetName("Type");
306 for(unsigned int k = 0; k < 2; ++k)
307 typeScalars->InsertTuple1(k, type);
308 lineData->GetCellData()->AddArray(typeScalars);
309 typeScalars->Delete();
310 // Size
311 sizeScalars = ttkSimplexIdTypeArray::New();
312 sizeScalars->SetName("RegionSize");
313 for(unsigned int k = 0; k < 2; ++k)
314 sizeScalars->InsertTuple1(k, regionSize);
315 lineData->GetCellData()->AddArray(sizeScalars);
316 sizeScalars->Delete();
317 // Span
318 spanScalars = vtkDoubleArray::New();
319 spanScalars->SetName("RegionSpan");
320 for(unsigned int k = 0; k < 2; ++k)
321 spanScalars->InsertTuple1(k, regionSpan);
322 lineData->GetCellData()->AddArray(spanScalars);
323 spanScalars->Delete();
324
325 app->AddInputData(lineData);
326 }
327
329 line->SetPoint1(point2.data());
330 }
331
332 // end: max
333 SimplexId upNodeVId;
334 if(treeType_ == TreeType::Split)
335 upNodeVId = tree_->getNode(a->getDownNodeId())->getVertexId();
336 else
337 upNodeVId = tree_->getNode(a->getUpNodeId())->getVertexId();
338
339 std::array<float, 3> pt{};
340 triangulation_->getVertexPoint(upNodeVId, pt[0], pt[1], pt[2]);
341 point2[0] = pt[0];
342 point2[1] = pt[1];
343 point2[2] = pt[2];
344 line->SetPoint2(point2.data());
345
346 if(!isCoincident(point1, point2.data())) {
347 line->Update();
348 vtkPolyData *lineData = line->GetOutput();
349
350 // Point data //
351 {
352 inputScalar
353 = skeletonScalars[i][barycenters_[static_cast<int>(treeType_)][i]
354 .size()];
355
356 scalars = vtkDoubleArray::New();
357 scalars->SetName(vtkInputScalars_->GetName());
358 for(unsigned int k = 0; k < 2; ++k)
359 scalars->InsertTuple1(k, inputScalar);
360 lineData->GetPointData()->AddArray(scalars);
361 scalars->Delete();
362 }
363
364 // Cell data //
365 // Identifier
366 identifierScalars = ttkSimplexIdTypeArray::New();
367 identifierScalars->SetName("SegmentationId");
368 for(unsigned int k = 0; k < 2; ++k)
369 identifierScalars->InsertTuple1(k, regionId);
370 lineData->GetCellData()->AddArray(identifierScalars);
371 identifierScalars->Delete();
372 // Type
373 typeScalars = vtkIntArray::New();
374 typeScalars->SetName("Type");
375 for(unsigned int k = 0; k < 2; ++k)
376 typeScalars->InsertTuple1(k, type);
377 lineData->GetCellData()->AddArray(typeScalars);
378 typeScalars->Delete();
379 // Size
380 sizeScalars = ttkSimplexIdTypeArray::New();
381 sizeScalars->SetName("RegionSize");
382 for(unsigned int k = 0; k < 2; ++k)
383 sizeScalars->InsertTuple1(k, regionSize);
384 lineData->GetCellData()->AddArray(sizeScalars);
385 sizeScalars->Delete();
386 // Span
387 spanScalars = vtkDoubleArray::New();
388 spanScalars->SetName("RegionSpan");
389 for(unsigned int k = 0; k < 2; ++k)
390 spanScalars->InsertTuple1(k, regionSpan);
391 lineData->GetCellData()->AddArray(spanScalars);
392 spanScalars->Delete();
393
394 app->AddInputData(lineData);
395 }
396 } else {
397 vtkNew<vtkLineSource> line{};
398
399 SimplexId downNodeVId
400 = tree_->getNode(a->getDownNodeId())->getVertexId();
401 triangulation_->getVertexPoint(
402 downNodeVId, point1[0], point1[1], point1[2]);
403 line->SetPoint1(point1);
404
405 SimplexId upNodeVId = tree_->getNode(a->getUpNodeId())->getVertexId();
406 std::array<float, 3> pt{};
407 triangulation_->getVertexPoint(upNodeVId, pt[0], pt[1], pt[2]);
408 point2[0] = pt[0];
409 point2[1] = pt[1];
410 point2[2] = pt[2];
411 line->SetPoint2(point2.data());
412
413 if(!isCoincident(point1, point2.data())) {
414 line->Update();
415 vtkPolyData *lineData = line->GetOutput();
416
417 // Point data //
418 {
419 inputScalar = skeletonScalars[i][0];
420
421 scalars = vtkDoubleArray::New();
422 scalars->SetName(vtkInputScalars_->GetName());
423 for(unsigned int k = 0; k < 2; ++k)
424 scalars->InsertTuple1(k, inputScalar);
425 lineData->GetPointData()->AddArray(scalars);
426 scalars->Delete();
427 }
428
429 // Cell data //
430 // Identifier
431 identifierScalars = ttkSimplexIdTypeArray::New();
432 identifierScalars->SetName("SegmentationId");
433 for(int k = 0; k < 2; ++k)
434 identifierScalars->InsertTuple1(k, regionId);
435 lineData->GetCellData()->AddArray(identifierScalars);
436 identifierScalars->Delete();
437 // Type
438 typeScalars = vtkIntArray::New();
439 typeScalars->SetName("Type");
440 for(unsigned int k = 0; k < 2; ++k)
441 typeScalars->InsertTuple1(k, type);
442 lineData->GetCellData()->AddArray(typeScalars);
443 typeScalars->Delete();
444 // Size
445 sizeScalars = ttkSimplexIdTypeArray::New();
446 sizeScalars->SetName("RegionSize");
447 for(unsigned int k = 0; k < 2; ++k)
448 sizeScalars->InsertTuple1(k, regionSize);
449 lineData->GetCellData()->AddArray(sizeScalars);
450 sizeScalars->Delete();
451 // Span
452 spanScalars = vtkDoubleArray::New();
453 spanScalars->SetName("RegionSpan");
454 for(unsigned int k = 0; k < 2; ++k)
455 spanScalars->InsertTuple1(k, regionSpan);
456 lineData->GetCellData()->AddArray(spanScalars);
457 spanScalars->Delete();
458
459 app->AddInputData(lineData);
460 }
461 }
462 } else {
463 // cout << " pruned _ :" <<
464 // tree_->getNode(a->getDownNodeId())->getVertexId() << " - "
465 //<< tree_->getNode(a->getUpNodeId())->getVertexId() << endl;
466 }
467 }
468
469 app->Update();
470 skeletonArcs_->ShallowCopy(app->GetOutput());
471}
472
474 const vector<double> &scalars,
475 vector<vector<double>> &skeletonScalars) const {
476 skeletonScalars.clear();
477 skeletonScalars.resize(tree_->getNumberOfSuperArcs());
478
479 SimplexId nodeId;
480 SimplexId vertexId;
481
482 double f;
483 double f0;
484 double f1;
485 double fmin;
486 double fmax;
487 SimplexId nodeMinId;
488 SimplexId nodeMaxId;
489 SimplexId nodeMinVId;
490 SimplexId nodeMaxVId;
491 const SuperArc *a;
492 for(SimplexId i = 0; i < (SimplexId)tree_->getNumberOfSuperArcs(); ++i) {
493 a = tree_->getSuperArc(i);
494
495 if(!a->isPruned()) {
496 if(treeType_ == TreeType::Split) {
497 nodeMinId = a->getUpNodeId();
498 nodeMaxId = a->getDownNodeId();
499 } else {
500 nodeMaxId = a->getUpNodeId();
501 nodeMinId = a->getDownNodeId();
502 }
503
504 nodeMaxVId = tree_->getNode(nodeMaxId)->getVertexId();
505 nodeMinVId = tree_->getNode(nodeMinId)->getVertexId();
506
507 fmax = scalars[nodeMaxVId];
508 fmin = scalars[nodeMinVId];
509
510 // init: min
511 f0 = fmin;
512
513 // iteration
514 for(SimplexId j = 0;
515 j < (SimplexId)samples_[static_cast<int>(treeType_)][i].size(); ++j) {
516 const vector<SimplexId> &sample
517 = samples_[static_cast<int>(treeType_)][i][j];
518
519 f = 0;
520 for(SimplexId k = 0; k < (SimplexId)sample.size(); ++k) {
521 nodeId = sample[k];
522 vertexId = nodeId;
523 f += scalars[vertexId];
524 }
525 if(sample.size()) {
526 f /= sample.size();
527
528 f1 = f;
529 // update the arc
530 skeletonScalars[i].push_back((f0 + f1) / 2);
531 f0 = f1;
532 }
533 }
534
535 // end: max
536 f1 = fmax;
537
538 // update the arc
539 skeletonScalars[i].push_back((f0 + f1) / 2);
540 }
541 }
542
543 return 0;
544}
545
547 vtkNew<vtkPoints> points{};
548 float point[3];
549
550 double scalar{};
551 vtkNew<vtkDoubleArray> scalars{};
552 scalars->SetName(vtkInputScalars_->GetName());
553
554 vtkNew<ttkSimplexIdTypeArray> nodeIdentifierScalars{};
555 nodeIdentifierScalars->SetName("NodeIdentifier");
556
557 vtkNew<ttkSimplexIdTypeArray> vertexIdentifierScalars{};
558 vertexIdentifierScalars->SetName(ttk::VertexScalarFieldName);
559
560 int type{};
561 vtkNew<vtkIntArray> nodeTypeScalars{};
562 nodeTypeScalars->SetName("CriticalType");
563
564 vtkNew<ttkSimplexIdTypeArray> regionSizeScalars{};
565 regionSizeScalars->SetName("RegionSize");
566
567 SimplexId identifier{};
568 for(unsigned i = 0; i < criticalPoints_.size(); ++i) {
569 SimplexId nodeId = criticalPoints_[i];
570 if(tree_->getNode(nodeId)->isHidden())
571 continue;
572 SimplexId vertexId = tree_->getNode(nodeId)->getVertexId();
573 CriticalType nodeType = getNodeType(nodeId);
574
575 if((nodeType == CriticalType::Local_minimum and showMin_)
576 or (nodeType == CriticalType::Local_maximum and showMax_)
577 or (nodeType == CriticalType::Saddle1 and showSaddle1_)
578 or (nodeType == CriticalType::Saddle2 and showSaddle2_)
579 or ((showSaddle1_ and showSaddle2_)
580 and (nodeType == CriticalType::Regular
581 or nodeType == CriticalType::Degenerate))) {
582 // Positions
583 triangulation_->getVertexPoint(vertexId, point[0], point[1], point[2]);
584 points->InsertPoint(identifier, point);
585
586 // Scalars
587 scalar = vertexScalars_[vertexId];
588 scalars->InsertTuple1(identifier, scalar);
589
590 // NodeIdentifier
591 nodeIdentifierScalars->InsertTuple1(identifier, nodeId);
592
593 // VertexIdentifier
594 vertexIdentifierScalars->InsertTuple1(identifier, vertexId);
595
596 // Type
597 type = static_cast<int>(nodeType);
598 nodeTypeScalars->InsertTuple1(identifier, type);
599
600 // RegionSize
601 SimplexId regionSize = 0;
602 if(nodeType == CriticalType::Local_maximum) {
603 const SimplexId arcId = tree_->getNode(nodeId)->getDownSuperArcId(0);
604 regionSize = tree_->getSuperArc(arcId)->getNumberOfRegularNodes() + 1;
605 } else if(nodeType == CriticalType::Local_minimum) {
606 const SimplexId arcId = tree_->getNode(nodeId)->getUpSuperArcId(0);
607 regionSize = tree_->getSuperArc(arcId)->getNumberOfRegularNodes() + 1;
608 }
609 regionSizeScalars->InsertTuple1(identifier, regionSize);
610
611 ++identifier;
612 }
613 }
614 skeletonNodes_->SetPoints(points);
615 skeletonNodes_->GetPointData()->AddArray(scalars);
616 skeletonNodes_->GetPointData()->AddArray(nodeIdentifierScalars);
617 skeletonNodes_->GetPointData()->AddArray(vertexIdentifierScalars);
618 skeletonNodes_->GetPointData()->AddArray(nodeTypeScalars);
619 skeletonNodes_->GetPointData()->AddArray(regionSizeScalars);
620}
621
623 return getNodeType(id, treeType_, tree_);
624}
625
628 int upDegree{};
629 int downDegree{};
630 if(type == TreeType::Join || type == TreeType::Contour) {
631 upDegree = tree->getNode(id)->getUpValence();
632 downDegree = tree->getNode(id)->getDownValence();
633 } else {
634 downDegree = tree->getNode(id)->getUpValence();
635 upDegree = tree->getNode(id)->getDownValence();
636 }
637 int degree = upDegree + downDegree;
638
639 // saddle point
640 if(degree > 1) {
641 if(upDegree == 2 && downDegree == 1)
642 return CriticalType::Saddle2;
643 else if(upDegree == 1 && downDegree == 2)
644 return CriticalType::Saddle1;
645 else if(upDegree == 1 && downDegree == 1)
646 return CriticalType::Regular;
647 else
648 return CriticalType::Degenerate;
649 }
650 // local extremum
651 else {
652 if(upDegree)
653 return CriticalType::Local_minimum;
654 else
655 return CriticalType::Local_maximum;
656 }
657}
658
660 vector<bool> isCriticalPoint(numberOfVertices_);
661
662 criticalPoints_.clear();
663 for(SimplexId i = 0; i < numberOfVertices_; ++i)
664 isCriticalPoint[i] = false;
665
666 // const int nbVert = triangulation_->getNumberOfOriginalVertices();
667
668 // looking for critical points
669 for(SimplexId i = 0; i < (SimplexId)tree_->getNumberOfSuperArcs(); ++i) {
670 auto a = tree_->getSuperArc(i);
671
672 if(!a->isPruned()) {
673 SimplexId upId = a->getUpNodeId();
674 SimplexId up_vId = tree_->getNode(upId)->getVertexId();
675 if(!isCriticalPoint[up_vId]) {
676 isCriticalPoint[up_vId] = true;
677 criticalPoints_.push_back(upId);
678 }
679
680 SimplexId downId = a->getDownNodeId();
681 SimplexId down_vId = tree_->getNode(downId)->getVertexId();
682 if(!isCriticalPoint[down_vId]) {
683 isCriticalPoint[down_vId] = true;
684 criticalPoints_.push_back(downId);
685 }
686 }
687 }
688 //{
689 // stringstream msg;
690 // msg << "[ttkContourForests] List of critical points :" << endl;
691 // for (unsigned int it = 0; it < criticalPoints_->size(); ++it)
692 // msg << "[ttkContourForests] NodeId:" << (*criticalPoints_)[it]
693 //<< ", VertexId:" << tree_->getNode(it)->getVertexId() << endl;
694 // dMsg(cout, msg.str(), advancedInfoMsg);
695 //}
696}
697
698int ttkContourForests::sample(unsigned int samplingLevel) {
699 samples_.resize(3);
700 samples_[static_cast<int>(treeType_)].resize(tree_->getNumberOfSuperArcs());
701 vector<vector<SimplexId>> sampleList(samplingLevel);
702
703 SuperArc *a;
704 for(SimplexId i = 0; i < (SimplexId)tree_->getNumberOfSuperArcs(); ++i) {
705 a = tree_->getSuperArc(i);
706
707 if(!a->isPruned()) {
708 for(unsigned int j = 0; j < samplingLevel; ++j) {
709 sampleList[j].clear();
710 }
711
712 double fmax, fmin;
713 SimplexId nodeMaxId, nodeMinId;
714 SimplexId nodeMaxVId, nodeMinVId;
715 double delta;
716 if(a->getNumberOfRegularNodes()) {
717 if(treeType_ == TreeType::Split) {
718 nodeMaxId = a->getDownNodeId();
719 nodeMinId = a->getUpNodeId();
720 } else {
721 nodeMaxId = a->getUpNodeId();
722 nodeMinId = a->getDownNodeId();
723 }
724
725 nodeMaxVId = tree_->getNode(nodeMaxId)->getVertexId();
726 nodeMinVId = tree_->getNode(nodeMinId)->getVertexId();
727
728 fmax = vertexScalars_[nodeMaxVId];
729 fmin = vertexScalars_[nodeMinVId];
730
731 delta = (fmax - fmin) / samplingLevel;
732
733 double f;
734 SimplexId nodeId;
735 SimplexId vertexId;
736 for(SimplexId j = 0; j < a->getNumberOfRegularNodes(); ++j) {
737 nodeId = a->getRegularNodeId(j);
738 if(a->isMasqued(j))
739 continue;
740 vertexId = nodeId;
741 f = vertexScalars_[vertexId];
742
743 for(unsigned int k = 0; k < samplingLevel; ++k) {
744 if(f <= (k + 1) * delta + fmin) {
745 sampleList[k].push_back(nodeId);
746 break;
747 }
748 }
749 }
750
751 // update the arc
752 for(SimplexId j = 0; j < (SimplexId)sampleList.size(); ++j)
753 samples_[static_cast<int>(treeType_)][i].push_back(sampleList[j]);
754 }
755 }
756 }
757
758 return 0;
759}
760
762 barycenters_.resize(3);
763 barycenters_[static_cast<int>(treeType_)].resize(
764 tree_->getNumberOfSuperArcs());
765 vector<float> barycenter(3);
766 SimplexId vertexId;
767
768 const SuperArc *a;
769 for(SimplexId i = 0; i < (SimplexId)tree_->getNumberOfSuperArcs(); ++i) {
770 a = tree_->getSuperArc(i);
771 if(!a->isPruned()) {
772 for(SimplexId j = 0;
773 j < (SimplexId)samples_[static_cast<int>(treeType_)][i].size(); ++j) {
774 vector<SimplexId> &sample = samples_[static_cast<int>(treeType_)][i][j];
775
776 for(unsigned int k = 0; k < 3; ++k)
777 barycenter[k] = 0;
778
779 for(SimplexId k = 0; k < (SimplexId)sample.size(); ++k) {
780 vertexId = sample[k];
781
782 float pt[3];
783 triangulation_->getVertexPoint(vertexId, pt[0], pt[1], pt[2]);
784 barycenter[0] += pt[0];
785 barycenter[1] += pt[1];
786 barycenter[2] += pt[2];
787 }
788 if(sample.size()) {
789 for(unsigned int k = 0; k < 3; ++k)
790 barycenter[k] /= sample.size();
791
792 // update the arc
793 unsigned int nbBar
794 = barycenters_[static_cast<int>(treeType_)][i].size();
795 barycenters_[static_cast<int>(treeType_)][i].resize(nbBar + 1);
796 barycenters_[static_cast<int>(treeType_)][i][nbBar].resize(3);
797
798 for(unsigned int k = 0; k < 3; ++k)
799 barycenters_[static_cast<int>(treeType_)][i][nbBar][k]
800 = barycenter[k];
801 }
802 }
803 }
804 }
805
806 return 0;
807}
808
809void ttkContourForests::computeSkeleton(unsigned int arcRes) {
810 sample(arcRes);
812}
813
814void ttkContourForests::smoothSkeleton(unsigned int skeletonSmoothing) {
815 for(unsigned int i = 0; i < skeletonSmoothing; i++) {
816 for(SimplexId j = 0; j < (SimplexId)tree_->getNumberOfSuperArcs(); j++) {
817 if(!tree_->getSuperArc(j)->isPruned()) {
818 smooth(j, !(treeType_ == TreeType::Split));
819 }
820 }
821 }
822}
823
824void ttkContourForests::smooth(const SimplexId idArc, bool order) {
825 int N = barycenters_[static_cast<int>(treeType_)][idArc].size();
826 if(N) {
827 // init //
828 vector<vector<double>> barycenterList(N);
829 for(unsigned int i = 0; i < barycenterList.size(); ++i)
830 barycenterList[i].resize(3);
831
832 SimplexId up_vId;
833 SimplexId down_vId;
834 if(order) {
835 up_vId = tree_->getNode(tree_->getSuperArc(idArc)->getUpNodeId())
836 ->getVertexId();
837 down_vId = tree_->getNode(tree_->getSuperArc(idArc)->getDownNodeId())
838 ->getVertexId();
839 } else {
840 down_vId = tree_->getNode(tree_->getSuperArc(idArc)->getUpNodeId())
841 ->getVertexId();
842 up_vId = tree_->getNode(tree_->getSuperArc(idArc)->getDownNodeId())
843 ->getVertexId();
844 }
845
846 std::array<float, 3> p0{};
847 std::array<float, 3> p1{};
848 triangulation_->getVertexPoint(down_vId, p0[0], p0[1], p0[2]);
849 triangulation_->getVertexPoint(up_vId, p1[0], p1[1], p1[2]);
850
851 // filtering //
852 if(N > 1) {
853 // first
854 for(unsigned int k = 0; k < 3; ++k)
855 barycenterList[0][k]
856 = (p0[k] + barycenters_[static_cast<int>(treeType_)][idArc][1][k])
857 * 0.5;
858
859 // main
860 for(int i = 1; i < N - 1; ++i) {
861 for(unsigned int k = 0; k < 3; ++k)
862 barycenterList[i][k]
863 = (barycenters_[static_cast<int>(treeType_)][idArc][i - 1][k]
864 + barycenters_[static_cast<int>(treeType_)][idArc][i + 1][k])
865 * 0.5;
866 }
867 // last
868 for(unsigned int k = 0; k < 3; ++k)
869 barycenterList[N - 1][k]
870 = (p1[k] + barycenters_[static_cast<int>(treeType_)][idArc][N - 1][k])
871 * 0.5;
872 } else {
873 for(unsigned int k = 0; k < 3; ++k)
874 barycenterList[0][k] = (p0[k] + p1[k]) * 0.5;
875 }
876
877 // copy //
878 for(int i = 0; i < N; ++i) {
879 for(unsigned int k = 0; k < 3; ++k)
880 barycenters_[static_cast<int>(treeType_)][idArc][i][k]
881 = barycenterList[i][k];
882 }
883 }
884}
885
887 Timer t;
888 computeSkeleton(arcResolution_);
889 smoothSkeleton(skeletonSmoothing_);
890
891 // nodes
892 if(showMin_ || showMax_ || showSaddle1_ || showSaddle2_)
894 else
895 skeletonNodes_->ShallowCopy(voidUnstructuredGrid_);
896
897 // arcs
898 if(showArc_)
900 else
901 skeletonArcs_->ShallowCopy(voidPolyData_);
902
903 // what is done is no longer to be done
904 toComputeSkeleton_ = false;
905
906 this->printMsg(
907 "Topological skeleton built", 1.0, t.getElapsedTime(), this->threadNumber_);
908 this->printMsg(std::vector<std::vector<std::string>>{
909 {"Arc resolution", std::to_string(arcResolution_)},
910 {"Smoothing", std::to_string(skeletonSmoothing_)}});
911}
912
913void ttkContourForests::getSegmentation(vtkDataSet *input) {
914 Timer t;
915
916 // field
917 SimplexId regionId{};
918 vtkNew<ttkSimplexIdTypeArray> scalarsRegionId{};
919 scalarsRegionId->SetName("SegmentationId");
920 scalarsRegionId->SetNumberOfTuples(vertexScalars_.size());
921
922 int regionType{};
923 vtkNew<vtkIntArray> scalarsRegionType{};
924 scalarsRegionType->SetName("RegionType");
925 scalarsRegionType->SetNumberOfTuples(vertexScalars_.size());
926
927 SimplexId regionSize{};
928 vtkNew<ttkSimplexIdTypeArray> scalarsRegionSize{};
929 scalarsRegionSize->SetName("RegionSize");
930 scalarsRegionSize->SetNumberOfTuples(vertexScalars_.size());
931
932 double regionSpan{};
933 vtkNew<vtkDoubleArray> scalarsRegionSpan{};
934 scalarsRegionSpan->SetName("RegionSpan");
935 scalarsRegionSpan->SetNumberOfTuples(vertexScalars_.size());
936
937 SimplexId currentZone{};
938
939 if(!segmentation_) {
940 segmentation_ = input->NewInstance();
941 segmentation_->ShallowCopy(input);
942 }
943
944 for(SimplexId i = 0; i < numberOfVertices_; i++) {
945 scalarsRegionId->SetTuple1(i, -1);
946 }
947
948 // nodes
949 for(SimplexId it = 0; it < (SimplexId)criticalPoints_.size(); ++it) {
950 SimplexId nodeId = criticalPoints_[it];
951 SimplexId vertexId = tree_->getNode(nodeId)->getVertexId();
952
953 // RegionType
954 regionType = -1;
955 scalarsRegionType->SetTuple1(vertexId, regionType);
956 }
957
958 // arcs
959 for(SimplexId i = 0; i < (SimplexId)tree_->getNumberOfSuperArcs(); ++i) {
960 auto a = tree_->getSuperArc(i);
961 if(a->isVisible()) {
962 SimplexId upNodeId = tree_->getSuperArc(i)->getUpNodeId();
963 CriticalType upNodeType = getNodeType(upNodeId);
964 SimplexId upVertex = tree_->getNode(upNodeId)->getVertexId();
965 float coordUp[3];
966 triangulation_->getVertexPoint(
967 upVertex, coordUp[0], coordUp[1], coordUp[2]);
968
969 SimplexId downNodeId = tree_->getSuperArc(i)->getDownNodeId();
970 CriticalType downNodeType = getNodeType(downNodeId);
971 SimplexId downVertex = tree_->getNode(downNodeId)->getVertexId();
972 float coordDown[3];
973 triangulation_->getVertexPoint(
974 downVertex, coordDown[0], coordDown[1], coordDown[2]);
975
976 regionSize = tree_->getNumberOfVisibleRegularNode(i);
977 regionSpan = Geometry::distance(coordUp, coordDown);
978 regionId = currentZone++;
979 // regionId = i;
980
981 // cout << "arc : " << tree_->printArc(i);
982 // cout << " span : " << regionSpan;
983 // cout << " coords : ";
984 // cout << coordDown[0] << ",";
985 // cout << coordDown[1] << ",";
986 // cout << coordDown[2] << " || ";
987 // cout << coordUp[0] << ",";
988 // cout << coordUp[1] << ",";
989 // cout << coordUp[2] << endl;
990
991 scalarsRegionId->SetTuple1(
992 tree_->getNode(downNodeId)->getVertexId(), regionId);
993 scalarsRegionId->SetTuple1(
994 tree_->getNode(upNodeId)->getVertexId(), regionId);
995
996 scalarsRegionSize->SetTuple1(
997 tree_->getNode(downNodeId)->getVertexId(), regionSize);
998 scalarsRegionSize->SetTuple1(
999 tree_->getNode(upNodeId)->getVertexId(), regionSize);
1000
1001 scalarsRegionSpan->SetTuple1(
1002 tree_->getNode(downNodeId)->getVertexId(), regionSpan);
1003 scalarsRegionSpan->SetTuple1(
1004 tree_->getNode(upNodeId)->getVertexId(), regionSpan);
1005
1006 for(SimplexId j = 0; j < tree_->getSuperArc(i)->getNumberOfRegularNodes();
1007 ++j) {
1008 SimplexId nodeId = tree_->getSuperArc(i)->getRegularNodeId(j);
1009 SimplexId vertexId = nodeId;
1010 // cout << vertexId << ", ";
1011 if(tree_->getSuperArc(i)->isMasqued(j)) {
1012 // cout << vertexId << ", ";
1013 continue;
1014 }
1015
1016 // cout << vertexId << ", ";
1017 scalarsRegionId->SetTuple1(vertexId, regionId);
1018 scalarsRegionSize->SetTuple1(vertexId, regionSize);
1019 scalarsRegionSpan->SetTuple1(vertexId, regionSpan);
1020 }
1021 // cout << endl;
1022
1023 // RegionType
1024 if((upNodeType == CriticalType::Local_minimum
1025 && downNodeType == CriticalType::Local_maximum)
1026 || (upNodeType == CriticalType::Local_minimum
1027 || downNodeType == CriticalType::Local_minimum))
1028 regionType = static_cast<int>(ArcType::Min_arc);
1029 else if(upNodeType == CriticalType::Local_maximum
1030 || downNodeType == CriticalType::Local_maximum)
1031 regionType = static_cast<int>(ArcType::Max_arc);
1032 else if(upNodeType == CriticalType::Saddle1
1033 && downNodeType == CriticalType::Saddle1)
1034 regionType = static_cast<int>(ArcType::Saddle1_arc);
1035 else if(upNodeType == CriticalType::Saddle2
1036 && downNodeType == CriticalType::Saddle2)
1037 regionType = static_cast<int>(ArcType::Saddle2_arc);
1038 else
1039 regionType = static_cast<int>(ArcType::Saddle1_saddle2_arc);
1040
1041 for(SimplexId j = 0; j < tree_->getSuperArc(i)->getNumberOfRegularNodes();
1042 ++j) {
1043 SimplexId nodeId = tree_->getSuperArc(i)->getRegularNodeId(j);
1044 if(tree_->getSuperArc(i)->isMasqued(j)) {
1045 // Ignore masqued ones
1046 continue;
1047 }
1048 SimplexId vertexId = nodeId;
1049 scalarsRegionType->SetTuple1(vertexId, regionType);
1050 }
1051 }
1052 }
1053
1054 // output
1055 segmentation_->GetPointData()->AddArray(scalarsRegionId);
1056 segmentation_->GetPointData()->AddArray(scalarsRegionType);
1057 segmentation_->GetPointData()->AddArray(scalarsRegionSize);
1058 segmentation_->GetPointData()->AddArray(scalarsRegionSpan);
1059
1060 this->printMsg("Topological segmentation built", 1.0, t.getElapsedTime(),
1061 this->threadNumber_);
1062 this->printMsg(std::vector<std::vector<std::string>>{
1063 {"Region type", std::to_string(scalarsRegionType->GetNumberOfTuples())},
1064 {"Segmentation Id", std::to_string(scalarsRegionId->GetNumberOfTuples())}});
1065
1066 toComputeSegmentation_ = false;
1067}
1068
1070 // sequential params
1071 this->preconditionTriangulation(triangulation_);
1072 this->setVertexScalars(ttkUtils::GetVoidPointer(vtkInputScalars_));
1073 this->setVertexSoSoffsets(vertexSoSoffsets_);
1074
1075 this->setTreeType(treeType_);
1076 // parallel params
1077 this->setLessPartition(lessPartition_);
1079 this->setPartitionNum(partitionNum_);
1080 // simplification params
1081 this->setSimplificationMethod(simplificationType_);
1082 this->setSimplificationThreshold(simplificationThreshold_);
1083 // build
1084 ttkVtkTemplateMacro(vtkInputScalars_->GetDataType(),
1085 triangulation_->getType(),
1086 (this->build<VTK_TT, TTK_TT *>(
1087 static_cast<TTK_TT *>(triangulation_->getData()))));
1088
1089 // what is done is no longer to be done
1090 toComputeContourTree_ = false;
1091}
1092
1094 // polymorphic tree
1095 switch(treeType_) {
1096 case TreeType::Join:
1097 tree_ = this->getJoinTree();
1098 break;
1099 case TreeType::Split:
1100 tree_ = this->getSplitTree();
1101 break;
1102 case TreeType::JoinAndSplit:
1103 tree_ = this->getJoinTree();
1104 tree_ = this->getSplitTree();
1105 break;
1106 case TreeType::Contour:
1107 tree_ = this;
1108 break;
1109 }
1110
1112
1113 toUpdateTree_ = false;
1114}
1115
1116int ttkContourForests::RequestData(vtkInformation *ttkNotUsed(request),
1117 vtkInformationVector **inputVector,
1118 vtkInformationVector *outputVector) {
1119 vtkWarningMacro(
1120 "DEPRECATED This plugin will be removed in a future release, please use "
1121 "FTM instead for contour trees and FTR for Reeb graphs.");
1122
1123 const auto input = vtkDataSet::GetData(inputVector[0]);
1124 auto outputSkeletonNodes = vtkPolyData::GetData(outputVector, 0);
1125 auto outputSkeletonArcs = vtkPolyData::GetData(outputVector, 1);
1126 auto outputSegmentation = vtkDataSet::GetData(outputVector, 2);
1127
1128#ifndef TTK_ENABLE_KAMIKAZE
1129 if(!input) {
1130 this->printErr("Input pointer is NULL.");
1131 return -1;
1132 }
1133
1134 if(!outputSkeletonNodes || !outputSkeletonArcs || !outputSegmentation) {
1135 this->printErr("Output pointer is NULL.");
1136 return -1;
1137 }
1138
1139 if(!input->GetNumberOfPoints()) {
1140 this->printErr("Input has no point.");
1141 return -1;
1142 }
1143#endif
1144
1145 triangulation_ = ttkAlgorithm::GetTriangulation(input);
1146
1147#ifndef TTK_ENABLE_KAMIKAZE
1148 if(!triangulation_) {
1149 this->printErr("Input triangulation is NULL.");
1150 return -1;
1151 }
1152#endif
1153
1154 varyingMesh_ = false;
1155 if(triangulation_->isEmpty())
1156 varyingMesh_ = true;
1157
1158 // init
1159 if(varyingMesh_ || !numberOfVertices_) {
1160 numberOfVertices_ = input->GetNumberOfPoints();
1161 }
1162
1163 if(varyingMesh_) {
1164 segmentation_ = input->NewInstance();
1165 if(segmentation_ != nullptr) {
1166 segmentation_->ShallowCopy(input);
1167 }
1168 }
1169
1170#ifndef TTK_ENABLE_KAMIKAZE
1171 if(!input->GetPointData()) {
1172 this->printErr("Input has no point data.");
1173 return -2;
1174 }
1175#endif
1176
1177 // scalars
1178 vtkInputScalars_ = this->GetInputArrayToProcess(0, input);
1179
1180#ifndef TTK_ENABLE_KAMIKAZE
1181 if(!vtkInputScalars_) {
1182 this->printErr("Input scalar is NULL.");
1183 return -2;
1184 }
1185#endif
1186
1187 varyingDataValues_ = (vtkInputScalars_->GetMTime() > GetMTime());
1188 if(input->GetPointData()) {
1189
1190 vertexScalars_.resize(numberOfVertices_);
1191 for(SimplexId j = 0; j < numberOfVertices_; ++j) {
1192 vertexScalars_[j] = vtkInputScalars_->GetTuple1(j);
1193 }
1194 }
1195
1196 auto result
1197 = std::minmax_element(vertexScalars_.begin(), vertexScalars_.end());
1198 double scalarMin = *result.first;
1199 double scalarMax = *result.second;
1200 deltaScalar_ = (scalarMax - scalarMin);
1201
1202 // offsets
1203 if(varyingMesh_ || varyingDataValues_ || vertexSoSoffsets_ == nullptr) {
1204
1205 const auto offsets
1206 = this->GetOrderArray(input, 0, 1, ForceInputOffsetScalarField);
1207
1208 if(offsets != nullptr) {
1209 vertexSoSoffsets_
1210 = static_cast<SimplexId *>(ttkUtils::GetVoidPointer(offsets));
1211 }
1212
1213 toUpdateVertexSoSoffsets_ = false;
1214 }
1215
1216 if(varyingMesh_ || varyingDataValues_ || !isLoaded_) {
1217 this->printMsg("Convenient data storage loaded", debug::Priority::DETAIL);
1218 this->printMsg(
1219 std::vector<std::vector<std::string>>{
1220 {"#Tuples", std::to_string(vertexScalars_.size())},
1221 {"#Vertices", std::to_string(numberOfVertices_)},
1222 {"Min", std::to_string(scalarMin)},
1223 {"Max", std::to_string(scalarMax)},
1224 },
1225 debug::Priority::DETAIL);
1226 }
1227
1228 this->printMsg("Launching computation for field `"
1229 + std::string{vtkInputScalars_->GetName()} + "'...");
1230
1231 isLoaded_ = true;
1232
1233 if(simplificationType_ == 0) {
1234 simplificationThreshold_ = simplificationThresholdBuffer_ * deltaScalar_;
1235 } else if(simplificationType_ == 1) {
1236 double coord0[3], coord1[3], spanTotal;
1237 double *bounds = input->GetBounds();
1238 coord0[0] = bounds[0];
1239 coord1[0] = bounds[1];
1240 coord0[1] = bounds[2];
1241 coord1[1] = bounds[3];
1242 coord0[2] = bounds[4];
1243 coord1[2] = bounds[5];
1244 spanTotal = Geometry::distance(coord0, coord1);
1245 simplificationThreshold_ = simplificationThresholdBuffer_ * spanTotal;
1246 } else if(simplificationType_ == 2) {
1247 simplificationThreshold_
1248 = simplificationThresholdBuffer_ * triangulation_->getNumberOfVertices();
1249 }
1250
1251 // ContourForestsTree //
1252 if(varyingMesh_ || varyingDataValues_ || toComputeContourTree_) {
1253 clearTree();
1254 getTree();
1255 updateTree();
1256 }
1257
1258 // Skeleton //
1259 if(varyingMesh_ || varyingDataValues_ || toComputeSkeleton_) {
1260#ifndef TTK_ENABLE_KAMIKAZE
1261 if(tree_ == nullptr) {
1262 this->printErr("MergeTree pointer is NULL.");
1263 return -2;
1264 }
1265#endif // TTK_ENABLE_KAMIKAZE
1266 clearSkeleton();
1267 getSkeleton();
1268 }
1269
1270 // Segmentation //
1271 if(varyingMesh_ || varyingDataValues_ || toComputeSegmentation_)
1272 getSegmentation(input);
1273
1274 // Output //
1275
1276 // skeleton
1277 outputSkeletonNodes->ShallowCopy(skeletonNodes_);
1278 outputSkeletonArcs->ShallowCopy(skeletonArcs_);
1279 // segmentation
1280 outputSegmentation->ShallowCopy(segmentation_);
1281
1282 return 1;
1283}
#define ttkNotUsed(x)
Mark function/method parameters that are not used in the function body at all.
Definition: BaseClass.h:47
static vtkInformationIntegerKey * SAME_DATA_TYPE_AS_INPUT_PORT()
vtkDataArray * GetOrderArray(vtkDataSet *const inputData, const int scalarArrayIdx, const int orderArrayIdx=0, const bool enforceOrderArrayIdx=false)
ttk::Triangulation * GetTriangulation(vtkDataSet *dataSet)
TTK VTK-filter that efficiently computes the contour tree of scalar data and more (data segmentation,...
int FillInputPortInformation(int port, vtkInformation *info) override
void SetTreeType(int tree)
int FillOutputPortInformation(int port, vtkInformation *info) override
void ShowSaddle2(bool state)
void ShowArc(bool state)
void smoothSkeleton(unsigned int skeletonSmoothing)
int getSkeletonScalars(const std::vector< double > &scalars, std::vector< std::vector< double > > &skeletonScalars) const
void SetSkeletonSmoothing(double skeletonSmooth)
void Modified() override
void getSegmentation(vtkDataSet *input)
void SetPartitionNumber(int partitionNum)
bool isCoincident(float p1[], double p2[])
void SetForceInputOffsetScalarField(bool onOff)
void ShowMin(bool state)
void SetSimplificationThreshold(double simplificationThreshold)
void ShowSaddle1(bool state)
void SetSimplificationType(int type)
void ShowMax(bool state)
void SetArcResolution(int arcResolution)
int sample(unsigned int samplingLevel)
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
void smooth(const ttk::SimplexId idArc, bool order)
void computeSkeleton(unsigned int arcRes)
void SetLessPartition(bool l)
ttk::CriticalType getNodeType(ttk::SimplexId id)
static void * GetVoidPointer(vtkDataArray *array, vtkIdType start=0)
Definition: ttkUtils.cpp:225
int threadNumber_
Definition: BaseClass.h:95
int printMsg(const std::string &msg, const debug::Priority &priority=debug::Priority::INFO, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cout) const
Definition: Debug.h:118
int printErr(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition: Debug.h:149
double getElapsedTime()
Definition: Timer.h:15
AbstractTriangulation * getData()
Triangulation::Type getType() const
int getVertexPoint(const SimplexId &vertexId, float &x, float &y, float &z) const override
bool isEmpty() const override
SimplexId getNumberOfVertices() const override
void setLessPartition(bool l)
int setThreadNumber(const int nbThread) override
void setTreeType(const int &local_treeType)
Definition: MergeTree.h:115
void preconditionTriangulation(AbstractTriangulation *const m, const bool preproc=true)
Definition: MergeTree.h:133
void setSimplificationMethod(const int &local_simplifyMethod)
Definition: MergeTree.h:119
idSuperArc getNumberOfSuperArcs() const
Definition: MergeTree.h:183
void setSimplificationThreshold(const double &local_simplificationThreshold)
Definition: MergeTree.h:124
SimplexId getNumberOfVisibleRegularNode(const idSuperArc &sa)
Definition: MergeTree.h:212
void setVertexScalars(scalarType *vals)
Definition: MergeTree.h:159
const std::vector< SuperArc > & getSuperArc() const
Definition: MergeTree.h:197
void setVertexSoSoffsets(const SimplexId *const offsets)
Definition: MergeTree.h:175
Node * getNode(const idNode &nodeId)
Definition: MergeTree.h:244
SimplexId getVertexId() const
idSuperArc getDownSuperArcId(const idSuperArc &neighborId) const
idSuperArc getDownValence() const
bool isHidden() const
idSuperArc getUpSuperArcId(const idSuperArc &neighborId) const
idSuperArc getUpValence() const
bool isMasqued(const SimplexId &v) const
const idNode & getUpNodeId() const
SimplexId getNumberOfRegularNodes()
const SimplexId & getRegularNodeId(const SimplexId &idx)
const idNode & getDownNodeId() const
T distance(const T *p0, const T *p1, const int &dimension=3)
Definition: Geometry.cpp:344
The Topology ToolKit.
CriticalType
default value for critical index
Definition: DataTypes.h:80
const char VertexScalarFieldName[]
default name for vertex scalar field
Definition: DataTypes.h:35
int SimplexId
Identifier type for simplices of any dimension.
Definition: DataTypes.h:22
vtkStandardNewMacro(ttkContourForests)
vtkIntArray ttkSimplexIdTypeArray
Definition: ttkMacros.h:11
#define ttkVtkTemplateMacro(dataType, triangulationType, call)
Definition: ttkMacros.h:69