TTK
Loading...
Searching...
No Matches
ttkDistanceField.cpp
Go to the documentation of this file.
1#include <vtkDoubleArray.h>
2#include <vtkFloatArray.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
10#include <ttkDistanceField.h>
11#include <ttkMacros.h>
12#include <ttkUtils.h>
13
15
17 this->SetNumberOfInputPorts(2);
18 this->SetNumberOfOutputPorts(1);
19}
20
21int ttkDistanceField::FillInputPortInformation(int port, vtkInformation *info) {
22 if(port == 0) {
23 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkDataSet");
24 return 1;
25 } else if(port == 1) {
26 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkPointSet");
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
41int ttkDistanceField::RequestData(vtkInformation *ttkNotUsed(request),
42 vtkInformationVector **inputVector,
43 vtkInformationVector *outputVector) {
44 ttk::Timer const globalTimer;
45
46 vtkDataSet *domain = vtkDataSet::GetData(inputVector[0]);
47 vtkPointSet *sources = vtkPointSet::GetData(inputVector[1]);
48 vtkDataSet *output = vtkDataSet::GetData(outputVector);
49
50 auto triangulation = ttkAlgorithm::GetTriangulation(domain);
51 if(triangulation == nullptr) {
52 this->printErr("Wrong triangulation.");
53 return 0;
54 }
55
56 this->preconditionTriangulation(triangulation);
57
58 std::vector<ttk::SimplexId> idSpareStorage{};
59 auto *identifiers = this->GetIdentifierArrayPtr(ForceInputVertexScalarField,
61 sources, idSpareStorage);
62 if(identifiers == nullptr) {
63 printErr("Wrong identifiers.");
64 return 0;
65 }
66
67 const int numberOfPointsInDomain = domain->GetNumberOfPoints();
68 if(numberOfPointsInDomain == 0) {
69 printErr("Domain has no points.");
70 return 0;
71 }
72
73 const int numberOfPointsInSources = sources->GetNumberOfPoints();
74 if(numberOfPointsInSources == 0) {
75 printErr("Sources have no points.");
76 return 0;
77 }
78
79 vtkNew<ttkSimplexIdTypeArray> origin{};
80 origin->SetNumberOfComponents(1);
81 origin->SetNumberOfTuples(numberOfPointsInDomain);
82 origin->SetName(ttk::VertexScalarFieldName);
83
84 vtkNew<ttkSimplexIdTypeArray> seg{};
85 seg->SetNumberOfComponents(1);
86 seg->SetNumberOfTuples(numberOfPointsInDomain);
87 seg->SetName("SeedIdentifier");
88
89 this->setVertexNumber(numberOfPointsInDomain);
90 this->setSourceNumber(numberOfPointsInSources);
94
95 vtkSmartPointer<vtkDataArray> distanceScalars{};
96 int ret{};
97 switch(static_cast<DistanceType>(OutputScalarFieldType)) {
98 case DistanceType::Float:
99 distanceScalars = vtkFloatArray::New();
100 distanceScalars->SetNumberOfComponents(1);
101 distanceScalars->SetNumberOfTuples(numberOfPointsInDomain);
102 distanceScalars->SetName(OutputScalarFieldName.data());
103
105 ttkUtils::GetVoidPointer(distanceScalars));
107 triangulation->getType(), (ret = this->execute<float, TTK_TT>(
108 (TTK_TT *)triangulation->getData())));
109 break;
110
111 case DistanceType::Double:
112 distanceScalars = vtkDoubleArray::New();
113 distanceScalars->SetNumberOfComponents(1);
114 distanceScalars->SetNumberOfTuples(numberOfPointsInDomain);
115 distanceScalars->SetName(OutputScalarFieldName.data());
116
118 ttkUtils::GetVoidPointer(distanceScalars));
119
121 triangulation->getType(), (ret = this->execute<double, TTK_TT>(
122 (TTK_TT *)triangulation->getData())));
123 break;
124
125 default:
126 printErr("Invalid scalar field type.");
127 return 0;
128 break;
129 }
130
131 // something wrong in baseCode
132 if(ret != 0) {
133 printErr("DistanceField.execute() error code : " + std::to_string(ret));
134 return 0;
135 }
136
137 // update result
138 output->ShallowCopy(domain);
139 output->GetPointData()->AddArray(distanceScalars);
140 output->GetPointData()->AddArray(origin);
141 output->GetPointData()->AddArray(seg);
142 distanceScalars->Delete();
143
144 return 1;
145}
#define ttkTemplateMacro(triangulationType, call)
#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)
TTK VTK-filter for distance field computations.
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
int FillInputPortInformation(int port, vtkInformation *info) override
int FillOutputPortInformation(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 setOutputIdentifiers(void *data)
void setVertexNumber(SimplexId vertexNumber)
void setVertexIdentifierScalarFieldPointer(SimplexId *const data)
void setOutputScalarFieldPointer(void *data)
void setOutputSegmentation(void *data)
int preconditionTriangulation(AbstractTriangulation *triangulation)
void setSourceNumber(SimplexId sourceNumber)
const char VertexScalarFieldName[]
default name for vertex scalar field
Definition DataTypes.h:35
vtkStandardNewMacro(ttkDistanceField)