1#include <vtkCellData.h>
3#include <vtkDoubleArray.h>
4#include <vtkFloatArray.h>
5#include <vtkIdTypeArray.h>
6#include <vtkInformation.h>
8#include <vtkPointData.h>
9#include <vtkSmartPointer.h>
18 SetNumberOfInputPorts(1);
19 SetNumberOfOutputPorts(1);
23 vtkInformation *info) {
25 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(),
"vtkDataSet");
32 vtkInformation *info) {
34 info->Set(vtkDataObject::DATA_TYPE_NAME(),
"vtkUnstructuredGrid");
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) {
55 std::cout <<
"Chosen approx" << std::endl;
56 double *range = inputScalarsArray->GetRange(0);
63 status = this->
execute(CTDiagram, inputScalars, inputScalarsArray->GetMTime(),
64 inputOrder, triangulation);
68 this->
printErr(
"PersistenceDiagram::execute() error code: "
69 + std::to_string(status));
73 if(CTDiagram.empty()) {
78 vtkNew<vtkUnstructuredGrid>
const vtu{};
82 triangulation->getDimensionality(), this->ShowInsideDomain);
84 outputCTPersistenceDiagram->ShallowCopy(vtu);
87 this->
printMsg(
"Clearing DiscreteGradient cache...");
95 vtkInformationVector **inputVector,
96 vtkInformationVector *outputVector) {
98 vtkDataSet *input = vtkDataSet::GetData(inputVector[0]);
99 vtkUnstructuredGrid *outputCTPersistenceDiagram
100 = vtkUnstructuredGrid::GetData(outputVector, 0);
103#ifndef TTK_ENABLE_KAMIKAZE
105 this->
printErr(
"Wrong triangulation");
112 vtkDataArray *inputScalars = this->GetInputArrayToProcess(0, inputVector);
113#ifndef TTK_ENABLE_KAMIKAZE
115 this->
printErr(
"Wrong input scalars");
121 input, 0, triangulation,
false, 1, ForceInputOffsetScalarField);
123#ifndef TTK_ENABLE_KAMIKAZE
125 this->
printErr(
"Wrong input offsets");
128 if(offsetField->GetDataType() != VTK_INT
129 and offsetField->GetDataType() != VTK_ID_TYPE) {
130 this->
printErr(
"Input offset field type not supported");
135 vtkNew<ttkSimplexIdTypeArray> outputOffsets{};
136 outputOffsets->SetNumberOfComponents(1);
137 outputOffsets->SetNumberOfTuples(inputScalars->GetNumberOfTuples());
138 outputOffsets->SetName(
"outputOffsets");
140 vtkNew<vtkIntArray> outputMonotonyOffsets{};
141 outputMonotonyOffsets->SetNumberOfComponents(1);
142 outputMonotonyOffsets->SetNumberOfTuples(inputScalars->GetNumberOfTuples());
143 outputMonotonyOffsets->SetName(
"outputMonotonyffsets");
144 outputMonotonyOffsets->FillComponent(0, 0);
148 outputScalars->SetNumberOfComponents(1);
149 outputScalars->SetNumberOfTuples(inputScalars->GetNumberOfTuples());
150 outputScalars->DeepCopy(inputScalars);
151 outputScalars->SetName(
"Cropped");
155 inputScalars->GetDataType(), triangulation->
getType(),
156 status = this->dispatch(
157 outputCTPersistenceDiagram, inputScalars,
163 static_cast<TTK_TT *
>(triangulation->
getData())));
166 outputCTPersistenceDiagram->GetFieldData()->ShallowCopy(
167 input->GetFieldData());
#define ttkNotUsed(x)
Mark function/method parameters that are not used in the function body at all.
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)
int printErr(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
void setOutputMonotonyOffsets(void *data)
void preconditionTriangulation(AbstractTriangulation *triangulation)
void setDeltaApproximate(double data)
void setOutputScalars(void *data)
@ DISCRETE_MORSE_SANDWICH
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)
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)