TTK
Loading...
Searching...
No Matches
ttkRangePolygon.cpp
Go to the documentation of this file.
1#include <ttkRangePolygon.h>
2
3#include <vtkCleanPolyData.h>
4#include <vtkDataSetSurfaceFilter.h>
5#include <vtkDataSetTriangleFilter.h>
6#include <vtkFeatureEdges.h>
7#include <vtkInformation.h>
8#include <vtkPointData.h>
9#include <vtkUnstructuredGrid.h>
10
11#include <ttkMacros.h>
12#include <ttkUtils.h>
13
14using namespace std;
15using namespace ttk;
16
18
20
21 this->SetNumberOfInputPorts(1);
22 this->SetNumberOfOutputPorts(1);
23
24 setDebugMsgPrefix("RangePolygon");
25
26 vtkWarningMacro("`TTK RangePolygon' is now deprecated. Please use instead "
27 "`Poly Line Source' followed by `Resample With Dataset'.");
28}
29
31
32int ttkRangePolygon::FillInputPortInformation(int port, vtkInformation *info) {
33 if(port == 0)
34 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkUnstructuredGrid");
35 else
36 return 0;
37
38 return 1;
39}
40
41int ttkRangePolygon::FillOutputPortInformation(int port, vtkInformation *info) {
42 if(port == 0)
43 info->Set(vtkDataObject::DATA_TYPE_NAME(), "vtkUnstructuredGrid");
44 else
45 return 0;
46
47 return 1;
48}
49
50int ttkRangePolygon::RequestData(vtkInformation *ttkNotUsed(request),
51 vtkInformationVector **inputVector,
52 vtkInformationVector *outputVector) {
53
54 vtkUnstructuredGrid *input = vtkUnstructuredGrid::GetData(inputVector[0]);
55 vtkUnstructuredGrid *output = vtkUnstructuredGrid::GetData(outputVector);
56
57 if(input->GetNumberOfCells()) {
58
59 if(input->GetCell(0)->GetCellDimension() == 0) {
60 processPoints(input, output);
61 } else {
62 processTriangles(input, output);
63 }
64 } else {
65 processPoints(input, output);
66 }
67
69
70 return 1;
71}
72
73int ttkRangePolygon::processPoints(vtkUnstructuredGrid *input,
74 vtkUnstructuredGrid *output) {
75
76 Timer t;
77
79 output->SetPoints(pointSet);
80
81 output->GetPoints()->ShallowCopy(input->GetPoints());
82 output->GetPointData()->ShallowCopy(input->GetPointData());
83
84 vtkSmartPointer<vtkCellArray> const edgeArray
87 idList->SetNumberOfIds(2);
88
89 for(SimplexId i = 0; i < input->GetNumberOfPoints(); i++) {
90 if(i) {
91 idList->SetId(0, i - 1);
92 idList->SetId(1, i);
93
94 edgeArray->InsertNextCell(idList);
95 }
96 }
97 if(ClosedLoop) {
98 idList->SetId(0, input->GetNumberOfPoints() - 1);
99 idList->SetId(1, 0);
100 edgeArray->InsertNextCell(idList);
101 }
102
103 output->SetCells(VTK_LINE, edgeArray);
104
105 printMsg(std::to_string(output->GetNumberOfCells()) + " edges extracted", 1,
107
108 return 0;
109}
110
111int ttkRangePolygon::processTriangles(vtkUnstructuredGrid *input,
112 vtkUnstructuredGrid *output) {
113
114 Timer t;
115
118
119 surfaceMaker->SetInputData(input);
120
121 vtkSmartPointer<vtkCleanPolyData> const surfaceCleaner
123
124 surfaceCleaner->SetInputConnection(surfaceMaker->GetOutputPort());
125
126 vtkSmartPointer<vtkFeatureEdges> const featureEdges
128
129 featureEdges->SetBoundaryEdges(true);
130 featureEdges->SetManifoldEdges(false);
131 featureEdges->SetFeatureEdges(false);
132 featureEdges->SetNonManifoldEdges(false);
133 featureEdges->SetColoring(false);
134 featureEdges->SetInputConnection(surfaceCleaner->GetOutputPort());
135
138
139 triangleMaker->SetInputConnection(featureEdges->GetOutputPort());
140 triangleMaker->Update();
141
142 output->ShallowCopy(triangleMaker->GetOutput());
143
144 if(NumberOfIterations > 0) {
145
146 // set up the triangulation
147 Triangulation *triangulation = ttkAlgorithm::GetTriangulation(output);
148
149 ScalarFieldSmoother smoother;
150 smoother.setDimensionNumber(3);
151 smoother.setInputDataPointer(output->GetPoints()->GetVoidPointer(0));
152 smoother.setOutputDataPointer(output->GetPoints()->GetVoidPointer(0));
153 smoother.preconditionTriangulation(triangulation);
154
155 switch(output->GetPoints()->GetDataType()) {
156 vtkTemplateMacro(
157 smoother.smooth<VTK_TT>(triangulation, NumberOfIterations));
158 }
159
160 for(int i = 0; i < output->GetPointData()->GetNumberOfArrays(); i++) {
161 vtkDataArray *field = output->GetPointData()->GetArray(i);
162
163 smoother.setDebugLevel(0);
164
165 smoother.setDimensionNumber(field->GetNumberOfComponents());
166 smoother.setInputDataPointer(field->GetVoidPointer(0));
167
168 smoother.setOutputDataPointer(field->GetVoidPointer(0));
169 smoother.preconditionTriangulation(triangulation);
170
171 switch(field->GetDataType()) {
172 vtkTemplateMacro(
173 smoother.smooth<VTK_TT>(triangulation, NumberOfIterations));
174 }
175 }
176 }
177
178 printMsg(std::to_string(output->GetNumberOfCells()) + " edges extracted", 1,
180
181 return 0;
182}
#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 which produces a valid range polygon for fiber surface extraction.
int FillOutputPortInformation(int port, vtkInformation *info) override
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
~ttkRangePolygon() override
int FillInputPortInformation(int port, vtkInformation *info) override
void setDebugMsgPrefix(const std::string &prefix)
Definition Debug.h:364
virtual int setDebugLevel(const int &debugLevel)
Definition Debug.cpp:147
TTK processing package for scalar field smoothing.
int preconditionTriangulation(AbstractTriangulation *triangulation)
void setDimensionNumber(const int &dimensionNumber)
void setOutputDataPointer(void *data)
int smooth(const triangulationType *triangulation, const int &numberOfIterations) const
void setInputDataPointer(void *data)
double getElapsedTime()
Definition Timer.h:15
Triangulation is a class that provides time and memory efficient traversal methods on triangulations ...
The Topology ToolKit.
int SimplexId
Identifier type for simplices of any dimension.
Definition DataTypes.h:22
vtkStandardNewMacro(ttkRangePolygon)
printMsg(debug::output::BOLD+" | | | | | . \\ | | (__| | / __/| |_| / __/|__ _|"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)