TTK
Loading...
Searching...
No Matches
ttkScalarFieldSmoother.cpp
Go to the documentation of this file.
2
3#include <ttkMacros.h>
4#include <ttkUtils.h>
5
6#include <vtkDataArray.h>
7#include <vtkDataSet.h>
8#include <vtkInformation.h>
9#include <vtkObjectFactory.h>
10#include <vtkPointData.h>
11
12using namespace std;
13using namespace ttk;
14
16
18 this->SetNumberOfInputPorts(1);
19 this->SetNumberOfOutputPorts(1);
20}
21
23
25 vtkInformation *info) {
26 if(port == 0) {
27 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkDataSet");
28 return 1;
29 }
30 return 0;
31}
32
34 vtkInformation *info) {
35 if(port == 0) {
37 return 1;
38 }
39 return 0;
40}
41
43 vtkInformationVector **inputVector,
44 vtkInformationVector *outputVector) {
45
46 auto input = vtkDataSet::GetData(inputVector[0]);
47 auto output = vtkDataSet::GetData(outputVector);
48
49 Triangulation *triangulation = ttkAlgorithm::GetTriangulation(input);
50
51 int const keepGoing = checkEmptyMPIInput<Triangulation>(triangulation);
52 if(keepGoing < 2) {
53 return keepGoing;
54 }
55
56 this->preconditionTriangulation(triangulation);
57
58 vtkDataArray *inputScalarField = this->GetInputArrayToProcess(0, inputVector);
59 if(!inputScalarField)
60 return 0;
61
62 if(inputScalarField->GetNumberOfComponents() != 1) {
63 printErr("Invalid scalar field ("
64 + std::to_string(inputScalarField->GetNumberOfComponents())
65 + " components)");
66 return 0;
67 }
68
69 vtkDataArray *inputMaskField = ttkAlgorithm::GetOptionalArray(
70 ForceInputMaskScalarField, 1, ttk::MaskScalarFieldName, input);
71
72 // preparing the output
73 vtkSmartPointer<vtkDataArray> const outputArray
74 = vtkSmartPointer<vtkDataArray>::Take(inputScalarField->NewInstance());
75 outputArray->SetName(inputScalarField->GetName());
76 outputArray->SetNumberOfComponents(1);
77 outputArray->SetNumberOfTuples(inputScalarField->GetNumberOfTuples());
78
79 // This filter copies the input into a new data-set (smoothed)
80 // let's use shallow copies, in order to only duplicate point positions
81 // (before and after). the rest is not changed, pointers are sufficient.
82 output->ShallowCopy(input);
83 output->GetPointData()->AddArray(outputArray);
84
85 printMsg("Starting computation...");
87 {{" Scalar Array", inputScalarField->GetName()},
88 {" Mask Array", inputMaskField ? inputMaskField->GetName() : "None"},
89 {" #iterations", std::to_string(NumberOfIterations)}});
90
91 const auto inputMaskPtr
92 = (inputMaskField) ? ttkUtils::GetPointer<char>(inputMaskField) : nullptr;
93
94 this->setDimensionNumber(inputScalarField->GetNumberOfComponents());
95 this->setInputDataPointer(ttkUtils::GetVoidPointer(inputScalarField));
97 this->setMaskDataPointer(inputMaskPtr);
98
99 // calling the smoothing package
101 inputScalarField->GetDataType(), triangulation->getType(),
102 (smooth<T0, T1>(
103 static_cast<const T1 *>(triangulation->getData()), NumberOfIterations)));
104
105 return 1;
106}
#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 * GetOptionalArray(const bool &enforceArrayIndex, const int &arrayIndex, const std::string &arrayName, vtkDataSet *const inputData, const int &inputPort=0)
TTK VTK-filter for scalar field smoothing.
int FillOutputPortInformation(int port, vtkInformation *info) override
int FillInputPortInformation(int port, vtkInformation *info) override
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
~ttkScalarFieldSmoother() 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 setMaskDataPointer(const char *const mask)
int preconditionTriangulation(AbstractTriangulation *triangulation)
void setDimensionNumber(const int &dimensionNumber)
void setOutputDataPointer(void *data)
void setInputDataPointer(void *data)
Triangulation is a class that provides time and memory efficient traversal methods on triangulations ...
AbstractTriangulation * getData()
Triangulation::Type getType() const
std::string to_string(__int128)
Definition ripserpy.cpp:99
The Topology ToolKit.
const char MaskScalarFieldName[]
default name for mask scalar field
Definition DataTypes.h:32
#define ttkTypeMacroAT(group0, group1, call)
Definition ttkMacros.h:217
vtkStandardNewMacro(ttkScalarFieldSmoother)
printMsg(debug::output::BOLD+" | | | | | . \\ | | (__| | / __/| |_| / __/|__ _|"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)