TTK
Loading...
Searching...
No Matches
ttkTopologicalSimplification.cpp
Go to the documentation of this file.
1#include <vtkDataArray.h>
2#include <vtkDataSet.h>
3#include <vtkIdTypeArray.h>
4#include <vtkInformation.h>
5#include <vtkIntArray.h>
6#include <vtkNew.h>
7#include <vtkPointData.h>
8#include <vtkPointSet.h>
9#include <vtkSmartPointer.h>
10
11#include <ttkMacros.h>
13#include <ttkUtils.h>
14
16
18 this->SetNumberOfInputPorts(2);
19 this->SetNumberOfOutputPorts(1);
20}
21
23 int port, vtkInformation *info) {
24 if(port == 0) {
25 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkDataSet");
26 return 1;
27 } else if(port == 1) {
28 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkPointSet");
29 return 1;
30 }
31 return 0;
32}
33
35 int port, vtkInformation *info) {
36 if(port == 0) {
38 return 1;
39 }
40 return 0;
41}
42
44 vtkInformation *ttkNotUsed(request),
45 vtkInformationVector **inputVector,
46 vtkInformationVector *outputVector) {
47
48 using ttk::SimplexId;
49
50 // Warning: this needs to be done before the preconditioning.
51 if(!this->UseLTS) {
53 }
54
55 const auto domain = vtkDataSet::GetData(inputVector[0]);
56 const auto constraints = vtkPointSet::GetData(inputVector[1]);
57 if(!domain || !constraints)
58 return !this->printErr("Unable to retrieve required input data objects.");
59
60 auto output = vtkDataSet::GetData(outputVector);
61
62 // triangulation
63 auto triangulation = ttkAlgorithm::GetTriangulation(domain);
64
65 if(!triangulation) {
66 this->printErr("Input triangulation pointer is NULL.");
67 return -1;
68 }
69
70 this->preconditionTriangulation(triangulation);
71
72 if(triangulation->isEmpty()) {
73 this->printErr("Triangulation allocation problem.");
74 return -1;
75 }
76
77 if(!domain) {
78 this->printErr("Input pointer is NULL.");
79 return -1;
80 }
81
82 const auto numberOfVertices = domain->GetNumberOfPoints();
83 if(numberOfVertices <= 0) {
84 this->printErr("Domain has no points.");
85 return -5;
86 }
87
88 // domain scalar field
89 const auto inputScalars = this->GetInputArrayToProcess(0, domain);
90 if(!inputScalars) {
91 this->printErr("Input scalar field pointer is null.");
92 return -3;
93 }
94
95 // domain offset field
96 const auto inputOrder
97 = this->GetOrderArray(domain, 0, 2, ForceInputOffsetScalarField);
98 if(!inputOrder) {
99 this->printErr("Wrong input offset scalar field.");
100 return -1;
101 }
102
103 // create output arrays
104 auto outputScalars
105 = vtkSmartPointer<vtkDataArray>::Take(inputScalars->NewInstance());
106 outputScalars->DeepCopy(inputScalars);
107 auto outputOrder
108 = vtkSmartPointer<vtkDataArray>::Take(inputOrder->NewInstance());
109 outputOrder->DeepCopy(inputOrder);
110
111 // constraint identifier field
112 int numberOfConstraints = constraints->GetNumberOfPoints();
113
114 std::vector<ttk::SimplexId> idSpareStorage{};
115 auto identifiers = this->GetIdentifierArrayPtr(ForceInputVertexScalarField, 1,
117 constraints, idSpareStorage);
118
119 int ret{};
120 switch(inputScalars->GetDataType()) {
121 vtkTemplateMacro(ret = this->execute(
122 ttkUtils::GetPointer<VTK_TT>(inputScalars),
123 ttkUtils::GetPointer<VTK_TT>(outputScalars), identifiers,
124 ttkUtils::GetPointer<SimplexId>(inputOrder),
125 ttkUtils::GetPointer<SimplexId>(outputOrder),
126 numberOfConstraints, this->AddPerturbation,
127 *triangulation->getData()));
128 }
129
130 // something wrong in baseCode
131 if(ret) {
132 this->printErr("TopologicalSimplification.execute() error code: "
133 + std::to_string(ret));
134 return -12;
135 }
136
137 output->ShallowCopy(domain);
138 output->GetPointData()->AddArray(outputOrder);
139 output->GetPointData()->AddArray(outputScalars);
140
141 return 1;
142}
#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::SimplexId * GetIdentifierArrayPtr(const bool &enforceArrayIndex, const int &arrayIndex, const std::string &arrayName, vtkDataSet *const inputData, std::vector< ttk::SimplexId > &spareStorage, const int inputPort=0, const bool printErr=true)
vtkDataArray * GetOrderArray(vtkDataSet *const inputData, const int scalarArrayIdx, const int orderArrayIdx=0, const bool enforceOrderArrayIdx=false)
ttk::Triangulation * GetTriangulation(vtkDataSet *dataSet)
TTK VTK-filter for the topological simplification of scalar data.
int FillOutputPortInformation(int port, vtkInformation *info) override
int FillInputPortInformation(int port, vtkInformation *info) override
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
int printErr(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition: Debug.h:149
int execute(const dataType *const inputScalars, dataType *const outputScalars, const SimplexId *const identifiers, const SimplexId *const inputOffsets, SimplexId *const offsets, const SimplexId constraintNumber, const bool addPerturbation, const triangulationType &triangulation)
int preconditionTriangulation(AbstractTriangulation *triangulation)
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(ttkTopologicalSimplification)