TTK
Loading...
Searching...
No Matches
ttkMandatoryCriticalPoints.cpp
Go to the documentation of this file.
1#include <ttkMacros.h>
3#include <ttkUtils.h>
4
5#include <vtkCellData.h>
6#include <vtkDataObject.h>
7#include <vtkDoubleArray.h>
8#include <vtkInformation.h>
9#include <vtkIntArray.h>
10#include <vtkObjectFactory.h>
11#include <vtkPointData.h>
12#include <vtkUnstructuredGrid.h>
13
14#include <array>
15
17
19 SetNumberOfInputPorts(1);
20 SetNumberOfOutputPorts(6);
21}
22
24 vtkInformation *info) {
25 if(port == 0) {
26 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkDataSet");
27 return 1;
28 }
29 return 0;
30}
31
33 int port, vtkInformation *info) {
34 if(port >= 0 && port < 4) {
36 return 1;
37 } else if(port == 4 || port == 5) {
38 info->Set(vtkDataObject::DATA_TYPE_NAME(), "vtkUnstructuredGrid");
39 return 1;
40 }
41 return 0;
42}
43
44static void buildVtkTree(
45 vtkUnstructuredGrid *outputTree,
46 const ttk::Graph &graph,
47 const std::vector<double> &xCoord,
48 const std::vector<double> &yCoord,
49 const std::vector<int> &mdtTreePointComponentId,
50 const std::vector<ttk::MandatoryCriticalPoints::PointType> &mdtTreePointType,
51 const std::vector<double> &mdtTreePointLowInterval,
52 const std::vector<double> &mdtTreePointUpInterval,
53 const std::vector<int> &mdtTreeEdgeSwitchable) {
54
55 const int numberOfPoints = graph.getNumberOfVertices();
56 const int numberOfEdges = graph.getNumberOfEdges();
57
58 /* Point data */
59 vtkNew<vtkIntArray> outputTreePointType{};
60 outputTreePointType->SetName("Type");
61 outputTreePointType->SetNumberOfTuples(numberOfPoints);
62 outputTree->GetPointData()->AddArray(outputTreePointType);
63
64 vtkNew<vtkDoubleArray> outputTreePointLowInterval{};
65 outputTreePointLowInterval->SetName("LowInterval");
66 outputTreePointLowInterval->SetNumberOfTuples(numberOfPoints);
67 outputTree->GetPointData()->AddArray(outputTreePointLowInterval);
68
69 vtkNew<vtkDoubleArray> outputTreePointUpInterval{};
70 outputTreePointUpInterval->SetName("UpInterval");
71 outputTreePointUpInterval->SetNumberOfTuples(numberOfPoints);
72 outputTree->GetPointData()->AddArray(outputTreePointUpInterval);
73
74 vtkNew<vtkIntArray> outputTreePointComponentId{};
75 outputTreePointComponentId->SetName("ComponentId");
76 outputTreePointComponentId->SetNumberOfTuples(numberOfPoints);
77 outputTree->GetPointData()->AddArray(outputTreePointComponentId);
78
79 for(int i = 0; i < numberOfPoints; i++) {
80 outputTreePointType->SetTuple1(i, static_cast<int>(mdtTreePointType[i]));
81 outputTreePointLowInterval->SetTuple1(i, mdtTreePointLowInterval[i]);
82 outputTreePointUpInterval->SetTuple1(i, mdtTreePointUpInterval[i]);
83 outputTreePointComponentId->SetTuple1(i, mdtTreePointComponentId[i]);
84 }
85
86 /* Cell data */
87 vtkNew<vtkIntArray> outputTreeEdgeSwitchable{};
88 outputTreeEdgeSwitchable->SetName("Switchable");
89 outputTreeEdgeSwitchable->SetNumberOfTuples(numberOfEdges);
90 outputTree->GetCellData()->AddArray(outputTreeEdgeSwitchable);
91 for(int i = 0; i < numberOfEdges; i++) {
92 outputTreeEdgeSwitchable->SetTuple1(i, mdtTreeEdgeSwitchable[i]);
93 }
94
95 /* New Objects */
96 // Points
97 vtkNew<vtkPoints> mdtTreePoints{};
98 outputTree->SetPoints(mdtTreePoints);
99 for(int i = 0; i < numberOfPoints; i++) {
100 mdtTreePoints->InsertNextPoint(xCoord[i], yCoord[i], 0.0);
101 }
102
103 // Edges (cells)
104 outputTree->Allocate(numberOfEdges);
105 for(int i = 0; i < numberOfEdges; i++) {
106 std::array<vtkIdType, 2> edge{graph.getEdge(i).getVertexIdx().first,
107 graph.getEdge(i).getVertexIdx().second};
108 outputTree->InsertNextCell(VTK_LINE, 2, edge.data());
109 }
110}
111
113 vtkInformation *ttkNotUsed(request),
114 vtkInformationVector **inputVector,
115 vtkInformationVector *outputVector) {
116
117 const auto input = vtkDataSet::GetData(inputVector[0]);
118 auto outputMinimum = vtkDataSet::GetData(outputVector, 0);
119 auto outputJoinSaddle = vtkDataSet::GetData(outputVector, 1);
120 auto outputSplitSaddle = vtkDataSet::GetData(outputVector, 2);
121 auto outputMaximum = vtkDataSet::GetData(outputVector, 3);
122 auto outputJoinTree = vtkUnstructuredGrid::GetData(outputVector, 4);
123 auto outputSplitTree = vtkUnstructuredGrid::GetData(outputVector, 5);
124
125 if(input == nullptr) {
126 return 0;
127 }
128
129 // Use a pointer-base copy for the input data
130 outputMinimum->ShallowCopy(input);
131 outputJoinSaddle->ShallowCopy(input);
132 outputSplitSaddle->ShallowCopy(input);
133 outputMaximum->ShallowCopy(input);
134
135 // Input data arrays
136 const auto inputLowerBoundField = this->GetInputArrayToProcess(0, input);
137 const auto inputUpperBoundField = this->GetInputArrayToProcess(1, input);
138
139 if(inputLowerBoundField == nullptr || inputUpperBoundField == nullptr) {
140 return 0;
141 }
142
143 this->printMsg("Using `" + std::string{inputLowerBoundField->GetName()}
144 + "' as lower bound...");
145 this->printMsg("Using `" + std::string{inputUpperBoundField->GetName()}
146 + "' as upper bound...");
147
148 // Initialize the triangulation object with the input data set
149 auto triangulation = ttkAlgorithm::GetTriangulation(input);
150
151 if(triangulation == nullptr)
152 return -1;
153
154 this->preconditionTriangulation(triangulation);
155
156 bool hasChangedConnectivity = false;
157
158 if(triangulation->isEmpty()) {
159 Modified();
160 hasChangedConnectivity = true;
161 }
162
163 // Allocate the memory for the output scalar field
164 vtkNew<vtkIntArray> outputMandatoryMinimum{};
165 outputMandatoryMinimum->SetNumberOfTuples(input->GetNumberOfPoints());
166 outputMandatoryMinimum->SetName("MinimumComponents");
167
168 vtkNew<vtkIntArray> outputMandatoryJoinSaddle{};
169 outputMandatoryJoinSaddle->SetNumberOfTuples(input->GetNumberOfPoints());
170 outputMandatoryJoinSaddle->SetName("JoinSaddleComponents");
171
172 vtkNew<vtkIntArray> outputMandatorySplitSaddle{};
173 outputMandatorySplitSaddle->SetNumberOfTuples(input->GetNumberOfPoints());
174 outputMandatorySplitSaddle->SetName("SplitSaddleComponents");
175
176 vtkNew<vtkIntArray> outputMandatoryMaximum{};
177 outputMandatoryMaximum->SetNumberOfTuples(input->GetNumberOfPoints());
178 outputMandatoryMaximum->SetName("MaximumComponents");
179
180 // Add array in the output data set
181 outputMinimum->GetPointData()->AddArray(outputMandatoryMinimum);
182 outputJoinSaddle->GetPointData()->AddArray(outputMandatoryJoinSaddle);
183 outputSplitSaddle->GetPointData()->AddArray(outputMandatorySplitSaddle);
184 outputMaximum->GetPointData()->AddArray(outputMandatoryMaximum);
185
186 // Reset the baseCode object
187 if(hasChangedConnectivity)
188 this->flush();
189
190 // Set the number of vertex
191 this->setVertexNumber(input->GetNumberOfPoints());
192 // Set the coordinates of each vertex
193 std::array<double, 3> point{};
194 for(int i = 0; i < input->GetNumberOfPoints(); i++) {
195 input->GetPoint(i, point.data());
196 this->setVertexPosition(i, point.data());
197 }
198 // Set the void pointers to the upper and lower bound fields
200 ttkUtils::GetVoidPointer(inputLowerBoundField));
202 ttkUtils::GetVoidPointer(inputUpperBoundField));
203 // Set the output data pointers
205 ttkUtils::GetVoidPointer(outputMandatoryMinimum));
207 ttkUtils::GetVoidPointer(outputMandatoryJoinSaddle));
209 ttkUtils::GetVoidPointer(outputMandatorySplitSaddle));
211 ttkUtils::GetVoidPointer(outputMandatoryMaximum));
212 // Set the triangulation object
213 // Set the offsets
214 this->setSoSoffsets();
215 // Simplification threshold
216
217 this->setSimplificationThreshold(simplificationThreshold_);
218
219 // Execute process
220 ttkVtkTemplateMacro(inputUpperBoundField->GetDataType(),
221 triangulation->getType(),
222 (this->execute<VTK_TT, TTK_TT>(
223 *static_cast<TTK_TT *>(triangulation->getData()))));
224
225 computeMinimumOutput_ = true;
226 computeJoinSaddleOutput_ = true;
227 computeSplitSaddleOutput_ = true;
228 computeMaximumOutput_ = true;
229
230 // Simplification
231 if(simplify_) {
232 // Simplify
233 this->simplifyJoinTree();
235 this->simplifySplitTree();
237 // Recompute the outputs
238 computeMinimumOutput_ = true;
239 computeJoinSaddleOutput_ = true;
240 computeSplitSaddleOutput_ = true;
241 computeMaximumOutput_ = true;
242 }
243 // Outputs
244 if(computeMinimumOutput_) {
245 if(outputAllMinimumComponents_)
246 this->computeAllMinima();
247 else
248 this->computeMinimum(outputMinimumComponentId_);
249 computeMinimumOutput_ = false;
250 }
251 if(computeJoinSaddleOutput_) {
252 if(outputAllJoinSaddleComponents_)
253 this->computeAllJoinSaddle(*triangulation);
254 else
255 this->computeJoinSaddle(outputJoinSaddleComponentId_, *triangulation);
256 computeJoinSaddleOutput_ = false;
257 }
258 if(computeSplitSaddleOutput_) {
259 if(outputAllSplitSaddleComponents_)
260 this->computeAllSplitSaddle(*triangulation);
261 else
262 this->computeSplitSaddle(outputSplitSaddleComponentId_, *triangulation);
263 computeSplitSaddleOutput_ = false;
264 }
265 if(computeMaximumOutput_) {
266 if(outputAllMaximumComponents_)
267 this->computeAllMaxima();
268 else
269 this->computeMaximum(outputMaximumComponentId_);
270 computeMaximumOutput_ = false;
271 }
272
273 buildVtkTree(outputJoinTree, mdtJoinTree_, mdtJoinTreePointXCoord_,
277 buildVtkTree(outputSplitTree, mdtSplitTree_, mdtSplitTreePointXCoord_,
281
282 return 1;
283}
#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)
TTK VTK-filter for the computation of mandatory critical points in uncertain scalar data.
int FillOutputPortInformation(int port, vtkInformation *info) override
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
int FillInputPortInformation(int port, vtkInformation *info) override
static void * GetVoidPointer(vtkDataArray *array, vtkIdType start=0)
Definition ttkUtils.cpp:226
const std::pair< int, int > & getVertexIdx() const
const Edge & getEdge(const int idx) const
int getNumberOfVertices() const
int getNumberOfEdges() const
int computeJoinSaddle(const int &id, const triangulationType &triangulation, const bool &reset=true)
std::vector< double > mdtJoinTreePointXCoord_
std::vector< double > mdtSplitTreePointLowInterval_
std::vector< double > mdtJoinTreePointLowInterval_
std::vector< PointType > mdtJoinTreePointType_
int computeAllSplitSaddle(const triangulationType &triangulation)
int computeSplitSaddle(const int &id, const triangulationType &triangulation, const bool &reset=true)
int computeMaximum(const int &id, const bool &reset=true, const bool &ttkNotUsed(parallel)=true)
std::vector< PointType > mdtSplitTreePointType_
int setVertexPosition(const int &i, const double point[3])
int setVertexNumber(const int &vertexNumber)
std::vector< double > mdtJoinTreePointUpInterval_
std::vector< int > mdtJoinTreePointComponentId_
int setSoSoffsets(int *offsets=nullptr)
std::vector< int > mdtSplitTreePointComponentId_
int computeMinimum(const int &id, const bool &reset=true, const bool &ttkNotUsed(parallel)=true)
std::vector< double > mdtSplitTreePointXCoord_
int setSimplificationThreshold(double normalizedThreshold)
std::vector< double > mdtSplitTreePointYCoord_
std::vector< double > mdtJoinTreePointYCoord_
std::vector< double > mdtSplitTreePointUpInterval_
void preconditionTriangulation(AbstractTriangulation *const triangulation)
int computeAllJoinSaddle(const triangulationType &triangulation)
#define ttkVtkTemplateMacro(dataType, triangulationType, call)
Definition ttkMacros.h:69
vtkStandardNewMacro(ttkMandatoryCriticalPoints)
printMsg(debug::output::BOLD+" | | | | | . \\ | | (__| | / __/| |_| / __/|__ _|"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)