TTK
Loading...
Searching...
No Matches
ttkIdentifyByScalarField.cpp
Go to the documentation of this file.
1#include <numeric>
2
4#include <ttkMacros.h>
5#include <ttkUtils.h>
6
7#include <vtkCellData.h>
8#include <vtkDataSet.h>
9#include <vtkIdTypeArray.h>
10#include <vtkInformation.h>
11#include <vtkIntArray.h>
12#include <vtkPointData.h>
13
14using namespace std;
15using namespace ttk;
16
18
20 this->setDebugMsgPrefix("IdentifyByScalarField");
21 SetNumberOfInputPorts(1);
22 SetNumberOfOutputPorts(1);
23}
24
26 vtkInformation *info) {
27 if(port == 0) {
28 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkDataSet");
29 return 1;
30 }
31
32 return 0;
33}
34
36 vtkInformation *info) {
37 if(port == 0) {
39 return 1;
40 }
41
42 return 0;
43}
44
45template <typename VTK_TT>
46int ttkIdentifyByScalarField::dispatch(vector<SimplexId> &inputIds) {
47 VTK_TT *arr = static_cast<VTK_TT *>(ttkUtils::GetVoidPointer(inputScalars_));
48
49 auto greater_cmp = [=](int a, int b) { return arr[a] > arr[b]; };
50 auto lower_cmp = [=](int a, int b) { return arr[a] < arr[b]; };
51
52 if(IncreasingOrder) {
53 std::sort(inputIds.begin(), inputIds.end(), lower_cmp);
54 } else {
55 std::sort(inputIds.begin(), inputIds.end(), greater_cmp);
56 }
57
58 return 0;
59}
60
62 vtkInformationVector **inputVector,
63 vtkInformationVector *outputVector) {
64 vtkDataSet *input = vtkDataSet::GetData(inputVector[0]);
65 vtkDataSet *output = vtkDataSet::GetData(outputVector);
66
67 inputScalars_ = this->GetInputArrayToProcess(0, inputVector);
68 if(!inputScalars_)
69 return 0;
70
71 int inputArrayAssociation = this->GetInputArrayAssociation(0, inputVector);
72
73 if(inputArrayAssociation > 1 || inputArrayAssociation < 0) {
74 printErr("input array has to be cell data or point data.");
75 return 0;
76 }
77
78 ttk::Timer t;
79 this->printMsg("Computing Identifiers", 0, t.getElapsedTime(), 1,
81
82 const SimplexId numberOfValues = inputArrayAssociation == 0
83 ? input->GetNumberOfPoints()
84 : input->GetNumberOfCells();
85
86 vector<SimplexId> inputIds(numberOfValues);
87 std::iota(inputIds.begin(), inputIds.end(), 0);
88 switch(inputScalars_->GetDataType()) {
89 vtkTemplateMacro(dispatch<VTK_TT>(inputIds));
90 }
91
92 this->printMsg("Computing Identifiers", 1, t.getElapsedTime(), 1);
93
96 ids->SetNumberOfComponents(1);
97 ids->SetNumberOfTuples(numberOfValues);
98 ids->SetName(inputArrayAssociation == 0 ? "PointScalarFieldName"
99 : "CellScalarFieldName");
100
101 SimplexId *outputIds
102 = static_cast<SimplexId *>(ttkUtils::GetVoidPointer(ids));
103
104 for(int i = 0; i < numberOfValues; ++i)
105 outputIds[inputIds[i]] = i;
106 if(StartByOne) {
107 for(int i = 0; i < numberOfValues; ++i)
108 outputIds[i] += 1;
109 }
110
111 output->ShallowCopy(input);
112
113 if(inputArrayAssociation == 0) {
114 output->GetPointData()->AddArray(ids);
115 } else {
116 output->GetCellData()->AddArray(ids);
117 }
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 computes a new scalar array based on a sorting of the input array.
int dispatch(std::vector< ttk::SimplexId > &inputIds)
int FillInputPortInformation(int port, vtkInformation *info) override
int FillOutputPortInformation(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:225
void setDebugMsgPrefix(const std::string &prefix)
Definition: Debug.h:364
int printMsg(const std::string &msg, const debug::Priority &priority=debug::Priority::INFO, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cout) const
Definition: Debug.h:118
int printErr(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition: Debug.h:149
double getElapsedTime()
Definition: Timer.h:15
The Topology ToolKit.
int SimplexId
Identifier type for simplices of any dimension.
Definition: DataTypes.h:22
vtkStandardNewMacro(ttkIdentifyByScalarField)