TTK
Loading...
Searching...
No Matches
ttkPointSetToCurve.cpp
Go to the documentation of this file.
2
3#include <vtkPointData.h>
4#include <vtkPointSet.h>
5#include <vtkUnstructuredGrid.h>
6
7#include <vtkInformation.h>
8#include <vtkInformationVector.h>
9
10#include <ttkUtils.h>
11
12#include <array>
13#include <map>
14#include <set>
15
17
19 this->setDebugMsgPrefix("PointSetToCurve");
20
21 this->SetNumberOfInputPorts(1);
22 this->SetNumberOfOutputPorts(1);
23}
24
26 vtkInformation *info) {
27 if(port == 0) {
28 info->Set(vtkDataObject::DATA_TYPE_NAME(), "vtkPointSet");
29 return 1;
30 }
31 return 0;
32}
33
35 vtkInformation *info) {
36 if(port == 0) {
37 info->Set(vtkDataObject::DATA_TYPE_NAME(), "vtkUnstructuredGrid");
38 return 1;
39 }
40 return 0;
41}
42
43template <typename VTK_TT>
45 std::vector<std::pair<vtkIdType, double>> &storage,
46 const VTK_TT *const values,
47 const size_t nvalues) {
48
49 for(size_t i = 0; i < nvalues; ++i) {
50 storage.emplace_back(i, static_cast<double>(values[i]));
51 }
52}
53
54int ttkPointSetToCurve::RequestData(vtkInformation *ttkNotUsed(request),
55 vtkInformationVector **inputVector,
56 vtkInformationVector *outputVector) {
57 const auto input = vtkPointSet::GetData(inputVector[0]);
58 auto output = vtkUnstructuredGrid::GetData(outputVector);
59
60 if(input == nullptr || output == nullptr) {
61 this->printErr("Null input data, aborting");
62 return 0;
63 }
64
65 // ordering array
66 const auto oa = this->GetInputArrayToProcess(0, inputVector);
67
68 if(oa == nullptr) {
69 this->printErr("Cannot find the required data array");
70 return 0;
71 }
72
73 const auto nvalues = oa->GetNumberOfTuples();
74
75 // store point index <-> ordering value in vector
76 std::vector<std::pair<vtkIdType, double>> orderedValues{};
77
78 switch(oa->GetDataType()) {
79 vtkTemplateMacro(
80 dispatch(orderedValues,
81 static_cast<VTK_TT *>(ttkUtils::GetVoidPointer(oa)), nvalues));
82 }
83
84 // compare two pairs of index/value according to their values
85 const auto cmp
86 = [](const std::pair<vtkIdType, double> &a,
87 const std::pair<vtkIdType, double> &b) { return a.second < b.second; };
88
89 // sort the vector of indices/values in ascending order
91 this->threadNumber_, orderedValues.begin(), orderedValues.end(), cmp);
92
93 // shallow-copy input into output
94 output->ShallowCopy(input);
95
96 // output cell array
97 vtkNew<vtkCellArray> cells{};
98
99 for(size_t i = 1; i < orderedValues.size(); ++i) {
100 std::array<vtkIdType, 2> linePoints{
101 orderedValues[i - 1].first, orderedValues[i].first};
102 cells->InsertNextCell(2, linePoints.data());
103 }
104
105 if(CloseCurve) {
106 std::array<vtkIdType, 2> linePoints{
107 orderedValues.back().first, orderedValues.front().first};
108 cells->InsertNextCell(2, linePoints.data());
109 }
110
111 output->SetCells(VTK_LINE, cells);
112
113 return 1;
114}
#define ttkNotUsed(x)
Mark function/method parameters that are not used in the function body at all.
Definition BaseClass.h:47
#define TTK_PSORT(NTHREADS,...)
Parallel sort macro.
Definition OpenMP.h:46
TTK VTK-filter that reads a Cinema Spec D Database.
int FillInputPortInformation(int port, vtkInformation *info) override
void dispatch(std::vector< std::pair< vtkIdType, double > > &storage, const VTK_TT *const values, const size_t nvalues)
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:226
void setDebugMsgPrefix(const std::string &prefix)
Definition Debug.h:364
int printErr(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition Debug.h:149
vtkStandardNewMacro(ttkPointSetToCurve)