TTK
Loading...
Searching...
No Matches
ttkPointMerger.cpp
Go to the documentation of this file.
1#include <vtkCellData.h>
2#include <vtkDoubleArray.h>
3#include <vtkGenericCell.h>
4#include <vtkInformation.h>
5#include <vtkInformationVector.h>
6#include <vtkNew.h>
7#include <vtkObjectFactory.h>
8#include <vtkPointData.h>
9#include <vtkUnstructuredGrid.h>
10
11#include <Geometry.h>
12#include <Triangulation.h>
13#include <ttkPointMerger.h>
14
15#include <array>
16
18
20 this->setDebugMsgPrefix("PointMerger");
21 this->SetNumberOfInputPorts(1);
22 this->SetNumberOfOutputPorts(1);
23}
24
25int ttkPointMerger::FillInputPortInformation(int port, vtkInformation *info) {
26 if(port == 0) {
27 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkUnstructuredGrid");
28 return 1;
29 }
30 return 0;
31}
32
33int ttkPointMerger::FillOutputPortInformation(int port, vtkInformation *info) {
34 if(port == 0) {
36 return 1;
37 }
38 return 0;
39}
40
41int ttkPointMerger::RequestData(vtkInformation *ttkNotUsed(request),
42 vtkInformationVector **inputVector,
43 vtkInformationVector *outputVector) {
44
45 using ttk::SimplexId;
46 ttk::Timer t;
47
48 auto input = vtkUnstructuredGrid::GetData(inputVector[0]);
49 auto output = vtkUnstructuredGrid::GetData(outputVector);
50
51 auto triangulation = ttkAlgorithm::GetTriangulation(input);
52
53 if(triangulation == nullptr) {
54 return -1;
55 }
56
57 if(BoundaryOnly) {
58 triangulation->preconditionBoundaryVertices();
59 }
60
61 SimplexId const vertexNumber = input->GetNumberOfPoints();
62 std::vector<SimplexId> candidateVertices;
63
64 if(BoundaryOnly) {
65 for(SimplexId i = 0; i < vertexNumber; i++) {
66 if(triangulation->isVertexOnBoundary(i)) {
67 candidateVertices.push_back(i);
68 }
69 }
70 } else {
71 candidateVertices.resize(vertexNumber);
72 for(SimplexId i = 0; i < vertexNumber; i++)
73 candidateVertices[i] = i;
74 }
75
76 std::vector<std::vector<SimplexId>> closePoints(vertexNumber);
77
78 this->printMsg("Computing pointwise distances with "
79 + std::to_string(candidateVertices.size()) + " candidates");
80
81#ifdef TTK_ENABLE_OPENMP
82#pragma omp parallel for num_threads(threadNumber_)
83#endif
84 for(SimplexId i = 0; i < (SimplexId)candidateVertices.size(); i++) {
85 double distance = -1;
86 std::array<double, 3> p0{};
87
88 SimplexId const vertexId0 = candidateVertices[i];
89
90 input->GetPoint(vertexId0, p0.data());
91 for(SimplexId j = 0; j < (SimplexId)candidateVertices.size(); j++) {
92 if(i != j) {
93 std::array<double, 3> p1{};
94 SimplexId const vertexId1 = candidateVertices[j];
95 input->GetPoint(vertexId1, p1.data());
96 distance = ttk::Geometry::distance(p0.data(), p1.data());
97 if(distance < DistanceThreshold) {
98 closePoints[vertexId0].push_back(vertexId1);
99 }
100 }
101 }
102 }
103
104 this->printMsg("Merging points...");
105
106 std::vector<double> minMergeDistance(vertexNumber, -1);
107 std::vector<double> maxMergeDistance(vertexNumber, -1);
108 std::vector<SimplexId> mergeCount(vertexNumber, 0);
109 std::vector<SimplexId> mergeMap(vertexNumber, -1);
110
111 for(SimplexId i = 0; i < vertexNumber; i++) {
112 if(closePoints[i].size()) {
113
114 SimplexId targetVertexId = i;
115 do {
116 SimplexId nextTargetVertexId = targetVertexId;
117 for(SimplexId j = 0; j < (SimplexId)closePoints[targetVertexId].size();
118 j++) {
119 if(closePoints[targetVertexId][j] < nextTargetVertexId)
120 nextTargetVertexId = closePoints[targetVertexId][j];
121 }
122 targetVertexId = nextTargetVertexId;
123 mergeMap[i] = targetVertexId;
124 } while(mergeMap[targetVertexId] != targetVertexId);
125 mergeCount[targetVertexId]++;
126
127 if(i != targetVertexId) {
128 std::array<double, 3> p0{}, p1{};
129 input->GetPoint(i, p0.data());
130 input->GetPoint(targetVertexId, p1.data());
131 double const distance = ttk::Geometry::distance(p0.data(), p1.data());
132 if((minMergeDistance[targetVertexId] == -1)
133 || (distance < minMergeDistance[targetVertexId])) {
134 minMergeDistance[targetVertexId] = distance;
135 }
136 if((maxMergeDistance[targetVertexId] == -1)
137 || (distance > maxMergeDistance[targetVertexId])) {
138 maxMergeDistance[targetVertexId] = distance;
139 }
140 }
141 } else {
142 mergeMap[i] = i;
143 }
144 }
145
146 SimplexId vertexIdGen = 0;
147 std::vector<SimplexId> old2new(vertexNumber, 0);
148 std::vector<SimplexId> new2old;
149 for(SimplexId i = 0; i < vertexNumber; i++) {
150 if(mergeMap[i] == i) {
151 old2new[i] = vertexIdGen;
152 vertexIdGen++;
153 new2old.push_back(i);
154 } else {
155 old2new[i] = old2new[mergeMap[i]];
156 }
157 }
158
159 // now create the output
160 vtkNew<vtkPoints> pointSet{};
161 pointSet->SetNumberOfPoints(vertexIdGen);
162
163 vtkNew<vtkIdTypeArray> mergeCountArray{};
164 mergeCountArray->SetNumberOfTuples(vertexIdGen);
165 mergeCountArray->SetName("VertexMergeCount");
166
167 vtkNew<vtkDoubleArray> minDistanceArray{};
168 minDistanceArray->SetNumberOfTuples(vertexIdGen);
169 minDistanceArray->SetName("MinMergeDistance");
170
171 vtkNew<vtkDoubleArray> maxDistanceArray{};
172 maxDistanceArray->SetNumberOfTuples(vertexIdGen);
173 maxDistanceArray->SetName("MaxMergeDistance");
174
175 std::vector<vtkNew<vtkDoubleArray>> pointData{};
176 pointData.resize(input->GetPointData()->GetNumberOfArrays());
177 for(size_t i = 0; i < pointData.size(); i++) {
178 pointData[i]->SetName(input->GetPointData()->GetArray(i)->GetName());
179 pointData[i]->SetNumberOfComponents(
180 input->GetPointData()->GetArray(i)->GetNumberOfComponents());
181 pointData[i]->SetNumberOfTuples(vertexIdGen);
182 }
183
184#ifdef TTK_ENABLE_OPENMP
185#pragma omp parallel for num_threads(threadNumber_)
186#endif
187 for(SimplexId i = 0; i < vertexIdGen; i++) {
188 std::array<double, 3> p{};
189 input->GetPoint(new2old[i], p.data());
190 pointSet->SetPoint(i, p.data());
191
192 mergeCountArray->SetTuple1(i, mergeCount[new2old[i]]);
193 minDistanceArray->SetTuple1(i, minMergeDistance[new2old[i]]);
194 maxDistanceArray->SetTuple1(i, maxMergeDistance[new2old[i]]);
195
196 for(size_t j = 0; j < pointData.size(); j++) {
197 std::vector<double> data(pointData[j]->GetNumberOfComponents());
198 input->GetPointData()->GetArray(j)->GetTuple(new2old[i], data.data());
199 pointData[j]->SetTuple(i, data.data());
200 }
201 }
202
203 output->SetPoints(pointSet);
204 for(size_t i = 0; i < pointData.size(); i++) {
205 output->GetPointData()->AddArray(pointData[i]);
206 }
207 output->GetPointData()->AddArray(mergeCountArray);
208 output->GetPointData()->AddArray(minDistanceArray);
209 output->GetPointData()->AddArray(maxDistanceArray);
210
211 // now do the cells
212 output->Allocate();
213 std::vector<vtkNew<vtkDoubleArray>> cellData{};
214 cellData.resize(input->GetCellData()->GetNumberOfArrays());
215 for(size_t i = 0; i < cellData.size(); i++) {
216 cellData[i]->SetName(input->GetCellData()->GetArray(i)->GetName());
217 cellData[i]->SetNumberOfComponents(
218 input->GetCellData()->GetArray(i)->GetNumberOfComponents());
219 }
220
221 for(SimplexId i = 0; i < input->GetNumberOfCells(); i++) {
222 vtkNew<vtkGenericCell> c{};
223 input->GetCell(i, c);
224
225 std::vector<SimplexId> newVertexIds;
226 for(int j = 0; j < c->GetNumberOfPoints(); j++) {
227 SimplexId const vertexId = old2new[c->GetPointId(j)];
228 bool isIn = false;
229 for(SimplexId k = 0; k < (SimplexId)newVertexIds.size(); k++) {
230 if(newVertexIds[k] == vertexId) {
231 isIn = true;
232 break;
233 }
234 }
235 if(!isIn) {
236 newVertexIds.push_back(vertexId);
237 }
238 }
239
240 vtkNew<vtkIdList> idList{};
241 idList->SetNumberOfIds(newVertexIds.size());
242 for(size_t j = 0; j < newVertexIds.size(); j++) {
243 idList->SetId(j, newVertexIds[j]);
244 }
245
246 vtkIdType cellId = -1;
247
248 if((c->GetCellDimension() == 1) && (newVertexIds.size() == 2)) {
249 cellId = output->InsertNextCell(VTK_LINE, idList);
250 }
251 if(c->GetCellDimension() == 2) {
252 if(newVertexIds.size() > 2) {
253 cellId = output->InsertNextCell(VTK_POLYGON, idList);
254 }
255 }
256 if(c->GetCellDimension() == 3) {
257 if(newVertexIds.size() == 4) {
258 cellId = output->InsertNextCell(VTK_TETRA, idList);
259 }
260 if(newVertexIds.size() == 8) {
261 cellId = output->InsertNextCell(VTK_HEXAHEDRON, idList);
262 } else {
263 this->printWrn("Ill-defined cell type for cell #" + std::to_string(i)
264 + " (" + std::to_string(newVertexIds.size())
265 + " vertices)");
266 }
267 }
268
269 if(cellId != -1) {
270 // insert the cell data
271 for(size_t j = 0; j < cellData.size(); j++) {
272 std::vector<double> data(cellData[j]->GetNumberOfComponents());
273 input->GetCellData()->GetArray(j)->GetTuple(i, data.data());
274 cellData[j]->InsertNextTuple(data.data());
275 }
276 }
277 }
278 for(size_t i = 0; i < cellData.size(); i++)
279 output->GetCellData()->AddArray(cellData[i]);
280
281 this->printMsg(
282 "Merge performed", 1.0, t.getElapsedTime(), this->threadNumber_);
283
284 return 1;
285}
#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()
ttk::Triangulation * GetTriangulation(vtkDataSet *dataSet)
TTK VTK-filter for point merging.
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
int FillInputPortInformation(int port, vtkInformation *info) override
int FillOutputPortInformation(int port, vtkInformation *info) override
int printWrn(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition Debug.h:159
void setDebugMsgPrefix(const std::string &prefix)
Definition Debug.h:364
double getElapsedTime()
Definition Timer.h:15
T distance(const T *p0, const T *p1, const int &dimension=3)
Definition Geometry.cpp:362
int SimplexId
Identifier type for simplices of any dimension.
Definition DataTypes.h:22
vtkStandardNewMacro(ttkPointMerger)
printMsg(debug::output::BOLD+" | | | | | . \\ | | (__| | / __/| |_| / __/|__ _|"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)