TTK
Loading...
Searching...
No Matches
ttkDimensionReductionMetrics.cpp
Go to the documentation of this file.
2
3#include <vtkDoubleArray.h>
4#include <vtkInformation.h>
5#include <vtkObjectFactory.h>
6#include <vtkStringArray.h>
7#include <vtkTable.h>
8
9#include <regex>
10
12
14 this->SetNumberOfInputPorts(2);
15 this->SetNumberOfOutputPorts(1);
16}
17
19 int port, vtkInformation *info) {
20 if(port == 0 || port == 1) {
21 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkTable");
22 return 1;
23 }
24 return 0;
25}
26
28 int port, vtkInformation *info) {
29 if(port == 0) {
30 info->Set(vtkDataObject::DATA_TYPE_NAME(), "vtkTable");
31 return 1;
32 }
33 return 0;
34}
35
37 vtkInformation *ttkNotUsed(request),
38 vtkInformationVector **inputVector,
39 vtkInformationVector *outputVector) {
40
41 vtkTable *input = vtkTable::GetData(inputVector[0]);
42 vtkTable *representation = vtkTable::GetData(inputVector[1]);
43 if(!input || !representation)
44 return 0;
45
46 vtkTable *output = vtkTable::GetData(outputVector, 0);
47
48 /*** input points ***/
49 if(SelectInputFieldsWithRegexp) {
50 // select all input columns whose name is matching the regexp
51 InputScalarFields.clear();
52 const auto n = input->GetNumberOfColumns();
53 for(int i = 0; i < n; ++i) {
54 const auto &name = input->GetColumnName(i);
55 if(std::regex_match(name, std::regex(InputRegexpString))) {
56 InputScalarFields.emplace_back(name);
57 }
58 }
59 }
60 if(input->GetNumberOfRows() <= 0 || InputScalarFields.size() <= 0) {
61 this->printErr("Input matrix has invalid dimensions (rows: "
62 + std::to_string(input->GetNumberOfRows()) + ", columns: "
63 + std::to_string(InputScalarFields.size()) + ")");
64 return 0;
65 }
66 std::vector<vtkAbstractArray *> inputArrays;
67 inputArrays.reserve(InputScalarFields.size());
68 for(const auto &s : InputScalarFields)
69 inputArrays.push_back(input->GetColumnByName(s.data()));
70
71 const int numberOfPoints = input->GetNumberOfRows();
72 const int dimensionHigh = InputScalarFields.size();
73 std::vector<std::vector<double>> inputPoints(numberOfPoints);
74 for(int i = 0; i < numberOfPoints; ++i) {
75 for(int j = 0; j < dimensionHigh; ++j)
76 inputPoints[i].push_back(inputArrays[j]->GetVariantValue(i).ToDouble());
77 }
78
79 /*** representation points ***/
80 if(SelectRepresentationFieldsWithRegexp) {
81 // select all input columns whose name is matching the regexp
82 RepresentationScalarFields.clear();
83 const auto n = representation->GetNumberOfColumns();
84 for(int i = 0; i < n; ++i) {
85 const auto &name = representation->GetColumnName(i);
86 if(std::regex_match(name, std::regex(RepresentationRegexpString))) {
87 RepresentationScalarFields.emplace_back(name);
88 }
89 }
90 }
91 if(representation->GetNumberOfRows() != numberOfPoints
92 || RepresentationScalarFields.size() <= 0) {
93 this->printErr("Representation matrix has invalid dimensions (rows: "
94 + std::to_string(representation->GetNumberOfRows())
95 + ", columns: "
96 + std::to_string(RepresentationScalarFields.size()) + ")");
97 return 0;
98 }
99 std::vector<vtkAbstractArray *> representationArrays;
100 representationArrays.reserve(RepresentationScalarFields.size());
101 for(const auto &s : RepresentationScalarFields)
102 representationArrays.push_back(representation->GetColumnByName(s.data()));
103
104 const int dimensionLow = RepresentationScalarFields.size();
105 std::vector<std::vector<double>> representationPoints(numberOfPoints);
106 for(int i = 0; i < numberOfPoints; ++i) {
107 for(int j = 0; j < dimensionLow; ++j)
108 representationPoints[i].push_back(
109 representationArrays[j]->GetVariantValue(i).ToDouble());
110 }
111
112 execute(inputPoints, representationPoints);
113
114 auto metricNames = vtkSmartPointer<vtkStringArray>::New();
115 metricNames->SetName("Metric");
116 metricNames->InsertNextValue("0-Wasserstein");
117 metricNames->InsertNextValue("1-Wasserstein");
118 metricNames->InsertNextValue("Triplet accuracy");
119 metricNames->InsertNextValue("Linear correlation");
120 metricNames->InsertNextValue("Trustworthiness");
121 metricNames->InsertNextValue("Continuity");
122 metricNames->InsertNextValue("LCMC");
123 metricNames->InsertNextValue("MRRE (input)");
124 metricNames->InsertNextValue("MRRE (latent)");
125 metricNames->InsertNextValue("RMSE");
126 output->AddColumn(metricNames);
127
128 auto metricValues = vtkSmartPointer<vtkDoubleArray>::New();
129 metricValues->SetName("Value");
130 metricValues->InsertNextValue(w0_);
131 metricValues->InsertNextValue(w1_);
132 metricValues->InsertNextValue(ta_);
133 metricValues->InsertNextValue(lc_);
134 metricValues->InsertNextValue(trust_);
135 metricValues->InsertNextValue(cont_);
136 metricValues->InsertNextValue(lcmc_);
137 metricValues->InsertNextValue(mrreh_);
138 metricValues->InsertNextValue(mrrel_);
139 metricValues->InsertNextValue(rmse_);
140 output->AddColumn(metricValues);
141
142 auto metricRanges = vtkSmartPointer<vtkStringArray>::New();
143 metricRanges->SetName("Range");
144 metricRanges->InsertNextValue("[0,inf)");
145 metricRanges->InsertNextValue("[0,inf)");
146 metricRanges->InsertNextValue("[0,1]");
147 metricRanges->InsertNextValue("[-1,1]");
148 metricRanges->InsertNextValue("[0,1]");
149 metricRanges->InsertNextValue("[0,1]");
150 metricRanges->InsertNextValue("[0,1]");
151 metricRanges->InsertNextValue("[0,1]");
152 metricRanges->InsertNextValue("[0,1]");
153 metricRanges->InsertNextValue("[0,inf)");
154 output->AddColumn(metricRanges);
155
156 // return success
157 return 1;
158}
#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::DimensionReductionMetrics module.
int FillInputPortInformation(int port, vtkInformation *info) override
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
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
void execute(std::vector< std::vector< double > > const &input, std::vector< std::vector< double > > const &latent)
Main entry point.
vtkStandardNewMacro(ttkDimensionReductionMetrics)