TTK
Loading...
Searching...
No Matches
ttkLDistanceMatrix.cpp
Go to the documentation of this file.
2#include <ttkUtils.h>
3
4#include <vtkDataArray.h>
5#include <vtkDataSet.h>
6#include <vtkDoubleArray.h>
7#include <vtkInformation.h>
8#include <vtkInformationVector.h>
9#include <vtkMultiBlockDataSet.h>
10#include <vtkNew.h>
11#include <vtkObjectFactory.h>
12#include <vtkPointData.h>
13#include <vtkTable.h>
14
15#include <set>
16
18
20 SetNumberOfInputPorts(1);
21 SetNumberOfOutputPorts(1);
22}
23
25 vtkInformation *info) {
26 if(port == 0) {
27 info->Set(vtkDataObject::DATA_TYPE_NAME(), "vtkMultiBlockDataSet");
28 return 1;
29 }
30 return 0;
31}
32
34 vtkInformation *info) {
35 if(port == 0) {
36 info->Set(vtkDataObject::DATA_TYPE_NAME(), "vtkTable");
37 return 1;
38 }
39 return 0;
40}
41
42template <typename T>
44 std::vector<std::vector<double>> &distanceMatrix,
45 const std::vector<vtkDataSet *> &inputData,
46 const size_t nPoints) {
47
48 std::vector<const T *> inputPtrs(inputData.size());
49 for(size_t i = 0; i < inputData.size(); ++i) {
50 inputPtrs[i]
51 = ttkUtils::GetPointer<T>(this->GetInputArrayToProcess(0, inputData[i]));
52 }
53 return this->execute(distanceMatrix, inputPtrs, nPoints);
54}
55
56int ttkLDistanceMatrix::RequestData(vtkInformation * /*request*/,
57 vtkInformationVector **inputVector,
58 vtkInformationVector *outputVector) {
59 ttk::Timer tm{};
60
61 const auto blocks = vtkMultiBlockDataSet::GetData(inputVector[0], 0);
62
63 if(blocks == nullptr) {
64 return 0;
65 }
66
67 const size_t nInputs{blocks->GetNumberOfBlocks()};
68
69 // Get input data
70 std::vector<vtkDataSet *> inputData(nInputs);
71
72 for(size_t i = 0; i < nInputs; ++i) {
73 inputData[i] = vtkDataSet::SafeDownCast(blocks->GetBlock(i));
74 }
75
76 // sanity check: data non null and same data size
77 {
78 std::set<vtkIdType> sizes{};
79 for(size_t i = 0; i < nInputs; ++i) {
80 if(inputData[i] == nullptr) {
81 this->printErr("One input block is not a vtkDataSet");
82 return 0;
83 }
84 sizes.emplace(inputData[i]->GetNumberOfPoints());
85 }
86 if(sizes.size() > 1) {
87 this->printErr("Input blocks do not have the same number of points");
88 return 0;
89 }
90 }
91
92 // Get output
93 auto DistTable = vtkTable::GetData(outputVector);
94
95 std::vector<std::vector<double>> distMatrix{};
96
97 const auto firstField = this->GetInputArrayToProcess(0, inputData[0]);
98 const auto dataType = firstField->GetDataType();
99 const size_t nPoints = firstField->GetNumberOfTuples();
100
101 switch(dataType) {
102 vtkTemplateMacro(this->dispatch<VTK_TT>(distMatrix, inputData, nPoints));
103 }
104
105 // zero-padd column name to keep Row Data columns ordered
106 const auto zeroPad
107 = [](std::string &colName, const size_t numberCols, const size_t colIdx) {
108 std::string const max{std::to_string(numberCols - 1)};
109 std::string const cur{std::to_string(colIdx)};
110 std::string const zer(max.size() - cur.size(), '0');
111 colName.append(zer).append(cur);
112 };
113
114 // copy distance matrix to output
115 for(size_t i = 0; i < nInputs; ++i) {
116 std::string name{"Dataset"};
117 zeroPad(name, distMatrix.size(), i);
118
119 vtkNew<vtkDoubleArray> col{};
120 col->SetNumberOfTuples(nInputs);
121 col->SetName(name.c_str());
122 for(size_t j = 0; j < nInputs; ++j) {
123 col->SetTuple1(j, distMatrix[i][j]);
124 }
125 DistTable->AddColumn(col);
126 }
127
128 // aggregate input field data
129 vtkNew<vtkFieldData> fd{};
130 fd->CopyStructure(inputData[0]->GetFieldData());
131 fd->SetNumberOfTuples(inputData.size());
132 for(size_t i = 0; i < inputData.size(); ++i) {
133 fd->SetTuple(i, 0, inputData[i]->GetFieldData());
134 }
135
136 // copy input field data to output row data
137 for(int i = 0; i < fd->GetNumberOfArrays(); ++i) {
138 DistTable->AddColumn(fd->GetAbstractArray(i));
139 }
140
141 this->printMsg("Complete (#datasets: " + std::to_string(nInputs)
142 + ", #points: " + std::to_string(nPoints) + ")",
143 1.0, tm.getElapsedTime(), this->threadNumber_);
144
145 return 1;
146}
Computes a distance matrix using LDistance between several input datasets with the same number of poi...
int FillInputPortInformation(int port, vtkInformation *info) override
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
int dispatch(std::vector< std::vector< double > > &distanceMatrix, const std::vector< vtkDataSet * > &inputData, const size_t nPoints)
int FillOutputPortInformation(int port, vtkInformation *info) override
int printErr(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition Debug.h:149
int execute(std::vector< TOut * > &output, const std::vector< const TIn * > &inputs, const size_t nPoints) const
vtkStandardNewMacro(ttkLDistanceMatrix)
printMsg(debug::output::BOLD+" | | | | | . \\ | | (__| | / __/| |_| / __/|__ _|"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)