TTK
Loading...
Searching...
No Matches
ttkMeshGraph.cpp
Go to the documentation of this file.
1#include <ttkMeshGraph.h>
2
3#include <vtkAbstractArray.h>
4#include <vtkCellArray.h>
5#include <vtkCellData.h>
6#include <vtkDataSetTriangleFilter.h>
7
8#include <vtkFloatArray.h>
9#include <vtkIdTypeArray.h>
10
11#include <vtkInformation.h>
12#include <vtkInformationVector.h>
13#include <vtkPointData.h>
14#include <vtkUnstructuredGrid.h>
15
16#include <vtkArrayDispatch.h>
17#include <vtkDataArrayAccessor.h>
18
19#include <ttkMacros.h>
20#include <ttkUtils.h>
21
23
25 this->SetNumberOfInputPorts(1);
26 this->SetNumberOfOutputPorts(1);
27}
28
30
31int ttkMeshGraph::FillInputPortInformation(int port, vtkInformation *info) {
32 if(port == 0) {
33 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkUnstructuredGrid");
34 return 1;
35 }
36 return 0;
37}
38
39int ttkMeshGraph::FillOutputPortInformation(int port, vtkInformation *info) {
40 if(port == 0) {
41 info->Set(vtkDataObject::DATA_TYPE_NAME(), "vtkUnstructuredGrid");
42 return 1;
43 }
44 return 0;
45}
46
47int ttkMeshGraph::RequestData(vtkInformation *ttkNotUsed(request),
48 vtkInformationVector **inputVector,
49 vtkInformationVector *outputVector) {
50 ttk::Timer t;
51
52 // ---------------------------------------------------------------------------
53 // Get Input
54 // ---------------------------------------------------------------------------
55 auto input = vtkUnstructuredGrid::GetData(inputVector[0]);
56 if(input == nullptr) {
57 return -1;
58 }
59
60 size_t const nInputPoints = input->GetNumberOfPoints();
61 size_t const nInputCells = input->GetNumberOfCells();
62
63 auto inputPointSizes = this->GetInputArrayToProcess(0, inputVector);
64 if(!inputPointSizes) {
65 this->printErr("Unable to retrieve point size array.");
66 return 0;
67 }
68
69 // ---------------------------------------------------------------------------
70 // Init Output
71 // ---------------------------------------------------------------------------
72
73 // Output Points
74 auto inputPoints = input->GetPoints();
75 auto outputPoints = vtkSmartPointer<vtkPoints>::New();
76
77 outputPoints->SetDataType(inputPoints->GetDataType());
78 auto nOutputPoints = this->computeNumberOfOutputPoints(
79 nInputPoints, nInputCells, this->GetUseQuadraticCells(),
80 this->GetSubdivisions());
81 outputPoints->SetNumberOfPoints(nOutputPoints);
82
83 // Output Topology
84 auto inputConnectivityArray = input->GetCells()->GetConnectivityArray();
85
86 auto nOutputCells = this->computeNumberOfOutputCells(
87 nInputCells, this->GetUseQuadraticCells());
88 auto outputTopologySize = this->computeOutputConnectivityArraySize(
89 nInputCells, this->GetUseQuadraticCells(), this->GetSubdivisions());
90
91 auto outputConnectivityArray = vtkSmartPointer<vtkDataArray>::Take(
92 inputConnectivityArray->NewInstance());
93 outputConnectivityArray->SetNumberOfValues(outputTopologySize);
94
95 auto outputOffsetArray = vtkSmartPointer<vtkDataArray>::Take(
96 inputConnectivityArray->NewInstance());
97 outputOffsetArray->SetNumberOfValues(nOutputCells + 1);
98
99 // ---------------------------------------------------------------------------
100 // Compute cells with base code
101 // ---------------------------------------------------------------------------
102
103 int status = 0;
104 if(this->GetUseQuadraticCells()) {
106 inputPoints->GetDataType(), inputPointSizes->GetDataType(),
107 outputConnectivityArray->GetDataType(),
108 (status = this->execute<T2, T0, T1>(
109 // Output
110 ttkUtils::GetPointer<T0>(outputPoints->GetData()),
111 ttkUtils::GetPointer<T2>(outputConnectivityArray),
112 ttkUtils::GetPointer<T2>(outputOffsetArray),
113
114 // Input
115 ttkUtils::GetPointer<T0>(inputPoints->GetData()),
116 ttkUtils::GetPointer<T2>(input->GetCells()->GetConnectivityArray()),
117 nInputPoints, nInputCells, ttkUtils::GetPointer<T1>(inputPointSizes),
118 this->GetSizeScale(), this->GetSizeAxis())));
119 } else {
121 inputPoints->GetDataType(), inputPointSizes->GetDataType(),
122 outputConnectivityArray->GetDataType(),
123 (status = this->execute2<T2, T0, T1>(
124 // Output
125 ttkUtils::GetPointer<T0>(outputPoints->GetData()),
126 ttkUtils::GetPointer<T2>(outputConnectivityArray),
127 ttkUtils::GetPointer<T2>(outputOffsetArray),
128
129 // Input
130 ttkUtils::GetPointer<T0>(inputPoints->GetData()),
131 ttkUtils::GetPointer<T2>(input->GetCells()->GetConnectivityArray()),
132 nInputPoints, nInputCells, this->GetSubdivisions(),
133 ttkUtils::GetPointer<T1>(inputPointSizes), this->GetSizeScale(),
134 this->GetSizeAxis())));
135 }
136 if(!status)
137 return 0;
138
139 // ---------------------------------------------------------------------------
140 // Generate meshed graph as vtkUnstructuredGrid
141 // ---------------------------------------------------------------------------
142
143 // Create new vtkUnstructuredGrid for meshed graph
144 auto meshedGraph = vtkUnstructuredGrid::GetData(outputVector);
145
146 meshedGraph->SetPoints(outputPoints);
147
148 auto outputCellArray = vtkSmartPointer<vtkCellArray>::New();
149 outputCellArray->SetData(outputOffsetArray, outputConnectivityArray);
150 meshedGraph->SetCells(
151 this->GetUseQuadraticCells() ? VTK_QUADRATIC_QUAD : VTK_POLYGON,
152 outputCellArray);
153
154 // Copy input point data to output point data
155 {
156 auto iPointData = input->GetPointData();
157 auto oPointData = meshedGraph->GetPointData();
158
159 for(int i = 0; i < iPointData->GetNumberOfArrays(); i++) {
160 auto iArray = iPointData->GetArray(i);
161 if(iArray->GetNumberOfComponents() > 1)
162 continue;
163
165 vtkDataArray::CreateDataArray(iArray->GetDataType()));
166 oArray->SetName(iArray->GetName());
167 oArray->SetNumberOfTuples(nOutputPoints);
168 oArray->SetNumberOfComponents(1);
169 oPointData->AddArray(oArray);
170
172 iArray->GetDataType(), inputConnectivityArray->GetDataType(),
173 (status = this->mapInputPointDataToOutputPointData<T0, T1>(
174 ttkUtils::GetPointer<T0>(oArray),
175
176 nInputPoints, nInputCells,
177 ttkUtils::GetPointer<T1>(inputConnectivityArray),
178 ttkUtils::GetPointer<T0>(iArray), this->GetUseQuadraticCells(),
179 this->GetSubdivisions())));
180
181 if(!status)
182 return 0;
183 }
184 }
185
186 // Copy input cell data to output cell data
187 {
188 auto iCellData = input->GetCellData();
189 auto oCellData = meshedGraph->GetCellData();
190
191 for(int i = 0; i < iCellData->GetNumberOfArrays(); i++) {
192 auto iArray = iCellData->GetArray(i);
193 if(iArray->GetNumberOfComponents() > 1)
194 continue;
195
197 vtkDataArray::CreateDataArray(iArray->GetDataType()));
198 oArray->SetName(iArray->GetName());
199 oArray->SetNumberOfTuples(nOutputCells);
200 oArray->SetNumberOfComponents(1);
201 oCellData->AddArray(oArray);
202
204 iArray->GetDataType(),
205 (status = this->mapInputCellDataToOutputCellData<T0>(
206 ttkUtils::GetPointer<T0>(oArray), nInputCells,
207 ttkUtils::GetPointer<T0>(iArray), this->GetUseQuadraticCells())));
208 if(!status)
209 return 0;
210 }
211 }
212
213 // -------------------------------------------------------------------------
214 // Print status
215 // -------------------------------------------------------------------------
217 this->printMsg("Complete", 1, t.getElapsedTime());
219
220 return 1;
221}
#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 generates a mesh for a graph.
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
virtual bool GetUseQuadraticCells()
int FillInputPortInformation(int port, vtkInformation *info) override
int FillOutputPortInformation(int port, vtkInformation *info) override
virtual int GetSubdivisions()
~ttkMeshGraph() override
int printErr(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition Debug.h:149
size_t computeOutputConnectivityArraySize(const size_t &nInputCells, const bool &useQuadraticCells, const size_t &nSubdivisions=0) const
Definition MeshGraph.h:85
size_t computeNumberOfOutputPoints(const size_t &nInputPoints, const size_t &nInputCells, const bool &useQuadraticCells, const size_t &nSubdivisions=0) const
Definition MeshGraph.h:54
size_t computeNumberOfOutputCells(const size_t &nInputCells, const bool &useQuadraticCells) const
Definition MeshGraph.h:71
double getElapsedTime()
Definition Timer.h:15
#define ttkTypeMacroRRI(group0, group1, group2, call)
Definition ttkMacros.h:425
#define ttkTypeMacroA(group, call)
Definition ttkMacros.h:234
#define ttkTypeMacroAI(group0, group1, call)
Definition ttkMacros.h:341
vtkStandardNewMacro(ttkMeshGraph)
printMsg(debug::output::BOLD+" | | | | | . \\ | | (__| | / __/| |_| / __/|__ _|"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)