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#include <vtkUnstructuredGrid.h>
11
12#include <ttkMacros.h>
14#include <ttkUtils.h>
15
17
19 this->SetNumberOfInputPorts(2);
20 this->SetNumberOfOutputPorts(1);
21}
22
24 int port, vtkInformation *info) {
25 if(port == 0) {
26 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkDataSet");
27 return 1;
28 } else if(port == 1) {
29 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkUnstructuredGrid");
30 return 1;
31 }
32 return 0;
33}
34
36 int port, vtkInformation *info) {
37 if(port == 0) {
39 return 1;
40 }
41 return 0;
42}
43
45 vtkInformation *ttkNotUsed(request),
46 vtkInformationVector **inputVector,
47 vtkInformationVector *outputVector) {
48
49 using ttk::SimplexId;
50
51 if(this->Method == 0) {
53 } else if(this->Method == 1) {
55 } else if(this->Method == 2) {
57 }
58
59 const auto domain = vtkDataSet::GetData(inputVector[0]);
60 const auto constraints = vtkUnstructuredGrid::GetData(inputVector[1]);
61 if(!domain || !constraints)
62 return !this->printErr("Unable to retrieve required input data objects.");
63
64 auto output = vtkDataSet::GetData(outputVector);
65
66 // triangulation
67 auto triangulation = ttkAlgorithm::GetTriangulation(domain);
68
69 if(!triangulation) {
70 this->printErr("Input triangulation pointer is NULL.");
71 return -1;
72 }
73
74 this->preconditionTriangulation(triangulation);
75
76 if(triangulation->isEmpty()) {
77 this->printErr("Triangulation allocation problem.");
78 return -1;
79 }
80
81 if(!domain) {
82 this->printErr("Input pointer is NULL.");
83 return -1;
84 }
85
86 const auto numberOfVertices = domain->GetNumberOfPoints();
87 if(numberOfVertices <= 0) {
88 this->printErr("Domain has no points.");
89 return -5;
90 }
91
92 // domain scalar field
93 const auto inputScalars = this->GetInputArrayToProcess(0, domain);
94 if(!inputScalars) {
95 this->printErr("Input scalar field pointer is null.");
96 return -3;
97 }
98
99 // domain offset field
100 const auto inputOrder = this->GetOrderArray(
101 domain, 0, triangulation, false, 2, ForceInputOffsetScalarField);
102 if(!inputOrder) {
103 this->printErr("Wrong input offset scalar field.");
104 return -1;
105 }
106
107 // Constraints
108 ttk::DiagramType constraintDiagram;
109 const ttk::Debug dbg;
110 VTUToDiagram(constraintDiagram, constraints, dbg);
111
112 // create output arrays
113 auto outputScalars
114 = vtkSmartPointer<vtkDataArray>::Take(inputScalars->NewInstance());
115 outputScalars->DeepCopy(inputScalars);
116 auto outputOrder
117 = vtkSmartPointer<vtkDataArray>::Take(inputOrder->NewInstance());
118 outputOrder->DeepCopy(inputOrder);
119
120 // constraint identifier field
121 int const numberOfConstraints = constraints->GetNumberOfPoints();
122
123 std::vector<ttk::SimplexId> idSpareStorage{};
124 auto identifiers = this->GetIdentifierArrayPtr(ForceInputVertexScalarField, 1,
126 constraints, idSpareStorage);
127
128 int ret{};
129 switch(inputScalars->GetDataType()) {
130 vtkTemplateMacro(ret = this->execute(
131 ttkUtils::GetPointer<VTK_TT>(inputScalars),
132 ttkUtils::GetPointer<VTK_TT>(outputScalars), identifiers,
133 ttkUtils::GetPointer<SimplexId>(inputOrder),
134 ttkUtils::GetPointer<SimplexId>(outputOrder),
135 numberOfConstraints, this->AddPerturbation,
136 *triangulation->getData(), constraintDiagram));
137 }
138
139 // something wrong in baseCode
140 if(ret) {
141 this->printErr("TopologicalSimplification.execute() error code: "
142 + std::to_string(ret));
143 return -12;
144 }
145
146 output->ShallowCopy(domain);
147 output->GetPointData()->AddArray(outputOrder);
148 output->GetPointData()->AddArray(outputScalars);
149
150 return 1;
151}
#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)
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 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
Minimalist debugging class.
Definition Debug.h:88
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(AbstractTriangulation *triangulation)
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, triangulationType &triangulation, const ttk::DiagramType &constraintDiagram={})
std::vector< PersistencePair > DiagramType
Persistence Diagram type as a vector of Persistence pairs.
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
int VTUToDiagram(ttk::DiagramType &diagram, vtkUnstructuredGrid *vtu, const ttk::Debug &dbg)
Converts a Persistence Diagram in the VTK Unstructured Grid format (as generated by the ttkPersistenc...
vtkStandardNewMacro(ttkTopologicalSimplification)