TTK
Loading...
Searching...
No Matches
ttkDelaunayRipsPersistenceDiagram.cpp
Go to the documentation of this file.
3
4#include <vtkCellData.h>
5#include <vtkDoubleArray.h>
6#include <vtkInformation.h>
7#include <vtkInformationVector.h>
8#include <vtkPointData.h>
9#include <vtkTable.h>
10
11#include <regex>
12
14
15using namespace ttk::rpd;
16
18 this->SetNumberOfInputPorts(1);
19 this->SetNumberOfOutputPorts(1);
20}
21
23 int port, vtkInformation *info) {
24 if(port == 0) {
25 info->Append(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkTable");
26 info->Append(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkPointSet");
27 return 1;
28 }
29 return 0;
30}
31
33 int port, vtkInformation *info) {
34 if(port == 0) {
35 info->Set(vtkDataObject::DATA_TYPE_NAME(), "vtkUnstructuredGrid");
36 return 1;
37 }
38 return 0;
39}
40
42 vtkInformation *ttkNotUsed(request),
43 vtkInformationVector **inputVector,
44 vtkInformationVector *outputVector) {
45
46 ttk::Timer tm{};
47
48 vtkInformation *info = inputVector[0]->GetInformationObject(0);
49 vtkDataObject *input = info->Get(vtkDataObject::DATA_OBJECT());
50 vtkUnstructuredGrid *outputPersistenceDiagram
51 = vtkUnstructuredGrid::GetData(outputVector);
52
53 if(!input)
54 return 0;
55
56 PointCloud points;
57 int numberOfPoints = 0;
58 int dimension = 0;
59
60 if(vtkTable *table = vtkTable::SafeDownCast(input)) {
61 if(SelectFieldsWithRegexp) {
62 // select all input columns whose name is matching the regexp
63 ScalarFields.clear();
64 const auto n = table->GetNumberOfColumns();
65 for(int i = 0; i < n; ++i) {
66 const auto &name = table->GetColumnName(i);
67 if(std::regex_match(name, std::regex(RegexpString))) {
68 ScalarFields.emplace_back(name);
69 }
70 }
71 }
72
73 if(table->GetNumberOfRows() <= 0 || ScalarFields.size() <= 1) {
74 this->printErr("Input matrix has invalid dimensions (rows: "
75 + std::to_string(table->GetNumberOfRows()) + ", columns: "
76 + std::to_string(ScalarFields.size()) + ")");
77 return 0;
78 }
79
80 std::vector<vtkAbstractArray *> arrays;
81 arrays.reserve(ScalarFields.size());
82 for(const auto &s : ScalarFields)
83 arrays.push_back(table->GetColumnByName(s.data()));
84
85 numberOfPoints = table->GetNumberOfRows();
86 dimension = ScalarFields.size();
87
88 points.resize(numberOfPoints);
89 for(int i = 0; i < numberOfPoints; ++i) {
90 for(int j = 0; j < dimension; ++j)
91 points[i].push_back(arrays[j]->GetVariantValue(i).ToDouble());
92 }
93 }
94
95 else if(vtkPointSet *pointset = vtkPointSet::SafeDownCast(input)) {
96 numberOfPoints = pointset->GetNumberOfPoints();
97 dimension = 3;
98 points.resize(numberOfPoints, std::vector<double>(3));
99 for(int i = 0; i < numberOfPoints; ++i)
100 pointset->GetPoint(i, points[i].data());
101 }
102
103 this->printMsg("Computing Delaunay-Rips persistence diagram", 1.0,
104 tm.getElapsedTime(), getThreadNumber());
105 this->printMsg("#dimensions: " + std::to_string(dimension)
106 + ", #points: " + std::to_string(numberOfPoints),
107 0.0, tm.getElapsedTime(), getThreadNumber());
108
110 if(this->execute(points, diagram) != 0)
111 return 0;
112
113 DiagramToVTU(outputPersistenceDiagram, diagram, inf);
114
115 this->printMsg("Complete", 1.0, tm.getElapsedTime(), getThreadNumber());
116
117 // shallow copy input Field Data
118 outputPersistenceDiagram->GetFieldData()->ShallowCopy(input->GetFieldData());
119
120 // return success
121 return 1;
122}
#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::DelaunayRipsPersistenceDiagram module.
int FillOutputPortInformation(int port, vtkInformation *info) override
int FillInputPortInformation(int port, vtkInformation *info) override
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
int getThreadNumber() const
Definition BaseClass.h:76
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 rpd::PointCloud &points, rpd::MultidimensionalDiagram &ph) const
Main entry point (without generators)
constexpr value_t inf
std::vector< std::vector< value_t > > PointCloud
std::vector< Diagram > MultidimensionalDiagram
vtkStandardNewMacro(ttkDelaunayRipsPersistenceDiagram)
int DiagramToVTU(vtkUnstructuredGrid *vtu, const ttk::DiagramType &diagram, vtkDataArray *const inputScalars, const ttk::Debug &dbg, const int dim, const bool embedInDomain)
Converts a Persistence Diagram in the ttk::DiagramType format to the VTK Unstructured Grid format (as...
printMsg(debug::output::BOLD+" | | | | | . \\ | | (__| | / __/| |_| / __/| (_) |"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)