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>
19 SetNumberOfInputPorts(1);
20 SetNumberOfOutputPorts(6);
24 vtkInformation *info) {
26 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(),
"vtkDataSet");
33 int port, vtkInformation *info) {
34 if(port >= 0 && port < 4) {
37 }
else if(port == 4 || port == 5) {
38 info->Set(vtkDataObject::DATA_TYPE_NAME(),
"vtkUnstructuredGrid");
44static void buildVtkTree(
45 vtkUnstructuredGrid *outputTree,
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) {
59 vtkNew<vtkIntArray> outputTreePointType{};
60 outputTreePointType->SetName(
"Type");
61 outputTreePointType->SetNumberOfTuples(numberOfPoints);
62 outputTree->GetPointData()->AddArray(outputTreePointType);
64 vtkNew<vtkDoubleArray> outputTreePointLowInterval{};
65 outputTreePointLowInterval->SetName(
"LowInterval");
66 outputTreePointLowInterval->SetNumberOfTuples(numberOfPoints);
67 outputTree->GetPointData()->AddArray(outputTreePointLowInterval);
69 vtkNew<vtkDoubleArray> outputTreePointUpInterval{};
70 outputTreePointUpInterval->SetName(
"UpInterval");
71 outputTreePointUpInterval->SetNumberOfTuples(numberOfPoints);
72 outputTree->GetPointData()->AddArray(outputTreePointUpInterval);
74 vtkNew<vtkIntArray> outputTreePointComponentId{};
75 outputTreePointComponentId->SetName(
"ComponentId");
76 outputTreePointComponentId->SetNumberOfTuples(numberOfPoints);
77 outputTree->GetPointData()->AddArray(outputTreePointComponentId);
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]);
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]);
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);
104 outputTree->Allocate(numberOfEdges);
105 for(
int i = 0; i < numberOfEdges; i++) {
108 outputTree->InsertNextCell(VTK_LINE, 2, edge.data());
114 vtkInformationVector **inputVector,
115 vtkInformationVector *outputVector) {
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);
125 if(input ==
nullptr) {
130 outputMinimum->ShallowCopy(input);
131 outputJoinSaddle->ShallowCopy(input);
132 outputSplitSaddle->ShallowCopy(input);
133 outputMaximum->ShallowCopy(input);
136 const auto inputLowerBoundField = this->GetInputArrayToProcess(0, input);
137 const auto inputUpperBoundField = this->GetInputArrayToProcess(1, input);
139 if(inputLowerBoundField ==
nullptr || inputUpperBoundField ==
nullptr) {
143 this->
printMsg(
"Using `" + std::string{inputLowerBoundField->GetName()}
144 +
"' as lower bound...");
145 this->
printMsg(
"Using `" + std::string{inputUpperBoundField->GetName()}
146 +
"' as upper bound...");
151 if(triangulation ==
nullptr)
156 bool hasChangedConnectivity =
false;
158 if(triangulation->isEmpty()) {
160 hasChangedConnectivity =
true;
164 vtkNew<vtkIntArray> outputMandatoryMinimum{};
165 outputMandatoryMinimum->SetNumberOfTuples(input->GetNumberOfPoints());
166 outputMandatoryMinimum->SetName(
"MinimumComponents");
168 vtkNew<vtkIntArray> outputMandatoryJoinSaddle{};
169 outputMandatoryJoinSaddle->SetNumberOfTuples(input->GetNumberOfPoints());
170 outputMandatoryJoinSaddle->SetName(
"JoinSaddleComponents");
172 vtkNew<vtkIntArray> outputMandatorySplitSaddle{};
173 outputMandatorySplitSaddle->SetNumberOfTuples(input->GetNumberOfPoints());
174 outputMandatorySplitSaddle->SetName(
"SplitSaddleComponents");
176 vtkNew<vtkIntArray> outputMandatoryMaximum{};
177 outputMandatoryMaximum->SetNumberOfTuples(input->GetNumberOfPoints());
178 outputMandatoryMaximum->SetName(
"MaximumComponents");
181 outputMinimum->GetPointData()->AddArray(outputMandatoryMinimum);
182 outputJoinSaddle->GetPointData()->AddArray(outputMandatoryJoinSaddle);
183 outputSplitSaddle->GetPointData()->AddArray(outputMandatorySplitSaddle);
184 outputMaximum->GetPointData()->AddArray(outputMandatoryMaximum);
187 if(hasChangedConnectivity)
193 std::array<double, 3> point{};
194 for(
int i = 0; i < input->GetNumberOfPoints(); i++) {
195 input->GetPoint(i, point.data());
221 triangulation->getType(),
222 (this->execute<VTK_TT, TTK_TT>(
223 *
static_cast<TTK_TT *
>(triangulation->getData()))));
225 computeMinimumOutput_ =
true;
226 computeJoinSaddleOutput_ =
true;
227 computeSplitSaddleOutput_ =
true;
228 computeMaximumOutput_ =
true;
238 computeMinimumOutput_ =
true;
239 computeJoinSaddleOutput_ =
true;
240 computeSplitSaddleOutput_ =
true;
241 computeMaximumOutput_ =
true;
244 if(computeMinimumOutput_) {
245 if(outputAllMinimumComponents_)
249 computeMinimumOutput_ =
false;
251 if(computeJoinSaddleOutput_) {
252 if(outputAllJoinSaddleComponents_)
256 computeJoinSaddleOutput_ =
false;
258 if(computeSplitSaddleOutput_) {
259 if(outputAllSplitSaddleComponents_)
263 computeSplitSaddleOutput_ =
false;
265 if(computeMaximumOutput_) {
266 if(outputAllMaximumComponents_)
270 computeMaximumOutput_ =
false;
#define ttkNotUsed(x)
Mark function/method parameters that are not used in the function body at all.
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
ttkMandatoryCriticalPoints()
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
int FillInputPortInformation(int port, vtkInformation *info) override
static void * GetVoidPointer(vtkDataArray *array, vtkIdType start=0)
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_
std::vector< int > mdtSplitTreeEdgeSwitchable_
int computeAllSplitSaddle(const triangulationType &triangulation)
int computeSplitSaddle(const int &id, const triangulationType &triangulation, const bool &reset=true)
int setLowerBoundFieldPointer(void *data)
int computeMaximum(const int &id, const bool &reset=true, const bool &ttkNotUsed(parallel)=true)
int setOutputJoinSaddleDataPointer(void *data)
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 buildJoinTreePlanarLayout()
int setSoSoffsets(int *offsets=nullptr)
int setUpperBoundFieldPointer(void *data)
std::vector< int > mdtSplitTreePointComponentId_
int setOutputSplitSaddleDataPointer(void *data)
std::vector< int > mdtJoinTreeEdgeSwitchable_
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_
int buildSplitTreePlanarLayout()
int setOutputMaximumDataPointer(void *data)
std::vector< double > mdtSplitTreePointUpInterval_
void preconditionTriangulation(AbstractTriangulation *const triangulation)
int setOutputMinimumDataPointer(void *data)
int computeAllJoinSaddle(const triangulationType &triangulation)
#define ttkVtkTemplateMacro(dataType, triangulationType, call)
vtkStandardNewMacro(ttkMandatoryCriticalPoints)
printMsg(debug::output::BOLD+" | | | | | . \\ | | (__| | / __/| |_| / __/|__ _|"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)