TTK
Loading...
Searching...
No Matches
ttkFiberSurface.cpp
Go to the documentation of this file.
1#include <vtkCellArray.h>
2#include <vtkCellData.h>
3#include <vtkDoubleArray.h>
4#include <vtkInformation.h>
5#include <vtkNew.h>
6#include <vtkObjectFactory.h>
7#include <vtkPointData.h>
8#include <vtkPolyData.h>
9#include <vtkUnstructuredGrid.h>
10
11#include <ttkFiberSurface.h>
12#include <ttkMacros.h>
13#include <ttkUtils.h>
14
16
18 this->SetNumberOfInputPorts(2);
19 this->SetNumberOfOutputPorts(1);
20}
21
22int ttkFiberSurface::FillInputPortInformation(int port, vtkInformation *info) {
23 if(port == 0) {
24 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkDataSet");
25 return 1;
26 } else if(port == 1) {
27 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkUnstructuredGrid");
28 return 1;
29 }
30 return 0;
31}
32
33int ttkFiberSurface::FillOutputPortInformation(int port, vtkInformation *info) {
34 if(port == 0) {
35 info->Set(vtkDataObject::DATA_TYPE_NAME(), "vtkPolyData");
36 return 1;
37 }
38 return 0;
39}
40
41template <typename VTK_T1, typename VTK_T2>
43
44#ifdef TTK_ENABLE_FIBER_SURFACE_WITH_RANGE_OCTREE
45 if(RangeOctree) {
46 ttkTemplateMacro(triangulation->getType(),
47 (this->buildOctree<VTK_T1, VTK_T2>(
48 static_cast<TTK_TT *>(triangulation->getData()))));
49 }
50#endif // TTK_ENABLE_FIBER_SURFACE_WITH_RANGE_OCTREE
51
52 ttkTemplateMacro(triangulation->getType(),
53 (this->computeSurface<VTK_T1, VTK_T2>(
54 static_cast<TTK_TT *>(triangulation->getData()))));
55 return 0;
56}
57
58int ttkFiberSurface::RequestData(vtkInformation *ttkNotUsed(request),
59 vtkInformationVector **inputVector,
60 vtkInformationVector *outputVector) {
61
62 using ttk::SimplexId;
63 ttk::Timer const t;
64
65 const auto input = vtkDataSet::GetData(inputVector[0]);
66 const auto polygon = vtkUnstructuredGrid::GetData(inputVector[1]);
67 auto output = vtkPolyData::GetData(outputVector);
68
69 const auto dataUfield = this->GetInputArrayToProcess(0, input);
70 const auto dataVfield = this->GetInputArrayToProcess(1, input);
71 const auto polygonUfield = this->GetInputArrayToProcess(2, polygon);
72 const auto polygonVfield = this->GetInputArrayToProcess(3, polygon);
73
74 if(dataUfield == nullptr || dataVfield == nullptr || polygonUfield == nullptr
75 || polygonVfield == nullptr) {
76 this->printErr("Could not find data array");
77 return -1;
78 }
79
80 if(!(input->GetDataObjectType() == VTK_UNSTRUCTURED_GRID
81 || input->GetDataObjectType() == VTK_IMAGE_DATA)) {
82 this->printErr("Unsupported VTK data structure");
83 return -5;
84 }
85
86 auto triangulation = ttkAlgorithm::GetTriangulation(input);
87 if(triangulation == nullptr) {
88 return -1;
89 }
90
91 this->preconditionTriangulation(triangulation);
92 outputVertexList_.clear();
93 this->setGlobalVertexList(&outputVertexList_);
94 this->setInputField(
95 ttkUtils::GetVoidPointer(dataUfield), ttkUtils::GetVoidPointer(dataVfield));
96 this->setPolygonEdgeNumber(polygon->GetNumberOfCells());
97 threadedTriangleList_.resize(polygon->GetNumberOfCells());
98 threadedVertexList_.resize(polygon->GetNumberOfCells());
99 this->setPolygon(&inputPolygon_);
100
101 this->setPointMerging(PointMerge);
102 this->setPointMergingThreshold(PointMergeDistanceThreshold);
103
104#ifdef TTK_ENABLE_FIBER_SURFACE_WITH_RANGE_OCTREE
105 if((!RangeOctree) || (dataUfield->GetMTime() > GetMTime())
106 || (dataVfield->GetMTime() > GetMTime())) {
107
108 this->printMsg("Resetting octree...");
109
110 this->flushOctree();
111 Modified();
112 }
113#endif
114
115 inputPolygon_.clear();
116
117 SimplexId const cellNumber = polygon->GetNumberOfCells();
118 vtkCellArray *connectivity = polygon->GetCells();
119
120 if(connectivity->GetData()->GetNumberOfTuples() < 3 * cellNumber) {
121 this->printErr("Error: ill-defined range polygon.");
122 return 0;
123 }
124
125#if !defined(_WIN32) || defined(_WIN32) && defined(VTK_USE_64BIT_IDS)
126 const long long int *cellArray = connectivity->GetData()->GetPointer(0);
127#else
128 int *pt = connectivity->GetPointer();
129 long long extra_pt = *pt;
130 const long long int *cellArray = &extra_pt;
131#endif
132
133 SimplexId vertexId0, vertexId1;
134 std::pair<std::pair<double, double>, std::pair<double, double>> rangeEdge;
135
136 for(SimplexId i = 0; i < cellNumber; i++) {
137
138 vertexId0 = cellArray[3 * i + 1];
139 vertexId1 = cellArray[3 * i + 2];
140
141 rangeEdge.first.first = polygonUfield->GetTuple1(vertexId0);
142 rangeEdge.first.second = polygonVfield->GetTuple1(vertexId0);
143
144 rangeEdge.second.first = polygonUfield->GetTuple1(vertexId1);
145 rangeEdge.second.second = polygonVfield->GetTuple1(vertexId1);
146
147 inputPolygon_.push_back(rangeEdge);
148 }
149
150 for(size_t i = 0; i < threadedTriangleList_.size(); i++) {
151 threadedTriangleList_[i].clear();
152 this->setTriangleList(i, &(threadedTriangleList_[i]));
153 threadedVertexList_[i].clear();
154 this->setVertexList(i, &(threadedVertexList_[i]));
155 }
156
157#ifndef TTK_ENABLE_DOUBLE_TEMPLATING
158 if(dataUfield->GetDataType() != dataVfield->GetDataType()) {
159 this->printErr(
160 "Scalar fields should have same input type. Use TTKPointDataConverter or "
161 "TTKArrayEditor to convert array types.");
162 return 0;
163 }
164 switch(dataUfield->GetDataType()) {
165 vtkTemplateMacro((dispatch<VTK_TT, VTK_TT>(triangulation)));
166 }
167#else
168 switch(vtkTemplate2PackMacro(
169 dataUfield->GetDataType(), dataVfield->GetDataType())) {
170 vtkTemplate2Macro((dispatch<VTK_T1, VTK_T2>(triangulation)));
171 }
172#endif // TTK_ENABLE_DOUBLE_TEMPLATING
173
174 // prepare the VTK output
175 // NOTE: right now, there is a copy of the output data. this is no good.
176 // to fix.
177
178 size_t triangleNumber = 0;
179
180 for(size_t i = 0; i < threadedTriangleList_.size(); i++) {
181 triangleNumber += threadedTriangleList_[i].size();
182 }
183
184 vtkNew<vtkPoints> outputVertexList{};
185 vtkNew<vtkDoubleArray> outputU{};
186 vtkNew<vtkDoubleArray> outputV{};
187 vtkNew<vtkDoubleArray> outputParameterization{};
188 vtkNew<vtkCellArray> outputTriangleList{};
189 vtkNew<ttkSimplexIdTypeArray> outputEdgeIds{};
190 vtkNew<ttkSimplexIdTypeArray> outputTetIds{};
191 vtkNew<ttkSimplexIdTypeArray> outputCaseIds{};
192
193 if(RangeCoordinates) {
194 outputU->SetName(dataUfield->GetName());
195 outputU->SetNumberOfTuples(outputVertexList_.size());
196
197 outputV->SetName(dataVfield->GetName());
198 outputV->SetNumberOfTuples(outputVertexList_.size());
199 }
200
201 if(EdgeParameterization) {
202 outputParameterization->SetName("EdgeParameterization");
203 outputParameterization->SetNumberOfTuples(outputVertexList_.size());
204 }
205
206 outputVertexList->SetNumberOfPoints(outputVertexList_.size());
207 output->SetPoints(outputVertexList);
208
209#ifdef TTK_ENABLE_OPENMP
210#pragma omp parallel for num_threads(threadNumber_)
211#endif
212 for(size_t i = 0; i < outputVertexList_.size(); i++) {
213 outputVertexList->SetPoint(i, outputVertexList_[i].p_[0],
214 outputVertexList_[i].p_[1],
215 outputVertexList_[i].p_[2]);
216 if(RangeCoordinates) {
217 outputU->SetTuple1(i, outputVertexList_[i].uv_.first);
218 outputV->SetTuple1(i, outputVertexList_[i].uv_.second);
219 }
220 if(EdgeParameterization) {
221 outputParameterization->SetTuple1(i, outputVertexList_[i].t_);
222 }
223 }
224 if(RangeCoordinates) {
225 output->GetPointData()->AddArray(outputU);
226 output->GetPointData()->AddArray(outputV);
227 } else {
228 output->GetPointData()->RemoveArray(dataUfield->GetName());
229 output->GetPointData()->RemoveArray(dataVfield->GetName());
230 }
231 if(EdgeParameterization) {
232 output->GetPointData()->AddArray(outputParameterization);
233 } else {
234 output->GetPointData()->RemoveArray("EdgeParameterization");
235 }
236
237 if(EdgeIds) {
238 outputEdgeIds->SetName("EdgeIds");
239 outputEdgeIds->SetNumberOfTuples(triangleNumber);
240 }
241
242 if(TetIds) {
243 outputTetIds->SetName("TetIds");
244 outputTetIds->SetNumberOfTuples(triangleNumber);
245 }
246
247 if(CaseIds) {
248 outputCaseIds->SetName("CaseIds");
249 outputCaseIds->SetNumberOfTuples(triangleNumber);
250 }
251
252 vtkNew<vtkIdList> idList{};
253 idList->SetNumberOfIds(3);
254
255 triangleNumber = 0;
256 for(size_t i = 0; i < threadedTriangleList_.size(); i++) {
257 for(size_t j = 0; j < threadedTriangleList_[i].size(); j++) {
258 for(int k = 0; k < 3; k++) {
259 idList->SetId(k, threadedTriangleList_[i][j].vertexIds_[k]);
260 }
261 outputTriangleList->InsertNextCell(idList);
262 if(EdgeIds) {
263 outputEdgeIds->SetTuple1(triangleNumber, i);
264 }
265 if(TetIds) {
266 outputTetIds->SetTuple1(
267 triangleNumber, threadedTriangleList_[i][j].tetId_);
268 }
269 if(CaseIds) {
270 outputCaseIds->SetTuple1(
271 triangleNumber, threadedTriangleList_[i][j].caseId_);
272 }
273 triangleNumber++;
274 }
275 }
276 output->SetPolys(outputTriangleList);
277 if(EdgeIds) {
278 output->GetCellData()->AddArray(outputEdgeIds);
279 } else {
280 output->GetCellData()->RemoveArray("EdgeIds");
281 }
282 if(TetIds) {
283 output->GetCellData()->AddArray(outputTetIds);
284 } else {
285 output->GetCellData()->RemoveArray("TetIds");
286 }
287 if(CaseIds) {
288 output->GetCellData()->AddArray(outputCaseIds);
289 } else {
290 output->GetCellData()->RemoveArray("CaseIds");
291 }
292
293 return 1;
294}
#define ttkTemplateMacro(triangulationType, call)
#define ttkNotUsed(x)
Mark function/method parameters that are not used in the function body at all.
Definition BaseClass.h:47
ttk::Triangulation * GetTriangulation(vtkDataSet *dataSet)
TTK VTK-filter that computes fiber surfaces.
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
int FillInputPortInformation(int port, vtkInformation *info) override
int FillOutputPortInformation(int port, vtkInformation *info) override
int dispatch(ttk::Triangulation *const triangulation)
static void * GetVoidPointer(vtkDataArray *array, vtkIdType start=0)
Definition ttkUtils.cpp:226
int printErr(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition Debug.h:149
int setGlobalVertexList(std::vector< Vertex > *globalList)
int setInputField(const void *uField, const void *vField)
int setPolygonEdgeNumber(const SimplexId &polygonEdgeNumber)
int setPointMerging(const bool &onOff)
int setPointMergingThreshold(const double &threshold)
int setPolygon(const std::vector< std::pair< std::pair< double, double >, std::pair< double, double > > > *polygon)
int setVertexList(const SimplexId &polygonEdgeId, std::vector< Vertex > *vertexList)
int setTriangleList(const SimplexId &polygonEdgeId, std::vector< Triangle > *triangleList)
void preconditionTriangulation(AbstractTriangulation *triangulation)
Triangulation is a class that provides time and memory efficient traversal methods on triangulations ...
AbstractTriangulation * getData()
Triangulation::Type getType() const
int SimplexId
Identifier type for simplices of any dimension.
Definition DataTypes.h:22
vtkStandardNewMacro(ttkFiberSurface)
printMsg(debug::output::BOLD+" | | | | | . \\ | | (__| | / __/| |_| / __/|__ _|"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)