TTK
Loading...
Searching...
No Matches
ttkMorphologicalOperators.cpp
Go to the documentation of this file.
2
3#include <vtkDataArray.h>
4#include <vtkDataSet.h>
5#include <vtkInformation.h>
6#include <vtkObjectFactory.h>
7#include <vtkPointData.h>
8#include <vtkSmartPointer.h>
9
10#include <ttkMacros.h>
11#include <ttkUtils.h>
12
14
16 this->SetNumberOfInputPorts(1);
17 this->SetNumberOfOutputPorts(1);
18}
19
21
23 vtkInformation *info) {
24 if(port == 0) {
25 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkDataSet");
26 return 1;
27 }
28 return 0;
29}
30
32 vtkInformation *info) {
33 if(port == 0) {
35 return 1;
36 }
37 return 0;
38}
39
41 vtkInformationVector **inputVector,
42 vtkInformationVector *outputVector) {
43 // get input and output
44 auto input = vtkDataSet::GetData(inputVector[0]);
45 auto output = vtkDataSet::GetData(outputVector);
46 if(!input || !output)
47 return 0;
48
49 // copy input to output
50 output->ShallowCopy(input);
51
52 // get input label array
53 auto inputLabels = this->GetInputArrayToProcess(0, inputVector);
54 if(!inputLabels)
55 return 0;
56 if(this->GetInputArrayAssociation(0, inputVector) != 0) {
57 this->printErr("Input labels needs to be a point data array.");
58 return 0;
59 }
60 if(inputLabels->GetNumberOfComponents() != 1) {
61 this->printErr("Input labels needs to be a scalar array.");
62 return 0;
63 }
64
65 // -------------------------------------------------------------------------
66 // Replace Variables in PivotLabel (e.g. {time[2]})
67 // -------------------------------------------------------------------------
68 double pivotLabel = 0;
69 {
70 std::string temp, errorMsg;
72 this->GetPivotLabel(), input->GetFieldData(), temp, errorMsg)) {
73 this->printErr(errorMsg);
74 return 0;
75 }
76
77 std::vector<double> values;
79
80 if(values.size() < 1) {
81 this->printErr("Unable to parse pivot label as double.");
82 return 0;
83 }
84 pivotLabel = values[0];
85 }
86
87 // create output labels
88 auto outputLabels
89 = vtkSmartPointer<vtkDataArray>::Take(inputLabels->NewInstance());
90 outputLabels->DeepCopy(inputLabels);
91 output->GetPointData()->AddArray(outputLabels);
92
93 // get triangulation
94 auto triangulation = ttkAlgorithm::GetTriangulation(input);
95 if(!triangulation)
96 return 0;
97
98 // precondition triangulation
99 this->preconditionTriangulation(triangulation);
100
101 int status = 0;
102
104 inputLabels->GetDataType(), triangulation->getType(),
105 (status = this->execute<VTK_TT, TTK_TT>(
106 // Output
107 static_cast<VTK_TT *>(ttkUtils::GetVoidPointer(outputLabels)),
108
109 // Input
110 this->Mode, this->Iterations, this->Grayscale,
111 static_cast<VTK_TT *>(ttkUtils::GetVoidPointer(inputLabels)), pivotLabel,
112 static_cast<TTK_TT *>(triangulation->getData()))));
113
114 if(!status)
115 return 0;
116
117 return 1;
118}
#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 that dilates or erodes a specified vertex label.
virtual std::string GetPivotLabel()
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
int FillInputPortInformation(int port, vtkInformation *info) override
~ttkMorphologicalOperators() override
int FillOutputPortInformation(int port, vtkInformation *info) override
static int replaceVariables(const std::string &iString, vtkFieldData *fieldData, std::string &oString, std::string &errorMsg)
Definition ttkUtils.cpp:73
static void * GetVoidPointer(vtkDataArray *array, vtkIdType start=0)
Definition ttkUtils.cpp:226
static int stringListToDoubleVector(const std::string &iString, std::vector< double > &v)
Definition ttkUtils.cpp:133
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(ttk::AbstractTriangulation *triangulation) const
#define ttkVtkTemplateMacro(dataType, triangulationType, call)
Definition ttkMacros.h:69
vtkStandardNewMacro(ttkMorphologicalOperators)