TTK
Loading...
Searching...
No Matches
ttkScalarFieldNormalizer.cpp
Go to the documentation of this file.
2
3#include <Geometry.h>
4
5#include <vtkDataArray.h>
6#include <vtkDataSet.h>
7#include <vtkInformation.h>
8#include <vtkObjectFactory.h>
9#include <vtkPointData.h>
10
11#include <ttkMacros.h>
12#include <ttkUtils.h>
13
14using namespace std;
15using namespace ttk;
16
18
20 // init
21 this->SetNumberOfInputPorts(1);
22 this->SetNumberOfOutputPorts(1);
23 setDebugMsgPrefix("ScalarFieldNormalizer");
24#ifdef TTK_ENABLE_MPI
25 hasMPISupport_ = true;
26#endif
27}
28
30
32 vtkInformation *info) {
33 if(port == 0) {
34 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkDataSet");
35 return 1;
36 }
37 return 0;
38}
39
41 vtkInformation *info) {
42 if(port == 0) {
44 return 1;
45 }
46 return 0;
47}
48
50 vtkDataArray *output) const {
51
52 if(!output)
53 return -1;
54 if(!input)
55 return -2;
56
57 double min = 0, max = 0;
58 for(SimplexId i = 0; i < input->GetNumberOfTuples(); i++) {
59
60 double const value = input->GetTuple1(i);
61
62 if((!i) || (value < min)) {
63 min = value;
64 }
65 if((!i) || (value > max)) {
66 max = value;
67 }
68 }
69
70#ifdef TTK_ENABLE_MPI
71 if(ttk::isRunningWithMPI()) {
72 // if we are built with MPI and are running with MPI, we need to calculate
73 // the min and max over all ranks. we do this by using MPI_Allreduce
74 MPI_Allreduce(MPI_IN_PLACE, &min, 1, MPI_DOUBLE, MPI_MIN, ttk::MPIcomm_);
75 MPI_Allreduce(MPI_IN_PLACE, &max, 1, MPI_DOUBLE, MPI_MAX, ttk::MPIcomm_);
76 }
77#endif
78
79 for(SimplexId i = 0; i < input->GetNumberOfTuples(); i++) {
80 double value = input->GetTuple1(i);
81
82 value = (value - min) / (max - min) + Geometry::powIntTen(-FLT_DIG);
83
84 output->SetTuple1(i, value);
85 }
86
87 return 0;
88}
89
91 vtkInformationVector **inputVector,
92 vtkInformationVector *outputVector) {
93
94 vtkDataSet *input = vtkDataSet::GetData(inputVector[0]);
95 vtkDataSet *output = vtkDataSet::GetData(outputVector);
96
97 // get input scalar field
98 vtkDataArray *inputArray = this->GetInputArrayToProcess(0, inputVector);
99
100 int const keepGoing
101 = ttkAlgorithm::checkEmptyMPIInput<vtkDataArray>(inputArray);
102 if(keepGoing < 2) {
103 return keepGoing;
104 }
105
106 vtkSmartPointer<vtkDataArray> const outputArray
107 = vtkSmartPointer<vtkDataArray>::Take(inputArray->NewInstance());
108 outputArray->SetName(inputArray->GetName());
109 outputArray->SetNumberOfComponents(1);
110 outputArray->SetNumberOfTuples(inputArray->GetNumberOfTuples());
111
112 // calling the executing package
113 normalize(inputArray, outputArray);
114
115 // prepare the output
116 output->ShallowCopy(input);
117 output->GetPointData()->AddArray(outputArray);
118
119 return 1;
120}
#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 normalizes an input scalar field.
int normalize(vtkDataArray *input, vtkDataArray *output) const
~ttkScalarFieldNormalizer() override
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
int FillOutputPortInformation(int port, vtkInformation *info) override
int FillInputPortInformation(int port, vtkInformation *info) override
void setDebugMsgPrefix(const std::string &prefix)
Definition Debug.h:364
T powIntTen(const int n)
Compute the nth power of ten.
Definition Geometry.h:403
The Topology ToolKit.
int SimplexId
Identifier type for simplices of any dimension.
Definition DataTypes.h:22
vtkStandardNewMacro(ttkScalarFieldNormalizer)