TTK
Loading...
Searching...
No Matches
ttkMorseSmaleComplex.cpp
Go to the documentation of this file.
1#include <ttkMacros.h>
3#include <ttkUtils.h>
4
5#include <vtkCellData.h>
6#include <vtkDataArray.h>
7#include <vtkDataObject.h>
8#include <vtkDataSet.h>
9#include <vtkDoubleArray.h>
10#include <vtkFloatArray.h>
11#include <vtkIdTypeArray.h>
12#include <vtkInformation.h>
13#include <vtkNew.h>
14#include <vtkPointData.h>
15#include <vtkPolyData.h>
16#include <vtkSignedCharArray.h>
17#include <vtkUnsignedCharArray.h>
18
20
22 this->setDebugMsgPrefix("MorseSmaleComplex");
23 SetNumberOfInputPorts(1);
24 SetNumberOfOutputPorts(4);
25}
26
28 vtkInformation *info) {
29 if(port == 0) {
30 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkDataSet");
31 return 1;
32 }
33 return 0;
34}
35
37 vtkInformation *info) {
38 if(port == 0 || port == 1 || port == 2) {
39 info->Set(vtkDataObject::DATA_TYPE_NAME(), "vtkPolyData");
40 return 1;
41 } else if(port == 3) {
43 return 1;
44 }
45 return 0;
46}
47
48template <typename vtkArrayType, typename vectorType>
49void setArray(vtkArrayType &vtkArray, vectorType &vector) {
50 ttkUtils::SetVoidArray(vtkArray, vector.data(), vector.size(), 1);
51}
52
53template <typename scalarType, typename triangulationType>
54int ttkMorseSmaleComplex::dispatch(vtkDataArray *const inputScalars,
55 vtkPolyData *const outputCriticalPoints,
56 vtkPolyData *const outputSeparatrices1,
57 vtkPolyData *const outputSeparatrices2,
58 const SimplexId *const inputOffsets,
59 const triangulationType &triangulation) {
60
61 const int dimensionality = triangulation.getDimensionality();
62 const auto scalars = ttkUtils::GetPointer<scalarType>(inputScalars);
63
64 const int ret = this->execute(
65 criticalPoints_, separatrices1_, separatrices2_, segmentations_, scalars,
66 inputScalars->GetMTime(), inputOffsets, triangulation);
67
68#ifndef TTK_ENABLE_KAMIKAZE
69 if(ret != 0) {
70 this->printErr("MorseSmaleComplex.execute() error");
71 return -1;
72 }
73#endif
74
75 // critical points
76 {
77 vtkNew<vtkPoints> points{};
78 vtkNew<vtkSignedCharArray> cellDimensions{};
79 vtkNew<ttkSimplexIdTypeArray> cellIds{};
80 vtkSmartPointer<vtkDataArray> const cellScalars{
81 inputScalars->NewInstance()};
82 vtkNew<vtkSignedCharArray> isOnBoundary{};
83 vtkNew<ttkSimplexIdTypeArray> PLVertexIdentifiers{};
84 vtkNew<ttkSimplexIdTypeArray> manifoldSizeScalars{};
85 const auto nPoints = criticalPoints_.points_.size();
86
87#ifndef TTK_ENABLE_KAMIKAZE
88 if(!points || !cellDimensions || !cellIds || !cellScalars || !isOnBoundary
89 || !PLVertexIdentifiers || !manifoldSizeScalars) {
90 this->printErr("Critical points vtkDataArray allocation problem.");
91 return -1;
92 }
93#endif
94
95 points->SetNumberOfPoints(nPoints);
96
97 cellDimensions->SetNumberOfComponents(1);
98 cellDimensions->SetName(ttk::MorseSmaleCellDimensionName);
99 setArray(cellDimensions, criticalPoints_.cellDimensions_);
100
101 cellIds->SetNumberOfComponents(1);
102 cellIds->SetName(ttk::MorseSmaleCellIdName);
103 setArray(cellIds, criticalPoints_.cellIds_);
104
105 cellScalars->SetNumberOfComponents(1);
106 cellScalars->SetName(inputScalars->GetName());
107 cellScalars->SetNumberOfTuples(nPoints);
108#ifdef TTK_ENABLE_OPENMP
109#pragma omp parallel for num_threads(this->threadNumber_)
110#endif // TTK_ENABLE_OPENMP
111 for(size_t i = 0; i < nPoints; ++i) {
112 points->SetPoint(i, criticalPoints_.points_[i].data());
113 cellScalars->SetTuple1(
114 i, scalars[criticalPoints_.PLVertexIdentifiers_[i]]);
115 }
116
117 isOnBoundary->SetNumberOfComponents(1);
118 isOnBoundary->SetName(ttk::MorseSmaleBoundaryName);
119 setArray(isOnBoundary, criticalPoints_.isOnBoundary_);
120
121 PLVertexIdentifiers->SetNumberOfComponents(1);
122 PLVertexIdentifiers->SetName(ttk::VertexScalarFieldName);
123 setArray(PLVertexIdentifiers, criticalPoints_.PLVertexIdentifiers_);
124
125 manifoldSizeScalars->SetNumberOfComponents(1);
126 manifoldSizeScalars->SetName(ttk::MorseSmaleManifoldSizeName);
128 criticalPoints_.manifoldSize_.resize(nPoints);
129 std::fill(criticalPoints_.manifoldSize_.begin(),
130 criticalPoints_.manifoldSize_.end(), -1);
131 }
132 setArray(manifoldSizeScalars, criticalPoints_.manifoldSize_);
133
134 ttkUtils::CellVertexFromPoints(outputCriticalPoints, points);
135
136 auto pointData = outputCriticalPoints->GetPointData();
137#ifndef TTK_ENABLE_KAMIKAZE
138 if(!pointData) {
139 this->printErr("outputCriticalPoints has no point data.");
140 return -1;
141 }
142#endif
143
144 pointData->SetScalars(cellDimensions);
145 pointData->AddArray(cellIds);
146 pointData->AddArray(cellScalars);
147 pointData->AddArray(isOnBoundary);
148 pointData->AddArray(PLVertexIdentifiers);
149 pointData->AddArray(manifoldSizeScalars);
150 }
151
152 // 1-separatrices
155
156 vtkNew<vtkFloatArray> pointsCoords{};
157 vtkNew<vtkSignedCharArray> smoothingMask{};
158 vtkNew<vtkSignedCharArray> cellDimensions{};
159 vtkNew<ttkSimplexIdTypeArray> cellIds{};
160 vtkNew<ttkSimplexIdTypeArray> sourceIds{};
161 vtkNew<ttkSimplexIdTypeArray> destinationIds{};
162 vtkNew<ttkSimplexIdTypeArray> separatrixIds{};
163 vtkNew<vtkSignedCharArray> separatrixTypes{};
164 vtkNew<vtkDoubleArray> separatrixFunctionMaxima{};
165 vtkNew<vtkDoubleArray> separatrixFunctionMinima{};
166 vtkNew<vtkDoubleArray> separatrixFunctionDiffs{};
167 vtkNew<vtkSignedCharArray> isOnBoundary{};
168
169#ifndef TTK_ENABLE_KAMIKAZE
170 if(!pointsCoords || !smoothingMask || !cellDimensions || !cellIds
171 || !sourceIds || !destinationIds || !separatrixIds || !separatrixTypes
172 || !separatrixFunctionMaxima || !separatrixFunctionMinima
173 || !separatrixFunctionDiffs || !isOnBoundary) {
174 this->printErr("1-separatrices vtkDataArray allocation problem.");
175 return -1;
176 }
177#endif
178
179 pointsCoords->SetNumberOfComponents(3);
180 setArray(pointsCoords, separatrices1_.pt.points_);
181
182 smoothingMask->SetNumberOfComponents(1);
183 smoothingMask->SetName(ttk::MaskScalarFieldName);
184 setArray(smoothingMask, separatrices1_.pt.smoothingMask_);
185
186 cellDimensions->SetNumberOfComponents(1);
187 cellDimensions->SetName(ttk::MorseSmaleCellDimensionName);
188 setArray(cellDimensions, separatrices1_.pt.cellDimensions_);
189
190 cellIds->SetNumberOfComponents(1);
191 cellIds->SetName(ttk::MorseSmaleCellIdName);
192 setArray(cellIds, separatrices1_.pt.cellIds_);
193
194 sourceIds->SetNumberOfComponents(1);
195 sourceIds->SetName(ttk::MorseSmaleSourceIdName);
196 setArray(sourceIds, separatrices1_.cl.sourceIds_);
197
198 destinationIds->SetNumberOfComponents(1);
199 destinationIds->SetName(ttk::MorseSmaleDestinationIdName);
200 setArray(destinationIds, separatrices1_.cl.destinationIds_);
201
202 separatrixIds->SetNumberOfComponents(1);
203 separatrixIds->SetName(ttk::MorseSmaleSeparatrixIdName);
204 setArray(separatrixIds, separatrices1_.cl.separatrixIds_);
205
206 separatrixTypes->SetNumberOfComponents(1);
207 separatrixTypes->SetName(ttk::MorseSmaleSeparatrixTypeName);
208 setArray(separatrixTypes, separatrices1_.cl.separatrixTypes_);
209
210 separatrixFunctionMaxima->SetNumberOfComponents(1);
211 separatrixFunctionMaxima->SetName(ttk::MorseSmaleSeparatrixMaximumName);
212 separatrixFunctionMaxima->SetNumberOfTuples(
213 separatrices1_.cl.numberOfCells_);
214
215 separatrixFunctionMinima->SetNumberOfComponents(1);
216 separatrixFunctionMinima->SetName(ttk::MorseSmaleSeparatrixMinimumName);
217 separatrixFunctionMinima->SetNumberOfTuples(
218 separatrices1_.cl.numberOfCells_);
219
220 separatrixFunctionDiffs->SetNumberOfComponents(1);
221 separatrixFunctionDiffs->SetName(ttk::MorseSmaleSeparatrixDifferenceName);
222 separatrixFunctionDiffs->SetNumberOfTuples(
223 separatrices1_.cl.numberOfCells_);
224
225#ifdef TTK_ENABLE_OPENMP
226#pragma omp parallel for num_threads(this->threadNumber_)
227#endif // TTK_ENABLE_OPENMP
228 for(SimplexId i = 0; i < separatrices1_.cl.numberOfCells_; ++i) {
229 const auto sepId = separatrices1_.cl.separatrixIds_[i];
230 // inputScalars->GetTuple1 not thread safe...
231 const auto min = scalars[separatrices1_.cl.sepFuncMinId_[sepId]];
232 const auto max = scalars[separatrices1_.cl.sepFuncMaxId_[sepId]];
233 separatrixFunctionMinima->SetTuple1(i, min);
234 separatrixFunctionMaxima->SetTuple1(i, max);
235 separatrixFunctionDiffs->SetTuple1(i, max - min);
236 }
237
238 isOnBoundary->SetNumberOfComponents(1);
239 isOnBoundary->SetName(ttk::MorseSmaleCriticalPointsOnBoundaryName);
240 setArray(isOnBoundary, separatrices1_.cl.isOnBoundary_);
241
242 vtkNew<ttkSimplexIdTypeArray> offsets{}, connectivity{};
243 offsets->SetNumberOfComponents(1);
244 offsets->SetNumberOfTuples(separatrices1_.cl.numberOfCells_ + 1);
245 connectivity->SetNumberOfComponents(1);
246 setArray(connectivity, separatrices1_.cl.connectivity_);
247
248#ifdef TTK_ENABLE_OPENMP
249#pragma omp parallel for num_threads(this->threadNumber_)
250#endif // TTK_ENABLE_OPENMP
251 for(SimplexId i = 0; i < separatrices1_.cl.numberOfCells_ + 1; ++i) {
252 offsets->SetTuple1(i, 2 * i);
253 }
254
255 vtkNew<vtkPoints> points{};
256 points->SetData(pointsCoords);
257 outputSeparatrices1->SetPoints(points);
258 vtkNew<vtkCellArray> cells{};
259#ifndef TTK_ENABLE_64BIT_IDS
260 cells->Use32BitStorage();
261#endif // TTK_ENABLE_64BIT_IDS
262 cells->SetData(offsets, connectivity);
263 outputSeparatrices1->SetLines(cells);
264
265 auto pointData = outputSeparatrices1->GetPointData();
266 auto cellData = outputSeparatrices1->GetCellData();
267
268#ifndef TTK_ENABLE_KAMIKAZE
269 if(!pointData || !cellData) {
270 this->printErr("outputSeparatrices1 has no point or no cell data.");
271 return -1;
272 }
273#endif
274
275 pointData->AddArray(smoothingMask);
276 pointData->AddArray(cellDimensions);
277 pointData->AddArray(cellIds);
278
279 cellData->AddArray(sourceIds);
280 cellData->AddArray(destinationIds);
281 cellData->AddArray(separatrixIds);
282 cellData->SetScalars(separatrixTypes);
283 cellData->AddArray(separatrixFunctionMaxima);
284 cellData->AddArray(separatrixFunctionMinima);
285 cellData->AddArray(separatrixFunctionDiffs);
286 cellData->AddArray(isOnBoundary);
287 }
288
289 // 2-separatrices
290 if(dimensionality == 3
292
293 vtkNew<vtkFloatArray> pointsCoords{};
294 vtkNew<ttkSimplexIdTypeArray> sourceIds{};
295 vtkNew<ttkSimplexIdTypeArray> separatrixIds{};
296 vtkNew<vtkSignedCharArray> separatrixTypes{};
297 vtkNew<vtkDoubleArray> separatrixFunctionMaxima{};
298 vtkNew<vtkDoubleArray> separatrixFunctionMinima{};
299 vtkNew<vtkDoubleArray> separatrixFunctionDiffs{};
300 vtkNew<vtkSignedCharArray> isOnBoundary{};
301
302#ifndef TTK_ENABLE_KAMIKAZE
303 if(!pointsCoords || !sourceIds || !separatrixIds || !separatrixTypes
304 || !separatrixFunctionMaxima || !separatrixFunctionMinima
305 || !separatrixFunctionDiffs || !isOnBoundary) {
306 this->printErr("2-separatrices vtkDataArray allocation problem.");
307 return -1;
308 }
309#endif
310
311 pointsCoords->SetNumberOfComponents(3);
312 setArray(pointsCoords, separatrices2_.pt.points_);
313
314 sourceIds->SetNumberOfComponents(1);
315 sourceIds->SetName(ttk::MorseSmaleSourceIdName);
316 setArray(sourceIds, separatrices2_.cl.sourceIds_);
317
318 separatrixIds->SetNumberOfComponents(1);
319 separatrixIds->SetName(ttk::MorseSmaleSeparatrixIdName);
320 setArray(separatrixIds, separatrices2_.cl.separatrixIds_);
321
322 separatrixTypes->SetNumberOfComponents(1);
323 separatrixTypes->SetName(ttk::MorseSmaleSeparatrixTypeName);
324 setArray(separatrixTypes, separatrices2_.cl.separatrixTypes_);
325
326 separatrixFunctionMaxima->SetNumberOfComponents(1);
327 separatrixFunctionMaxima->SetName(ttk::MorseSmaleSeparatrixMaximumName);
328 separatrixFunctionMaxima->SetNumberOfTuples(
329 separatrices2_.cl.numberOfCells_);
330
331 separatrixFunctionMinima->SetNumberOfComponents(1);
332 separatrixFunctionMinima->SetName(ttk::MorseSmaleSeparatrixMinimumName);
333 separatrixFunctionMinima->SetNumberOfTuples(
334 separatrices2_.cl.numberOfCells_);
335
336 separatrixFunctionDiffs->SetNumberOfComponents(1);
337 separatrixFunctionDiffs->SetName(ttk::MorseSmaleSeparatrixDifferenceName);
338 separatrixFunctionDiffs->SetNumberOfTuples(
339 separatrices2_.cl.numberOfCells_);
340
341#ifdef TTK_ENABLE_OPENMP
342#pragma omp parallel for num_threads(this->threadNumber_)
343#endif // TTK_ENABLE_OPENMP
344 for(SimplexId i = 0; i < separatrices2_.cl.numberOfCells_; ++i) {
345 const auto sepId = separatrices2_.cl.separatrixIds_[i];
346 // inputScalars->GetTuple1 not thread safe...
347 const auto min = scalars[separatrices2_.cl.sepFuncMinId_[sepId]];
348 const auto max = scalars[separatrices2_.cl.sepFuncMaxId_[sepId]];
349 separatrixFunctionMinima->SetTuple1(i, min);
350 separatrixFunctionMaxima->SetTuple1(i, max);
351 separatrixFunctionDiffs->SetTuple1(i, max - min);
352 }
353
354 isOnBoundary->SetNumberOfComponents(1);
355 isOnBoundary->SetName(ttk::MorseSmaleCriticalPointsOnBoundaryName);
356 setArray(isOnBoundary, separatrices2_.cl.isOnBoundary_);
357
358 vtkNew<ttkSimplexIdTypeArray> offsets{}, connectivity{};
359 offsets->SetNumberOfComponents(1);
360 setArray(offsets, separatrices2_.cl.offsets_);
361 connectivity->SetNumberOfComponents(1);
362 setArray(connectivity, separatrices2_.cl.connectivity_);
363
364 vtkNew<vtkPoints> points{};
365 points->SetData(pointsCoords);
366 outputSeparatrices2->SetPoints(points);
367 vtkNew<vtkCellArray> cells{};
368#ifndef TTK_ENABLE_64BIT_IDS
369 cells->Use32BitStorage();
370#endif // TTK_ENABLE_64BIT_IDS
371 cells->SetData(offsets, connectivity);
372 outputSeparatrices2->SetPolys(cells);
373
374 auto cellData = outputSeparatrices2->GetCellData();
375#ifndef TTK_ENABLE_KAMIKAZE
376 if(!cellData) {
377 this->printErr("outputSeparatrices2 has no cell data.");
378 return -1;
379 }
380#endif
381
382 cellData->AddArray(sourceIds);
383 cellData->AddArray(separatrixIds);
384 cellData->AddArray(separatrixTypes);
385 cellData->AddArray(separatrixFunctionMaxima);
386 cellData->AddArray(separatrixFunctionMinima);
387 cellData->AddArray(separatrixFunctionDiffs);
388 cellData->AddArray(isOnBoundary);
389 }
390
391 return ret;
392}
393
395 vtkInformationVector **inputVector,
396 vtkInformationVector *outputVector) {
397
398 const auto input = vtkDataSet::GetData(inputVector[0]);
399 auto outputCriticalPoints = vtkPolyData::GetData(outputVector, 0);
400 auto outputSeparatrices1 = vtkPolyData::GetData(outputVector, 1);
401 auto outputSeparatrices2 = vtkPolyData::GetData(outputVector, 2);
402 auto outputMorseComplexes = vtkDataSet::GetData(outputVector, 3);
403
404#ifndef TTK_ENABLE_KAMIKAZE
405 if(!input) {
406 this->printErr("Input pointer is NULL.");
407 return -1;
408 }
409 if(input->GetNumberOfPoints() == 0) {
410 this->printErr("Input has no point.");
411 return -1;
412 }
413 if(!outputCriticalPoints or !outputSeparatrices1 or !outputSeparatrices2
414 or !outputMorseComplexes) {
415 this->printErr("Output pointers are NULL.");
416 return -1;
417 }
418#endif
419
420 const auto triangulation = ttkAlgorithm::GetTriangulation(input);
421 if(triangulation == nullptr) {
422 this->printErr("Triangulation is null");
423 return 0;
424 }
425 this->preconditionTriangulation(triangulation);
426
427 const auto inputScalars = this->GetInputArrayToProcess(0, inputVector);
428
429#ifndef TTK_ENABLE_KAMIKAZE
430 if(inputScalars == nullptr) {
431 this->printErr("wrong scalars.");
432 return -1;
433 }
434#endif
435
436 auto inputOffsets = ttkAlgorithm::GetOrderArray(
437 input, 0, triangulation, false, 1, this->ForceInputOffsetScalarField);
438
439#ifndef TTK_ENABLE_KAMIKAZE
440 if(inputOffsets == nullptr) {
441 this->printErr("wrong offsets.");
442 return -1;
443 }
444 if(inputOffsets->GetDataType() != VTK_INT
445 and inputOffsets->GetDataType() != VTK_ID_TYPE) {
446 this->printErr("input offset field type not supported.");
447 return -1;
448 }
449#endif
450
451 this->printMsg("Launching computation on field `"
452 + std::string(inputScalars->GetName()) + "'...");
453
454 // morse complexes
455 const SimplexId numberOfVertices = triangulation->getNumberOfVertices();
456#ifndef TTK_ENABLE_KAMIKAZE
457 if(!numberOfVertices) {
458 this->printErr("Input has no vertices.");
459 return -1;
460 }
461#endif
462
463 vtkNew<ttkSimplexIdTypeArray> ascendingManifold{};
464 vtkNew<ttkSimplexIdTypeArray> descendingManifold{};
465 vtkNew<ttkSimplexIdTypeArray> morseSmaleManifold{};
466#ifndef TTK_ENABLE_KAMIKAZE
467 if(!ascendingManifold || !descendingManifold || !morseSmaleManifold) {
468 this->printErr("Manifold vtkDataArray allocation problem.");
469 return -1;
470 }
471#endif
472 ascendingManifold->SetNumberOfComponents(1);
473 ascendingManifold->SetNumberOfTuples(numberOfVertices);
474 ascendingManifold->SetName(ttk::MorseSmaleAscendingName);
475
476 descendingManifold->SetNumberOfComponents(1);
477 descendingManifold->SetNumberOfTuples(numberOfVertices);
478 descendingManifold->SetName(ttk::MorseSmaleDescendingName);
479
480 morseSmaleManifold->SetNumberOfComponents(1);
481 morseSmaleManifold->SetNumberOfTuples(numberOfVertices);
482 morseSmaleManifold->SetName(ttk::MorseSmaleManifoldName);
483
484 this->segmentations_ = {ttkUtils::GetPointer<SimplexId>(ascendingManifold),
485 ttkUtils::GetPointer<SimplexId>(descendingManifold),
486 ttkUtils::GetPointer<SimplexId>(morseSmaleManifold)};
487
491
492 int ret{};
493
495 inputScalars->GetDataType(), triangulation->getType(),
496 (ret = dispatch<VTK_TT, TTK_TT>(
497 inputScalars, outputCriticalPoints, outputSeparatrices1,
498 outputSeparatrices2, ttkUtils::GetPointer<SimplexId>(inputOffsets),
499 *static_cast<TTK_TT *>(triangulation->getData()))));
500
501 if(ret != 0) {
502 return -1;
503 }
504
505 outputMorseComplexes->ShallowCopy(input);
506 // morse complexes
508 vtkPointData *pointData = outputMorseComplexes->GetPointData();
509#ifndef TTK_ENABLE_KAMIKAZE
510 if(!pointData) {
511 this->printErr("outputMorseComplexes has no point data.");
512 return -1;
513 }
514#endif
515
516 if(ComputeDescendingSegmentation)
517 pointData->AddArray(descendingManifold);
519 pointData->AddArray(ascendingManifold);
522 pointData->AddArray(morseSmaleManifold);
523 }
524
525 return !ret;
526}
#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()
ttk::Triangulation * GetTriangulation(vtkDataSet *dataSet)
vtkDataArray * GetOrderArray(vtkDataSet *const inputData, const int scalarArrayIdx, ttk::Triangulation *triangulation, const bool getGlobalOrder=false, const int orderArrayIdx=0, const bool enforceOrderArrayIdx=false)
TTK VTK-filter that wraps the morseSmaleComplex processing package.
int dispatch(vtkDataArray *const inputScalars, vtkPolyData *const outputCriticalPoints, vtkPolyData *const outputSeparatrices1, vtkPolyData *const outputSeparatrices2, const SimplexId *const inputOffsets, const triangulationType &triangulation)
int FillInputPortInformation(int port, vtkInformation *info) override
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
int FillOutputPortInformation(int port, vtkInformation *info) override
static int CellVertexFromPoints(vtkDataSet *const dataSet, vtkPoints *const points)
Definition ttkUtils.cpp:328
static void SetVoidArray(vtkDataArray *array, void *data, vtkIdType size, int save)
Definition ttkUtils.cpp:280
void setDebugMsgPrefix(const std::string &prefix)
Definition Debug.h:364
int printErr(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition Debug.h:149
void setSaddleConnectorsPersistenceThreshold(const double threshold)
void preconditionTriangulation(AbstractTriangulation *const data)
void setReturnSaddleConnectors(const bool state)
int execute(OutputCriticalPoints &outCP, Output1Separatrices &outSeps1, Output2Separatrices &outSeps2, OutputManifold &outManifold, const dataType *const scalars, const size_t scalarsMTime, const SimplexId *const offsets, const triangulationType &triangulation)
const char MorseSmaleDestinationIdName[]
Definition DataTypes.h:53
const char MorseSmaleSeparatrixDifferenceName[]
Definition DataTypes.h:59
const char MorseSmaleManifoldName[]
Definition DataTypes.h:64
const char MorseSmaleSeparatrixMaximumName[]
Definition DataTypes.h:56
const char MorseSmaleSourceIdName[]
Definition DataTypes.h:52
const char MorseSmaleAscendingName[]
Definition DataTypes.h:62
const char MorseSmaleCellDimensionName[]
Definition DataTypes.h:48
const char MorseSmaleDescendingName[]
Definition DataTypes.h:63
const char MorseSmaleCriticalPointsOnBoundaryName[]
Definition DataTypes.h:61
const char MorseSmaleCellIdName[]
Definition DataTypes.h:49
const char MorseSmaleSeparatrixTypeName[]
Definition DataTypes.h:55
const char MorseSmaleManifoldSizeName[]
Definition DataTypes.h:51
const char MorseSmaleBoundaryName[]
Definition DataTypes.h:50
const char VertexScalarFieldName[]
default name for vertex scalar field
Definition DataTypes.h:35
const char MorseSmaleSeparatrixMinimumName[]
Definition DataTypes.h:57
const char MorseSmaleSeparatrixIdName[]
Definition DataTypes.h:54
const char MaskScalarFieldName[]
default name for mask scalar field
Definition DataTypes.h:32
std::vector< ttk::SimplexId > destinationIds_
std::vector< ttk::SimplexId > separatrixIds_
struct ttk::MorseSmaleComplex::Output1Separatrices::@1 cl
std::vector< ttk::SimplexId > connectivity_
struct ttk::MorseSmaleComplex::Output1Separatrices::@0 pt
struct ttk::MorseSmaleComplex::Output2Separatrices::@2 pt
struct ttk::MorseSmaleComplex::Output2Separatrices::@3 cl
std::vector< ttk::SimplexId > connectivity_
std::vector< ttk::SimplexId > separatrixIds_
std::vector< std::array< float, 3 > > points_
#define ttkVtkTemplateMacro(dataType, triangulationType, call)
Definition ttkMacros.h:69
vtkStandardNewMacro(ttkMorseSmaleComplex)
void setArray(vtkArrayType &vtkArray, vectorType &vector)
printMsg(debug::output::BOLD+" | | | | | . \\ | | (__| | / __/| |_| / __/|__ _|"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)