TTK
Loading...
Searching...
No Matches
ttkLDistance.cpp
Go to the documentation of this file.
1#include <vtkDataArray.h>
2#include <vtkDataSet.h>
3#include <vtkDoubleArray.h>
4#include <vtkInformation.h>
5#include <vtkObjectFactory.h>
6#include <vtkPointData.h>
7#include <vtkSmartPointer.h>
8
9#include "ttkLDistance.h"
10#include <ttkUtils.h>
11
13
15 this->SetNumberOfInputPorts(1);
16 this->SetNumberOfOutputPorts(1);
17}
18
19int ttkLDistance::FillInputPortInformation(int port, vtkInformation *info) {
20 if(port == 0) {
21 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkDataSet");
22 return 1;
23 }
24 return 0;
25}
26
27int ttkLDistance::FillOutputPortInformation(int port, vtkInformation *info) {
28 if(port == 0) {
30 return 1;
31 }
32 return 0;
33}
34
35int ttkLDistance::RequestData(vtkInformation *ttkNotUsed(request),
36 vtkInformationVector **inputVector,
37 vtkInformationVector *outputVector) {
38
39 const auto input = vtkDataSet::GetData(inputVector[0]);
40 auto output = vtkDataSet::GetData(outputVector);
41
42 // Test validity of datasets (must present the same number of points).
43#ifndef TTK_ENABLE_KAMIKAZE
44 if(!input) {
45 this->printErr("Input pointer is NULL.");
46 return 0;
47 }
48 if(!input->GetNumberOfPoints()) {
49 this->printErr("Input has no point.");
50 return 0;
51 }
52 if(!input->GetPointData()) {
53 this->printErr("Input has no point data.");
54 return 0;
55 }
56 if(!output) {
57 this->printErr("Output pointer is NULL.");
58 return 0;
59 }
60#endif
61
62 // Use a pointer-base copy for the input.
63 output->ShallowCopy(input);
64
65 const auto inputScalarField1 = this->GetInputArrayToProcess(0, input);
66 const auto inputScalarField2 = this->GetInputArrayToProcess(1, input);
67
68#ifndef TTK_ENABLE_KAMIKAZE
69 if(!inputScalarField1 || !inputScalarField2
70 || inputScalarField1->GetDataType() != inputScalarField2->GetDataType()) {
71 this->printErr("Input scalar fields are NULL or have different types.");
72 return 0;
73 }
74#endif
75
76 // Allocate memory for the output scalar field, based on the first input.
77 vtkSmartPointer<vtkDataArray> const outputScalarField{
78 inputScalarField1->NewInstance()};
79
80 ttk::SimplexId const numberOfPoints = input->GetNumberOfPoints();
81
82 outputScalarField->SetNumberOfTuples(numberOfPoints);
83 outputScalarField->SetName(DistanceFieldName.data());
84 output->GetPointData()->AddArray(outputScalarField);
85
86 // Calling the executing package.
87 switch(inputScalarField1->GetDataType()) {
88 vtkTemplateMacro(this->execute<VTK_TT>(
89 static_cast<VTK_TT *>(ttkUtils::GetVoidPointer(inputScalarField1)),
90 static_cast<VTK_TT *>(ttkUtils::GetVoidPointer(inputScalarField2)),
91 static_cast<VTK_TT *>(ttkUtils::GetVoidPointer(outputScalarField)),
92 DistanceType, numberOfPoints));
93 }
94
95 vtkNew<vtkDoubleArray> meanDistanceArray{};
96 const std::string arrayName = "L" + DistanceType + "-distance";
97 meanDistanceArray->SetName(arrayName.c_str());
98 meanDistanceArray->SetNumberOfTuples(1);
99 meanDistanceArray->SetTuple1(0, this->getResult());
100 output->GetFieldData()->AddArray(meanDistanceArray);
101
102 return 1;
103}
#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 VTK-filter that wraps the lDistance processing package.
int FillOutputPortInformation(int port, vtkInformation *info) override
int FillInputPortInformation(int port, vtkInformation *info) override
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) 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
double getResult()
Definition LDistance.h:51
int SimplexId
Identifier type for simplices of any dimension.
Definition DataTypes.h:22
vtkStandardNewMacro(ttkLDistance)