TTK
Loading...
Searching...
No Matches
ttkPersistenceDiagram.cpp
Go to the documentation of this file.
1#include <vtkCellData.h>
2#include <vtkDataSet.h>
3#include <vtkDoubleArray.h>
4#include <vtkFloatArray.h>
5#include <vtkIdTypeArray.h>
6#include <vtkInformation.h>
7#include <vtkNew.h>
8#include <vtkPointData.h>
9#include <vtkSmartPointer.h>
10
13#include <ttkUtils.h>
14
16
18 SetNumberOfInputPorts(1);
19 SetNumberOfOutputPorts(1);
20}
21
23 vtkInformation *info) {
24 if(port == 0) {
25 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkDataSet");
26 return 1;
27 }
28 return 0;
29}
30
32 vtkInformation *info) {
33 if(port == 0) {
34 info->Set(vtkDataObject::DATA_TYPE_NAME(), "vtkUnstructuredGrid");
35 return 1;
36 }
37 return 0;
38}
39
40template <typename scalarType, typename triangulationType>
41int ttkPersistenceDiagram::dispatch(
42 vtkUnstructuredGrid *outputCTPersistenceDiagram,
43 vtkDataArray *const inputScalarsArray,
44 const scalarType *const inputScalars,
45 scalarType *outputScalars,
46 SimplexId *outputOffsets,
47 int *outputMonotonyOffsets,
48 const SimplexId *const inputOrder,
49 const triangulationType *triangulation) {
50
51 int status{};
52 ttk::DiagramType CTDiagram{};
53
55 std::cout << "Chosen approx" << std::endl;
56 double *range = inputScalarsArray->GetRange(0);
57 this->setDeltaApproximate(range[1] - range[0]);
58 this->setOutputScalars(outputScalars);
59 this->setOutputOffsets(outputOffsets);
60 this->setOutputMonotonyOffsets(outputMonotonyOffsets);
61 }
62
63 status = this->execute(CTDiagram, inputScalars, inputScalarsArray->GetMTime(),
64 inputOrder, triangulation);
65
66 // something wrong in baseCode
67 if(status != 0) {
68 this->printErr("PersistenceDiagram::execute() error code: "
69 + std::to_string(status));
70 return 0;
71 }
72
73 if(CTDiagram.empty()) {
74 this->printErr("Empty diagram!");
75 return 0;
76 }
77
78 vtkNew<vtkUnstructuredGrid> const vtu{};
79
80 // convert CTDiagram to vtkUnstructuredGrid
81 DiagramToVTU(vtu, CTDiagram, inputScalarsArray, *this,
82 triangulation->getDimensionality(), this->ShowInsideDomain);
83
84 outputCTPersistenceDiagram->ShallowCopy(vtu);
85
86 if(this->ClearDGCache && this->BackEnd == BACKEND::DISCRETE_MORSE_SANDWICH) {
87 this->printMsg("Clearing DiscreteGradient cache...");
89 }
90
91 return 1;
92}
93
95 vtkInformationVector **inputVector,
96 vtkInformationVector *outputVector) {
97
98 vtkDataSet *input = vtkDataSet::GetData(inputVector[0]);
99 vtkUnstructuredGrid *outputCTPersistenceDiagram
100 = vtkUnstructuredGrid::GetData(outputVector, 0);
101
103#ifndef TTK_ENABLE_KAMIKAZE
104 if(!triangulation) {
105 this->printErr("Wrong triangulation");
106 return 0;
107 }
108#endif
109
110 this->preconditionTriangulation(triangulation);
111
112 vtkDataArray *inputScalars = this->GetInputArrayToProcess(0, inputVector);
113#ifndef TTK_ENABLE_KAMIKAZE
114 if(!inputScalars) {
115 this->printErr("Wrong input scalars");
116 return 0;
117 }
118#endif
119
120 vtkDataArray *offsetField = this->GetOrderArray(
121 input, 0, triangulation, false, 1, ForceInputOffsetScalarField);
122
123#ifndef TTK_ENABLE_KAMIKAZE
124 if(!offsetField) {
125 this->printErr("Wrong input offsets");
126 return 0;
127 }
128 if(offsetField->GetDataType() != VTK_INT
129 and offsetField->GetDataType() != VTK_ID_TYPE) {
130 this->printErr("Input offset field type not supported");
131 return 0;
132 }
133#endif
134
135 vtkNew<ttkSimplexIdTypeArray> outputOffsets{};
136 outputOffsets->SetNumberOfComponents(1);
137 outputOffsets->SetNumberOfTuples(inputScalars->GetNumberOfTuples());
138 outputOffsets->SetName("outputOffsets");
139
140 vtkNew<vtkIntArray> outputMonotonyOffsets{};
141 outputMonotonyOffsets->SetNumberOfComponents(1);
142 outputMonotonyOffsets->SetNumberOfTuples(inputScalars->GetNumberOfTuples());
143 outputMonotonyOffsets->SetName("outputMonotonyffsets");
144 outputMonotonyOffsets->FillComponent(0, 0);
145
146 vtkSmartPointer<vtkDataArray> const outputScalars
147 = vtkSmartPointer<vtkDataArray>::Take(inputScalars->NewInstance());
148 outputScalars->SetNumberOfComponents(1);
149 outputScalars->SetNumberOfTuples(inputScalars->GetNumberOfTuples());
150 outputScalars->DeepCopy(inputScalars);
151 outputScalars->SetName("Cropped");
152
153 int status{};
155 inputScalars->GetDataType(), triangulation->getType(),
156 status = this->dispatch(
157 outputCTPersistenceDiagram, inputScalars,
158 static_cast<VTK_TT *>(ttkUtils::GetVoidPointer(inputScalars)),
159 static_cast<VTK_TT *>(ttkUtils::GetVoidPointer(outputScalars)),
160 static_cast<SimplexId *>(ttkUtils::GetVoidPointer(outputOffsets)),
161 static_cast<int *>(ttkUtils::GetVoidPointer(outputMonotonyOffsets)),
162 static_cast<SimplexId *>(ttkUtils::GetVoidPointer(offsetField)),
163 static_cast<TTK_TT *>(triangulation->getData())));
164
165 // shallow copy input Field Data
166 outputCTPersistenceDiagram->GetFieldData()->ShallowCopy(
167 input->GetFieldData());
168
169 return status;
170}
#define ttkNotUsed(x)
Mark function/method parameters that are not used in the function body at all.
Definition BaseClass.h:47
ttk::Triangulation * GetTriangulation(vtkDataSet *dataSet)
vtkDataArray * GetOrderArray(vtkDataSet *const inputData, const int scalarArrayIdx, ttk::Triangulation *triangulation, const bool getGlobalOrder=false, const int orderArrayIdx=0, const bool enforceOrderArrayIdx=false)
TTK VTK-filter for the computation of persistence diagrams.
int FillInputPortInformation(int port, vtkInformation *info) override
int FillOutputPortInformation(int port, vtkInformation *info) override
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
static void * GetVoidPointer(vtkDataArray *array, vtkIdType start=0)
Definition ttkUtils.cpp:226
int printErr(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition Debug.h:149
void setOutputMonotonyOffsets(void *data)
void preconditionTriangulation(AbstractTriangulation *triangulation)
void setDeltaApproximate(double data)
void setOutputScalars(void *data)
int execute(std::vector< PersistencePair > &CTDiagram, const scalarType *inputScalars, const size_t scalarsMTime, const SimplexId *inputOffsets, const triangulationType *triangulation, const std::vector< bool > *updateMask=nullptr)
void setOutputOffsets(void *data)
Triangulation is a class that provides time and memory efficient traversal methods on triangulations ...
AbstractTriangulation * getData()
Triangulation::Type getType() const
static void clearCache(const AbstractTriangulation &triangulation)
std::vector< PersistencePair > DiagramType
Persistence Diagram type as a vector of Persistence pairs.
#define ttkVtkTemplateMacro(dataType, triangulationType, call)
Definition ttkMacros.h:69
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...
vtkStandardNewMacro(ttkPersistenceDiagram)
printMsg(debug::output::BOLD+" | | | | | . \\ | | (__| | / __/| |_| / __/|__ _|"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)