TTK
Loading...
Searching...
No Matches
ttkPlanarGraphLayout.cpp
Go to the documentation of this file.
2
3#include <ttkMacros.h>
4
5#include <vtkAbstractArray.h>
6#include <vtkCellArray.h>
7#include <vtkFloatArray.h>
8#include <vtkInformation.h>
9#include <vtkInformationVector.h>
10#include <vtkPointData.h>
11#include <vtkPolyData.h>
12#include <vtkSmartPointer.h>
13#include <vtkUnstructuredGrid.h>
14
15#include <FTMTreePPUtils.h>
16#include <ttkMergeTreeUtils.h>
18
20
22 this->SetNumberOfInputPorts(1);
23 this->SetNumberOfOutputPorts(1);
24}
26
28 vtkInformation *info) {
29 if(port == 0) {
30 info->Set(vtkAlgorithm::INPUT_IS_REPEATABLE(), 1);
31 info->Append(
32 vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkUnstructuredGrid");
33 info->Append(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkPolyData");
34 } else
35 return 0;
36 return 1;
37}
38
40 vtkInformation *info) {
41 if(port == 0)
43 else
44 return 0;
45 return 1;
46}
47
49 vtkInformation *ttkNotUsed(request),
50 vtkInformationVector **inputVector,
51 vtkInformationVector *outputVector) {
52 // Get input and output objects
53
54 auto input = vtkDataSet::GetData(inputVector[0]);
55 auto output = vtkDataSet::GetData(outputVector);
56
57 if(!input || !output)
58 return !this->printErr("Unable to retrieve input/output data objects.");
59
60 // Copy input to output
61 output->ShallowCopy(input);
62
63 size_t const nPoints = output->GetNumberOfPoints();
64 size_t const nEdges = output->GetNumberOfCells();
65
66 // Get input arrays
67 auto sequenceArray = this->GetInputArrayToProcess(0, inputVector);
68 if(this->GetUseSequences() && !sequenceArray)
69 return !this->printErr("Unable to retrieve sequence array.");
70
71 auto sizeArray = this->GetInputArrayToProcess(1, inputVector);
72 if(this->GetUseSizes() && !sizeArray)
73 return !this->printErr("Unable to retrieve size array.");
74
75 auto branchArray = this->GetInputArrayToProcess(2, inputVector);
76 if(this->GetUseBranches() && !branchArray)
77 return !this->printErr("Unable to retrieve branch array.");
78
79 auto levelArray = this->GetInputArrayToProcess(3, inputVector);
80 if(this->GetUseLevels() && !levelArray)
81 return !this->printErr("Unable to retrieve level array.");
82
83 // Initialize output array
84 auto outputArray = vtkSmartPointer<vtkFloatArray>::New();
85 outputArray->SetName(this->GetOutputArrayName().data());
86 outputArray->SetNumberOfComponents(2); // (x,y) position
87 outputArray->SetNumberOfValues(nPoints * 2);
88
89 vtkDataArray *cells{nullptr};
90 if(auto outputAsUG = vtkUnstructuredGrid::SafeDownCast(output))
91 cells = outputAsUG->GetCells()->GetConnectivityArray();
92 else if(auto outputAsPD = vtkPolyData::SafeDownCast(output))
93 cells = outputAsPD->GetLines()->GetConnectivityArray();
94
95 if(!cells)
96 return !this->printErr("Unable to retrieve connectivity array.");
97
98 int status = 1;
100 this->GetUseSequences() ? sequenceArray->GetDataType() : VTK_INT,
101 this->GetUseBranches() ? branchArray->GetDataType() : VTK_INT,
102 cells->GetDataType(),
103 (status = this->computeLayout<T0, T1, T2>(
104 // Output
105 ttkUtils::GetPointer<float>(outputArray),
106 // Input
107 ttkUtils::GetPointer<T2>(cells), nPoints, nEdges,
108 this->GetUseSequences() ? ttkUtils::GetPointer<T0>(sequenceArray)
109 : nullptr,
110 this->GetUseSizes() ? ttkUtils::GetPointer<float>(sizeArray) : nullptr,
111 this->GetUseBranches() ? ttkUtils::GetPointer<T1>(branchArray) : nullptr,
112 this->GetUseLevels() ? ttkUtils::GetPointer<T1>(levelArray) : nullptr)));
113
114 if(status != 1)
115 return 0;
116
117 // Add output field to output
118 output->GetPointData()->AddArray(outputArray);
119
120 return 1;
121}
122
123template <class dataType>
125 vtkUnstructuredGrid *treeNodes,
126 vtkUnstructuredGrid *treeArcs,
127 vtkUnstructuredGrid *output) {
128 auto mergeTree = ttk::ftm::makeTree<dataType>(treeNodes, treeArcs);
129 ttk::ftm::FTMTree_MT *tree = &(mergeTree.tree);
130
131 ttk::ftm::computePersistencePairs<dataType>(tree);
132
133 std::vector<std::vector<int>> treeNodeCorrMesh(1);
134 treeNodeCorrMesh[0] = std::vector<int>(tree->getNumberOfNodes());
135 for(unsigned int j = 0; j < tree->getNumberOfNodes(); ++j)
136 treeNodeCorrMesh[0][j] = j;
137
139 visuMaker.setPlanarLayout(true);
140 visuMaker.setOutputSegmentation(false);
141 visuMaker.setBranchDecompositionPlanarLayout(BranchDecompositionPlanarLayout);
142 visuMaker.setBranchSpacing(BranchSpacing);
143 visuMaker.setImportantPairs(ImportantPairs);
144 visuMaker.setMaximumImportantPairs(MaximumImportantPairs);
145 visuMaker.setMinimumImportantPairs(MinimumImportantPairs);
146 visuMaker.setImportantPairsSpacing(ImportantPairsSpacing);
147 visuMaker.setNonImportantPairsSpacing(NonImportantPairsSpacing);
148 visuMaker.setNonImportantPairsProximity(NonImportantPairsProximity);
149 visuMaker.setExcludeImportantPairsHigher(ExcludeImportantPairsHigher);
150 visuMaker.setExcludeImportantPairsLower(ExcludeImportantPairsLower);
151 visuMaker.setVtkOutputNode(output);
152 visuMaker.setVtkOutputArc(output);
153 visuMaker.setTreesNodes(treeNodes);
154 visuMaker.setTreesNodeCorrMesh(treeNodeCorrMesh);
155 visuMaker.setDebugLevel(this->debugLevel_);
156 visuMaker.copyPointData(treeNodes);
157 visuMaker.makeTreesOutput<dataType>(tree);
158
159 return 1;
160}
161
163 vtkInformation *ttkNotUsed(request),
164 vtkInformationVector **inputVector,
165 vtkInformationVector *outputVector) {
166
167 vtkUnstructuredGrid *treeNodes
168 = vtkUnstructuredGrid::GetData(inputVector[0], 0);
169 vtkUnstructuredGrid *treeArcs
170 = vtkUnstructuredGrid::GetData(inputVector[0], 1);
171 auto output = vtkUnstructuredGrid::GetData(outputVector);
172
173 return mergeTreePlanarLayoutCallTemplate<float>(treeNodes, treeArcs, output);
174}
175
176int ttkPlanarGraphLayout::RequestData(vtkInformation *request,
177 vtkInformationVector **inputVector,
178 vtkInformationVector *outputVector) {
179 if(not InputIsAMergeTree)
180 return planarGraphLayoutCall(request, inputVector, outputVector);
181 else
182 return mergeTreePlanarLayoutCall(request, inputVector, outputVector);
183}
#define ttkNotUsed(x)
Mark function/method parameters that are not used in the function body at all.
Definition BaseClass.h:47
static vtkInformationIntegerKey * SAME_DATA_TYPE_AS_INPUT_PORT()
void copyPointData(vtkUnstructuredGrid *treeNodes, std::vector< int > &nodeCorrT)
void setTreesNodeCorrMesh(std::vector< std::vector< int > > &nodeCorrMesh)
void makeTreesOutput(FTMTree_MT *tree1)
void setVtkOutputArc(vtkUnstructuredGrid *vtkArc)
void setVtkOutputNode(vtkUnstructuredGrid *vtkNode)
void setTreesNodes(std::vector< vtkUnstructuredGrid * > &nodes)
TTK VTK-filter that computes a planar graph layout.
int mergeTreePlanarLayoutCallTemplate(vtkUnstructuredGrid *treeNodes, vtkUnstructuredGrid *treeArcs, vtkUnstructuredGrid *output)
int mergeTreePlanarLayoutCall(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
virtual bool GetUseSizes()
virtual bool GetUseLevels()
int planarGraphLayoutCall(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
~ttkPlanarGraphLayout() override
virtual bool GetUseSequences()
int FillInputPortInformation(int port, vtkInformation *info) override
int FillOutputPortInformation(int port, vtkInformation *info) override
virtual std::string GetOutputArrayName()
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
virtual bool GetUseBranches()
int debugLevel_
Definition Debug.h:379
virtual int setDebugLevel(const int &debugLevel)
Definition Debug.cpp:147
int printErr(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition Debug.h:149
void setExcludeImportantPairsLower(std::string &d)
void setExcludeImportantPairsHigher(std::string &d)
idNode getNumberOfNodes() const
Definition FTMTree_MT.h:389
#define ttkTypeMacroAII(group0, group1, group2, call)
Definition ttkMacros.h:406
vtkStandardNewMacro(ttkPlanarGraphLayout)