1#include <vtkCellData.h>
2#include <vtkDataArray.h>
4#include <vtkDoubleArray.h>
5#include <vtkGenericCell.h>
6#include <vtkIdTypeArray.h>
7#include <vtkInformation.h>
8#include <vtkIntArray.h>
10#include <vtkPointData.h>
11#include <vtkPointSet.h>
12#include <vtkPolyData.h>
13#include <vtkSmartPointer.h>
14#include <vtkUnstructuredGrid.h>
25 this->SetNumberOfInputPorts(3);
26 this->SetNumberOfOutputPorts(1);
30 int port, vtkInformation *info) {
32 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(),
"vtkDataSet");
34 }
else if(port == 1) {
35 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(),
"vtkPolyData");
37 }
else if(port == 2) {
38 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(),
"vtkUnstructuredGrid");
45 int port, vtkInformation *info) {
53int ttkStableManifoldPersistence::AttachPersistence(vtkDataSet *output)
const {
60 auto ascendingManifoldArray
62 auto descendingManifoldArray
70 if((!ascendingManifoldArray) && (!descendingManifoldArray)
72 printErr(
"The input #0 is not a valid stable manifold.");
76 bool isSegmentation =
false;
78 if((ascendingManifoldArray) || (descendingManifoldArray))
79 isSegmentation =
true;
91 int const cellNumber = output->GetNumberOfCells();
93 persistenceArray->SetNumberOfTuples(cellNumber);
94 pairTypeArray->SetNumberOfTuples(cellNumber);
96#ifdef TTK_ENABLE_OPENMP
97#pragma omp parallel for num_threads(threadNumber_)
99 for(
int i = 0; i < cellNumber; i++) {
101 if((destinationArray) && (IsUnstable)) {
102 destinationArray->GetTuple(i, &cellId);
104 sourceArray->GetTuple(i, &cellId);
106 double const persistence = simplex2persistence_[(int)cellId];
107 double const pairType = simplex2pairType_[(int)cellId];
108 persistenceArray->SetTuple(i, &persistence);
109 pairTypeArray->SetTuple(i, &pairType);
112 output->GetCellData()->AddArray(persistenceArray);
113 output->GetCellData()->AddArray(pairTypeArray);
116 if((IsUnstable && descendingManifoldArray ==
nullptr)
117 || ascendingManifoldArray ==
nullptr) {
122 int const vertexNumber = output->GetNumberOfPoints();
124 persistenceArray->SetNumberOfTuples(vertexNumber);
125 pairTypeArray->SetNumberOfTuples(vertexNumber);
127#ifdef TTK_ENABLE_OPENMP
128#pragma omp parallel for num_threads(threadNumber_)
130 for(
int i = 0; i < vertexNumber; i++) {
132 double extremumId = -1;
134 descendingManifoldArray->GetTuple(i, &extremumId);
135 cellId = min2simplex_[(int)extremumId];
137 ascendingManifoldArray->GetTuple(i, &extremumId);
138 cellId = max2simplex_[(int)extremumId];
140 double const persistence = simplex2persistence_[cellId];
141 double const pairType = simplex2pairType_[cellId];
142 persistenceArray->SetTuple(i, &persistence);
143 pairTypeArray->SetTuple(i, &pairType);
146 output->GetPointData()->AddArray(persistenceArray);
147 output->GetPointData()->AddArray(pairTypeArray);
155int ttkStableManifoldPersistence::BuildSimplex2PersistenceMap(
156 vtkDataSet *stableManifold,
157 vtkPolyData *criticalPoints,
158 vtkUnstructuredGrid *persistenceDiagram) {
165 vtkDataArray *criticalPointVertexIdArray
167 vtkDataArray *criticalCellIdArray
169 vtkDataArray *criticalCellDimensionArray
170 = criticalPoints->GetPointData()->GetArray(
172 if((!criticalPointVertexIdArray) || (!criticalCellIdArray)
173 || (!criticalCellDimensionArray)) {
174 printErr(
"The input #1 is not a valid set of critical points.");
178 vtkDataArray *vertexIdArray
180 vtkDataArray *persistenceArray
182 vtkDataArray *persistencePairIdArray
183 = persistenceDiagram->GetCellData()->GetArray(
185 vtkDataArray *persistencePairTypeArray
188 if((!vertexIdArray) || (!persistenceArray) || (!persistencePairIdArray)
189 || (!persistencePairTypeArray)) {
190 printErr(
"The input #2 is not a valid persistence diagram.");
194 int const maximumVertexId = vertexIdArray->GetMaxNorm();
195 std::vector<double> vertex2persistence(maximumVertexId + 1, -1);
196 std::vector<int> vertex2pairType(maximumVertexId + 1, -1);
198 int const cellNumber = persistenceDiagram->GetNumberOfCells();
201 for(
int i = 0; i < cellNumber; i++) {
204 persistencePairIdArray->GetTuple(i, &pairId);
209 vtkNew<vtkGenericCell> c;
210 persistenceDiagram->GetCell(i, c);
211 int const pointId0 = c->GetPointId(0);
212 int const pointId1 = c->GetPointId(1);
215 persistenceArray->GetTuple(i, &persistence);
217 double pairType = -1;
218 persistencePairTypeArray->GetTuple(i, &pairType);
220 double vertexId0 = -1, vertexId1 = -1;
221 vertexIdArray->GetTuple(pointId0, &vertexId0);
222 vertexIdArray->GetTuple(pointId1, &vertexId1);
224 if(vertex2persistence[(
int)vertexId0] < persistence) {
225 vertex2persistence[(int)vertexId0] = persistence;
226 vertex2pairType[(int)vertexId0] = pairType;
228 if(vertex2persistence[(
int)vertexId1] <
persistence) {
229 vertex2persistence[(int)vertexId1] = persistence;
230 vertex2pairType[(int)vertexId1] = pairType;
235 int const maximumSimplexId = criticalCellIdArray->GetMaxNorm();
236 simplex2persistence_.resize(maximumSimplexId + 1, -1);
237 simplex2pairType_.resize(maximumSimplexId + 1, -1);
239 int const criticalPointNumber = criticalCellIdArray->GetNumberOfTuples();
240#ifdef TTK_ENABLE_OPENMP
241#pragma omp parallel for num_threads(threadNumber_)
243 for(
int i = 0; i < criticalPointNumber; i++) {
245 criticalCellIdArray->GetTuple(i, &cellId);
247 double vertexId = -1;
248 criticalPointVertexIdArray->GetTuple(i, &vertexId);
250 simplex2persistence_[(int)cellId] = vertex2persistence[(
int)vertexId];
251 simplex2pairType_[(int)cellId] = vertex2pairType[(
int)vertexId];
258 if(stableManifold->GetNumberOfCells()) {
259 int const cellType = stableManifold->GetCellType(0);
261 if((cellType == VTK_TETRA) || (cellType == VTK_VOXEL))
263 else if((cellType == VTK_TRIANGLE) || (cellType == VTK_PIXEL))
268 min2simplex_.clear();
269 max2simplex_.clear();
270 for(
int i = 0; i < criticalPointNumber; i++) {
272 double criticalSimplexDimension = -1;
273 criticalCellDimensionArray->GetTuple(i, &criticalSimplexDimension);
274 criticalCellIdArray->GetTuple(i, &cellId);
275 if(criticalSimplexDimension == 0) {
277 min2simplex_.push_back(cellId);
278 }
else if(criticalSimplexDimension == dimension) {
280 max2simplex_.push_back(cellId);
293 vtkInformationVector **inputVector,
294 vtkInformationVector *outputVector) {
298 const auto stableManifold = vtkDataSet::GetData(inputVector[0]);
299 const auto criticalPoints = vtkPolyData::GetData(inputVector[1]);
300 const auto persistenceDiagram = vtkUnstructuredGrid::GetData(inputVector[2]);
302 int ret = this->BuildSimplex2PersistenceMap(
303 stableManifold, criticalPoints, persistenceDiagram);
308 auto output = vtkDataSet::GetData(outputVector);
309 output->ShallowCopy(stableManifold);
311 ret = this->AttachPersistence(output);
#define ttkNotUsed(x)
Mark function/method parameters that are not used in the function body at all.
static vtkInformationIntegerKey * SAME_DATA_TYPE_AS_INPUT_PORT()
TTK VTK-filter for attaching to an input stable manifold (given by the Morse-Smale complex module) it...
int FillOutputPortInformation(int port, vtkInformation *info) override
int FillInputPortInformation(int port, vtkInformation *info) override
ttkStableManifoldPersistence()
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
void setDebugMsgPrefix(const std::string &prefix)
int printErr(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
const char MorseSmaleDestinationIdName[]
const char PersistenceName[]
const char MorseSmaleSourceIdName[]
const char MorseSmaleAscendingName[]
const char PersistencePairTypeName[]
const char MorseSmaleCellDimensionName[]
const char MorseSmaleDescendingName[]
const char MorseSmaleCellIdName[]
const char VertexScalarFieldName[]
default name for vertex scalar field
const char PersistencePairIdentifierName[]
vtkStandardNewMacro(ttkStableManifoldPersistence)
printMsg(debug::output::BOLD+" | | | | | . \\ | | (__| | / __/| |_| / __/|__ _|"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)