TTK
Loading...
Searching...
No Matches
ttkMatrixToHeatMap.cpp
Go to the documentation of this file.
2#include <vtkCellData.h>
3#include <vtkDataArray.h>
4#include <vtkDataObject.h>
5#include <vtkDoubleArray.h>
6#include <vtkFiltersCoreModule.h>
7#include <vtkIdTypeArray.h>
8#include <vtkInformation.h>
9#include <vtkInformationVector.h>
10#include <vtkObjectFactory.h>
11#include <vtkPolyData.h>
12#include <vtkTable.h>
13
14#include <array>
15#include <regex>
16
18
20 this->setDebugMsgPrefix("MatrixToHeatMap");
21 this->SetNumberOfInputPorts(1);
22 this->SetNumberOfOutputPorts(1);
23}
24
26 vtkInformation *info) {
27 if(port == 0) {
28 info->Set(vtkDataObject::DATA_TYPE_NAME(), "vtkTable");
29 return 1;
30 }
31 return 0;
32}
33
35 vtkInformation *info) {
36 if(port == 0) {
37 info->Set(vtkDataObject::DATA_TYPE_NAME(), "vtkPolyData");
38 return 1;
39 }
40 return 0;
41}
42
43int ttkMatrixToHeatMap::RequestData(vtkInformation * /*request*/,
44 vtkInformationVector **inputVector,
45 vtkInformationVector *outputVector) {
46 ttk::Timer tm{};
47
48 const auto input = vtkTable::GetData(inputVector[0]);
49 const auto output = vtkPolyData::GetData(outputVector);
50
51 if(SelectFieldsWithRegexp) {
52 // select all input columns whose name is matching the regexp
53 ScalarFields.clear();
54 const auto n = input->GetNumberOfColumns();
55 for(int i = 0; i < n; ++i) {
56 const auto &name = input->GetColumnName(i);
57 if(std::regex_match(name, std::regex(RegexpString))) {
58 ScalarFields.emplace_back(name);
59 }
60 }
61 } else {
62 // Remove manually selected arrays that are not present in the input table
63 ScalarFields.erase(
64 std::remove_if(ScalarFields.begin(), ScalarFields.end(),
65 [&input](const std::string &name) {
66 return not input->GetColumnByName(name.data());
67 }),
68 ScalarFields.end());
69 }
70
71 const auto nInputs = ScalarFields.size();
72 const auto nRows = static_cast<size_t>(input->GetNumberOfRows());
73 if(nInputs != nRows) {
74 this->printWrn("Distance matrix is not square");
75 }
76
77 // grid points
78 vtkNew<vtkPoints> points{};
79 points->SetNumberOfPoints((nInputs + 1) * (nRows + 1));
80 // grid cells connectivity arrays
81 vtkNew<vtkIdTypeArray> offsets{}, connectivity{};
82 offsets->SetNumberOfComponents(1);
83 offsets->SetNumberOfTuples(nInputs * nRows + 1);
84 connectivity->SetNumberOfComponents(1);
85 connectivity->SetNumberOfTuples(4 * nInputs * nRows);
86 // distance and proximity arrays
87 vtkNew<vtkDoubleArray> dist{}, prox{};
88 dist->SetNumberOfComponents(1);
89 dist->SetName("Distance");
90 dist->SetNumberOfTuples(nInputs * nRows);
91 prox->SetNumberOfComponents(1);
92 prox->SetName("Proximity");
93 prox->SetNumberOfTuples(nInputs * nRows);
94
95#ifdef TTK_ENABLE_OPENMP
96#pragma omp parallel for num_threads(this->threadNumber_)
97#endif // TTK_ENABLE_OPENMP
98 for(size_t i = 0; i < nInputs + 1; ++i) {
99 for(size_t j = 0; j < nRows + 1; ++j) {
100 points->SetPoint(i * (nRows + 1) + j, i, j, 0.0);
101 }
102 }
103
104#ifdef TTK_ENABLE_OPENMP
105#pragma omp parallel for num_threads(this->threadNumber_)
106#endif // TTK_ENABLE_OPENMP
107 for(size_t i = 0; i < nInputs; ++i) {
108 for(size_t j = 0; j < nRows; ++j) {
109 // build the cell connectivity and offsets array
110 const auto nptline{static_cast<vtkIdType>(nRows + 1)};
111 const auto curr{static_cast<vtkIdType>(i * nptline + j)};
112 const auto o = i * nRows + j;
113 connectivity->SetTuple1(4 * o + 0, curr);
114 connectivity->SetTuple1(4 * o + 1, curr + nptline);
115 connectivity->SetTuple1(4 * o + 2, curr + nptline + 1);
116 connectivity->SetTuple1(4 * o + 3, curr + 1);
117 offsets->SetTuple1(o, 4 * o);
118
119 // copy distance matrix to heat map cell data
120 const auto val = input->GetColumnByName(ScalarFields[i].data())
121 ->GetVariantValue(j)
122 .ToDouble();
123 dist->SetTuple1(i * nRows + j, val);
124 prox->SetTuple1(i * nRows + j, std::exp(-val));
125 }
126 }
127 offsets->SetTuple1(nInputs * nRows, connectivity->GetNumberOfTuples());
128
129 // gather arrays to make the PolyData
130 vtkNew<vtkCellArray> cells{};
131 cells->SetData(offsets, connectivity);
132 output->SetPoints(points);
133 output->SetPolys(cells);
134 output->GetCellData()->AddArray(dist);
135 output->GetCellData()->AddArray(prox);
136
137 this->printMsg("Complete (#inputs: " + std::to_string(nInputs) + ")", 1.0,
138 tm.getElapsedTime(), this->threadNumber_);
139
140 return 1;
141}
Generates a Heat Map from a Distance Matrix stored into a vtkTable.
int FillOutputPortInformation(int port, vtkInformation *info) override
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
int FillInputPortInformation(int port, vtkInformation *info) override
int printWrn(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition Debug.h:159
void setDebugMsgPrefix(const std::string &prefix)
Definition Debug.h:364
vtkStandardNewMacro(ttkMatrixToHeatMap)
printMsg(debug::output::BOLD+" | | | | | . \\ | | (__| | / __/| |_| / __/|__ _|"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)