TTK
Loading...
Searching...
No Matches
ttkDimensionReduction.cpp
Go to the documentation of this file.
2#include <ttkUtils.h>
3
4#include <vtkDoubleArray.h>
5#include <vtkInformation.h>
6#include <vtkIntArray.h>
7#include <vtkNew.h>
8#include <vtkObjectFactory.h>
9#include <vtkTable.h>
10
11#include <regex>
12
14
16 this->SetNumberOfInputPorts(1);
17 this->SetNumberOfOutputPorts(1);
18}
19
21 vtkInformation *info) {
22 if(port == 0) {
23 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkTable");
24 return 1;
25 }
26 return 0;
27}
28
30 vtkInformation *info) {
31 if(port == 0) {
32 info->Set(vtkDataObject::DATA_TYPE_NAME(), "vtkTable");
33 return 1;
34 }
35 return 0;
36}
37
39 vtkInformationVector **inputVector,
40 vtkInformationVector *outputVector) {
41
42 using ttk::SimplexId;
43
44 vtkTable *input = vtkTable::GetData(inputVector[0]);
45 vtkTable *output = vtkTable::GetData(outputVector);
46
47 // initialize method
49
50 if(SelectFieldsWithRegexp) {
51 // select all input columns whose name is matching the regexp
52 ScalarFields.clear();
53 const auto n = input->GetNumberOfColumns();
54 for(int i = 0; i < n; ++i) {
55 const auto &name = input->GetColumnName(i);
56 if(std::regex_match(name, std::regex(RegexpString))) {
57 ScalarFields.emplace_back(name);
58 }
59 }
60 }
61
62 const SimplexId numberOfRows = input->GetNumberOfRows();
63 const SimplexId numberOfColumns = ScalarFields.size();
64
65 if(numberOfRows <= 0 || numberOfColumns <= 0) {
66 this->printErr("Input matrix has invalid dimensions (rows: "
67 + std::to_string(numberOfRows)
68 + ", columns: " + std::to_string(numberOfColumns) + ")");
69 return 0;
70 }
71
72 std::vector<double> inputData;
73 std::vector<vtkAbstractArray *> arrays;
74 arrays.reserve(ScalarFields.size());
75 for(const auto &s : ScalarFields)
76 arrays.push_back(input->GetColumnByName(s.data()));
77 for(SimplexId i = 0; i < numberOfRows; ++i) {
78 for(auto arr : arrays)
79 inputData.push_back(arr->GetVariantValue(i).ToDouble());
80 }
81
82 outputData_.clear();
83
84 vtkNew<vtkIntArray> insertionTimeCol{};
85 int *insertionPtr = nullptr;
86 if(this->Method == DimensionReduction::METHOD::TOPOMAP) {
87 insertionTimeCol->SetNumberOfTuples(numberOfRows);
88 insertionTimeCol->SetName("InsertionTime");
89 insertionPtr = ttkUtils::GetPointer<int>(insertionTimeCol);
90 }
91
92 const int errorCode = this->execute(
93 this->outputData_, inputData, numberOfRows, numberOfColumns, insertionPtr);
94
95 if(!errorCode) {
96 if(KeepAllDataArrays)
97 output->ShallowCopy(input);
98
99 for(int i = 0; i < NumberOfComponents; ++i) {
100 std::string const s = "Component_" + std::to_string(i);
101 vtkNew<vtkDoubleArray> arr{};
102 ttkUtils::SetVoidArray(arr, outputData_[i].data(), numberOfRows, 1);
103 arr->SetName(s.data());
104 output->AddColumn(arr);
105 }
106 if(this->Method == DimensionReduction::METHOD::TOPOMAP) {
107 output->AddColumn(insertionTimeCol);
108 }
109 }
110
111 return 1;
112}
#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::DimensionReduction processing package.
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
int FillOutputPortInformation(int port, vtkInformation *info) override
int FillInputPortInformation(int port, vtkInformation *info) override
static void SetVoidArray(vtkDataArray *array, void *data, vtkIdType size, int save)
Definition ttkUtils.cpp:280
int printErr(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition Debug.h:149
void setInputMethod(METHOD method)
int execute(std::vector< std::vector< double > > &outputEmbedding, const std::vector< double > &inputMatrix, const int nRows, const int nColumns, int *insertionTimeForTopoMap=nullptr) const
int SimplexId
Identifier type for simplices of any dimension.
Definition DataTypes.h:22
vtkStandardNewMacro(ttkDimensionReduction)