TTK
Loading...
Searching...
No Matches
ttkTopologicalSimplificationByPersistence.cpp
Go to the documentation of this file.
2
3#include <vtkInformation.h>
4
5#include <vtkDataArray.h>
6#include <vtkDataSet.h>
7#include <vtkObjectFactory.h>
8#include <vtkPointData.h>
9#include <vtkSmartPointer.h>
10
11#include <ttkMacros.h>
12#include <ttkUtils.h>
13
15
18 this->setDebugMsgPrefix("PLTS");
19
20 this->SetNumberOfInputPorts(1);
21 this->SetNumberOfOutputPorts(1);
22}
23
26 = default;
27
29 int port, vtkInformation *info) {
30 if(port == 0) {
31 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkDataSet");
32 return 1;
33 }
34 return 0;
35}
36
38 int port, vtkInformation *info) {
39 if(port == 0) {
41 return 1;
42 }
43 return 0;
44}
45
47 vtkInformation *ttkNotUsed(request),
48 vtkInformationVector **inputVector,
49 vtkInformationVector *outputVector) {
50
51 // retrieve input data object
52 auto inputDataSet = vtkDataSet::GetData(inputVector[0]);
53 if(!inputDataSet)
54 return 0;
55
56 // retrieve input arrays
57 if(this->GetInputArrayAssociation(0, inputVector) != 0)
58 return !this->printErr("Input array needs to be a point data array.");
59
60 auto inputScalars = this->GetInputArrayToProcess(0, inputVector);
61 if(!inputScalars)
62 return !this->printErr("Unable to retrieve input array.");
63
64 if(inputScalars->GetNumberOfComponents() != 1)
65 return !this->printErr("Input array needs to be a scalar array.");
66
67 auto inputOrder = this->GetOrderArray(inputDataSet, 0);
68 if(!inputOrder)
69 return 0;
70
71 // create output arrays
72 auto outputScalars
73 = vtkSmartPointer<vtkDataArray>::Take(inputScalars->NewInstance());
74 outputScalars->DeepCopy(inputScalars);
75 auto outputOrder
76 = vtkSmartPointer<vtkDataArray>::Take(inputOrder->NewInstance());
77 outputOrder->DeepCopy(inputOrder);
78
79 // retrieve and precondition triangulation
80 auto triangulation = ttkAlgorithm::GetTriangulation(inputDataSet);
81 if(!triangulation)
82 return 0;
83 this->preconditionTriangulation(triangulation);
84
85 double persistenceThreshold = this->PersistenceThreshold;
86 if(!this->ThresholdIsAbsolute) {
87 double range[2];
88 outputScalars->GetRange(range);
89 persistenceThreshold = (range[1] - range[0]) * persistenceThreshold;
90 }
91
92 // perform simplification
93 int status = 0;
95 outputScalars->GetDataType(), triangulation->getType(),
96 (status = this->removeNonPersistentExtrema<VTK_TT, ttk::SimplexId, TTK_TT>(
97 ttkUtils::GetPointer<VTK_TT>(outputScalars),
98 ttkUtils::GetPointer<ttk::SimplexId>(outputOrder),
99
100 static_cast<TTK_TT *>(triangulation->getData()), persistenceThreshold,
101 this->ComputePerturbation, this->PairType)));
102
103 // On error cancel filter execution
104 if(status != 0)
105 return 0;
106
107 auto outputDataSet = vtkDataSet::GetData(outputVector, 0);
108 outputDataSet->ShallowCopy(inputDataSet);
109
110 auto outputDataSetPD = outputDataSet->GetPointData();
111 outputDataSetPD->AddArray(outputScalars);
112 outputDataSetPD->AddArray(outputOrder);
113
114 return 1;
115}
#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 computes a persistence-based simplification of a scalar field.
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
int FillInputPortInformation(int port, vtkInformation *info) override
int FillOutputPortInformation(int port, vtkInformation *info) override
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
int preconditionTriangulation(ttk::AbstractTriangulation *triangulation) const
int SimplexId
Identifier type for simplices of any dimension.
Definition: DataTypes.h:22
#define ttkVtkTemplateMacro(dataType, triangulationType, call)
Definition: ttkMacros.h:69
vtkStandardNewMacro(ttkTopologicalSimplificationByPersistence)