TTK
Loading...
Searching...
No Matches
ttkTopologicalCompression.cpp
Go to the documentation of this file.
2#include <ttkMacros.h>
3#include <ttkUtils.h>
4
5#include <vtkDataSet.h>
6#include <vtkDoubleArray.h>
7#include <vtkFloatArray.h>
8#include <vtkIdTypeArray.h>
9#include <vtkInformation.h>
10#include <vtkIntArray.h>
11#include <vtkNew.h>
12#include <vtkPointData.h>
13#include <vtkSignedCharArray.h>
14#include <vtkSmartPointer.h>
15
17
19 this->SetNumberOfInputPorts(1);
20 this->SetNumberOfOutputPorts(1);
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 vtkInformation *info) {
34 if(port == 0) {
36 return 1;
37 }
38 return 0;
39}
40
42 vtkInformationVector **inputVector,
43 vtkInformationVector *outputVector) {
44
45 auto input = vtkDataSet::GetData(inputVector[0]);
46 auto output = vtkDataSet::GetData(outputVector);
47
48#ifndef TTK_ENABLE_KAMIKAZE
49 if(input == nullptr) {
50 this->printErr("Input pointer is NULL.");
51 return -1;
52 }
53 if(input->GetNumberOfPoints() == 0) {
54 this->printErr("Input has no point.");
55 return -1;
56 }
57 if(input->GetPointData() == nullptr) {
58 this->printErr("Input has no point data.");
59 return -1;
60 }
61 if(output == nullptr) {
62 this->printErr("Output pointer is NULL.");
63 return -1;
64 }
65#endif
66
67 // Triangulate
68 auto triangulation = ttkAlgorithm::GetTriangulation(input);
69 if(triangulation == nullptr) {
70 return 0;
71 }
72 this->preconditionTriangulation(triangulation);
73
74 // use a pointer-base copy for the input data -- to adapt if your wrapper does
75 // not produce an output of the type of the input.
76 output->ShallowCopy(input);
77
78 // in the following, the target scalar field of the input is replaced in the
79 // variable 'output' with the result of the computation.
80 // if your wrapper produces an output of the same type of the input, you
81 // should proceed in the same way.
82 const auto inputScalarField = this->GetInputArrayToProcess(0, inputVector);
83
84#ifndef TTK_ENABLE_KAMIKAZE
85 if(inputScalarField == nullptr) {
86 this->printErr("Input scalar field pointer is NULL.");
87 return -1;
88 }
89#endif
90
91 const auto vertexNumber = inputScalarField->GetNumberOfTuples();
92
93 const auto inputOffsets
94 = this->GetOrderArray(input, 0, triangulation, false, 1, false);
95
96 // allocate the memory for the output scalar field
97 vtkSmartPointer<vtkDataArray> outputScalarField{};
98
99 switch(inputScalarField->GetDataType()) {
100 case VTK_CHAR:
101 outputScalarField = vtkSmartPointer<vtkSignedCharArray>::New();
102 break;
103 case VTK_DOUBLE:
104 outputScalarField = vtkSmartPointer<vtkDoubleArray>::New();
105 break;
106 case VTK_FLOAT:
107 outputScalarField = vtkSmartPointer<vtkFloatArray>::New();
108 break;
109 case VTK_INT:
110 outputScalarField = vtkSmartPointer<vtkIntArray>::New();
111 break;
112 case VTK_ID_TYPE:
113 outputScalarField = vtkSmartPointer<vtkIdTypeArray>::New();
114 break;
115 default:
116 this->printErr("Unsupported data type :(");
117 return -1;
118 }
119
120 outputScalarField->SetNumberOfTuples(vertexNumber);
121 outputScalarField->SetName(inputScalarField->GetName());
122
123 vtkNew<vtkIntArray> outputOffsetField{};
124 outputOffsetField->SetNumberOfTuples(vertexNumber);
125 outputOffsetField->SetName(this->GetOrderArrayName(inputScalarField).data());
126
127 // manage tolerance (relative % -> absolute)
128 std::array<double, 2> sfRange{};
129 inputScalarField->GetRange(sfRange.data());
130 this->relToAbsZFPTolerance(this->ZFPTolerance, sfRange);
131
132 // Call TopologicalCompression
134 inputScalarField->GetDataType(), triangulation->getType(),
135 this->execute(
136 static_cast<VTK_TT *>(ttkUtils::GetVoidPointer(inputScalarField)),
137 static_cast<ttk::SimplexId *>(ttkUtils::GetVoidPointer(inputOffsets)),
138 static_cast<VTK_TT *>(ttkUtils::GetVoidPointer(outputScalarField)),
139 *static_cast<TTK_TT *>(triangulation->getData())));
140
141 for(ttk::SimplexId i = 0; i < vertexNumber; ++i)
142 outputOffsetField->SetTuple1(i, this->compressedOffsets_[i]);
143
144 output->GetPointData()->AddArray(outputScalarField);
145 output->GetPointData()->AddArray(outputOffsetField);
146
147 return 1;
148}
#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)
vtkDataArray * GetOrderArray(vtkDataSet *const inputData, const int scalarArrayIdx, ttk::Triangulation *triangulation, const bool getGlobalOrder=false, const int orderArrayIdx=0, const bool enforceOrderArrayIdx=false)
static std::string GetOrderArrayName(vtkDataArray *const array)
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
int FillOutputPortInformation(int port, vtkInformation *info) override
int FillInputPortInformation(int port, vtkInformation *info) override
static void * GetVoidPointer(vtkDataArray *array, vtkIdType start=0)
Definition ttkUtils.cpp:226
int printErr(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition Debug.h:149
void relToAbsZFPTolerance(const double zfpRelTol, const std::array< double, 2 > &sfRange)
Switch from a relative (% of scalar field range) to an absolute ZFP tolerance.
void preconditionTriangulation(AbstractTriangulation *const triangulation)
int SimplexId
Identifier type for simplices of any dimension.
Definition DataTypes.h:22
#define ttkVtkTemplateMacro(dataType, triangulationType, call)
Definition ttkMacros.h:69
vtkStandardNewMacro(ttkTopologicalCompression)