TTK
Loading...
Searching...
No Matches
ttkDistanceMatrixDistortion.cpp
Go to the documentation of this file.
2
3#include <vtkInformation.h>
4
5#include <vtkDataArray.h>
6#include <vtkDataSet.h>
7#include <vtkDoubleArray.h>
8#include <vtkObjectFactory.h>
9#include <vtkPointData.h>
10#include <vtkSmartPointer.h>
11#include <vtkTable.h>
12#include <vtkVariantArray.h>
13
14#include <regex>
15
16#include <ttkMacros.h>
17#include <ttkUtils.h>
18
19// A VTK macro that enables the instantiation of this class via ::New()
20// You do not have to modify this
22
24 this->SetNumberOfInputPorts(2);
25 this->SetNumberOfOutputPorts(1);
26}
27
29 int port, vtkInformation *info) {
30 if(port == 0) {
31 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkTable");
32 return 1;
33 } else if(port == 1) {
34 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkTable");
35 return 1;
36 }
37 return 0;
38}
39
41 int port, vtkInformation *info) {
42 if(port == 0) {
43 info->Set(vtkDataObject::DATA_TYPE_NAME(), "vtkTable");
44 return 1;
45 }
46 return 0;
47}
48
49inline void fillWithInputColumns(vtkTable *input,
50 const std::string &regExpStr,
51 std::vector<std::string> &vectToFill) {
52 // Select all input columns whose name is matching the regexp.
53 vectToFill.clear();
54 const size_t n = input->GetNumberOfColumns();
55 for(size_t i = 0; i < n; ++i) {
56 const auto &name = input->GetColumnName(i);
57 if(std::regex_match(name, std::regex(regExpStr))) {
58 vectToFill.emplace_back(name);
59 }
60 }
61}
62
64 vtkInformation *ttkNotUsed(request),
65 vtkInformationVector **inputVector,
66 vtkInformationVector *outputVector) {
67
68 // Get input objects from input vector
69 vtkTable *inputHigh = vtkTable::GetData(inputVector[0]);
70 vtkTable *inputLow = vtkTable::GetData(inputVector[1]);
71
72 vtkTable *output = vtkTable::GetData(outputVector);
73
74 if(!inputLow || !inputHigh || !output)
75 return 0;
76
77 // Selecting the columns names if a regexp was given.
78 if(SelectFieldsWithRegexpHigh)
79 fillWithInputColumns(inputHigh, RegexpStringHigh, ScalarFieldsHigh);
80 if(SelectFieldsWithRegexpLow)
81 fillWithInputColumns(inputLow, RegexpStringLow, ScalarFieldsLow);
82
83 const size_t nRowsHigh = inputHigh->GetNumberOfRows();
84 const size_t nColsHigh = ScalarFieldsHigh.size();
85 const size_t nRowsLow = inputLow->GetNumberOfRows();
86 const size_t nColsLow = ScalarFieldsLow.size();
87 if(nRowsHigh == 0 || nRowsHigh != nColsHigh) {
88 this->printErr("High input matrix is not a valid square matrix (rows: "
89 + std::to_string(nRowsHigh)
90 + ", columns: " + std::to_string(nColsHigh) + ")");
91 return 0;
92 }
93
94 if(nRowsLow == 0 || nRowsLow != nColsLow) {
95 this->printErr("Low input matrix is not a valid square matrix (rows: "
96 + std::to_string(nRowsLow)
97 + ", columns: " + std::to_string(nColsLow) + ")");
98 return 0;
99 }
100 if(nRowsHigh != nRowsLow) {
101 this->printErr(
102 "High and low input matrices must have same size (rows(high): "
103 + std::to_string(nRowsHigh) + ", rows(low): " + std::to_string(nRowsLow)
104 + ")");
105 return 0;
106 }
107
108 int const n = nRowsHigh;
109 std::vector<double *> vectMatHigh(n),
110 vectMatLow(n); // No 2D vectors to avoid copy of data from the VTK layer.
111
112 // Getting the actual input columns.
113 std::vector<vtkDataArray *> arraysHigh{}, arraysLow{};
114 for(const auto &s : ScalarFieldsHigh)
115 arraysHigh.push_back(
116 vtkDoubleArray::SafeDownCast(inputHigh->GetColumnByName(s.data())));
117 for(const auto &s : ScalarFieldsLow)
118 arraysLow.push_back(
119 vtkDoubleArray::SafeDownCast(inputLow->GetColumnByName(s.data())));
120
121#ifdef TTK_ENABLE_OPENMP
122#pragma omp parallel for num_threads(this->threadNumber_)
123#endif // TTK_ENABLE_OPENMP
124 for(int i = 0; i < n; i++) {
125 vectMatHigh[i] = ttkUtils::GetPointer<double>(arraysHigh[i]);
126 vectMatLow[i] = ttkUtils::GetPointer<double>(arraysLow[i]);
127 }
128
129 double distortionValue = 0;
130 vtkNew<vtkDoubleArray> tmpCol{}, distortionValArray{};
131 tmpCol->SetNumberOfTuples(n);
132 this->printMsg("Starting computation of sim distortion value...");
133 this->execute(vectMatHigh, vectMatLow, distortionValue,
134 ttkUtils::GetPointer<double>(tmpCol));
135
136 tmpCol->SetName("SimValue");
137 // No deep copy, makes output->RowData points to the data of tmpCol.
138 output->AddColumn(tmpCol);
139 distortionValArray->SetName("DistortionValue");
140 distortionValArray->SetNumberOfTuples(1);
141 distortionValArray->SetTuple1(0, distortionValue);
142 output->GetFieldData()->AddArray(distortionValArray);
143
144 return 1;
145}
#define ttkNotUsed(x)
Mark function/method parameters that are not used in the function body at all.
Definition BaseClass.h:47
TTK VTK-filter that wraps the ttk::DistanceMatrixDistortion module.
int FillInputPortInformation(int port, vtkInformation *info) override
int FillOutputPortInformation(int port, vtkInformation *info) override
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) 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(const std::vector< double * > &highDistMatrix, const std::vector< double * > &lowDistMatrix, double &distortionValue, double *distortionVerticesValues) const
vtkStandardNewMacro(ttkDistanceMatrixDistortion)
void fillWithInputColumns(vtkTable *input, const std::string &regExpStr, std::vector< std::string > &vectToFill)
printMsg(debug::output::BOLD+" | | | | | . \\ | | (__| | / __/| |_| / __/|__ _|"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)