TTK
Loading...
Searching...
No Matches
ttkContourAroundPoint.cpp
Go to the documentation of this file.
2
3#include <ttkMacros.h>
4#include <ttkUtils.h>
5
6#include <vtkInformation.h>
7#include <vtkInformationVector.h>
8#include <vtkVersion.h>
9
10#include <vtkCellData.h>
11#include <vtkDataSet.h>
12#include <vtkUnstructuredGrid.h>
13
14#include <vtkFloatArray.h>
15#include <vtkIdTypeArray.h>
16
17#include <cassert>
18
20
21//----------------------------------------------------------------------------//
22
24 vtkInformation *info) {
25 switch(port) {
26 case 0:
27 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkDataSet");
28 return 1;
29 case 1:
30 case 2:
31 info->Set(
32 vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkUnstructuredGrid");
33 return 1;
34 }
35 return 0;
36}
37
39 vtkInformation *info) {
40 switch(port) {
41 case 0:
42 case 1:
43 info->Set(vtkDataObject::DATA_TYPE_NAME(), "vtkUnstructuredGrid");
44 return 1;
45 }
46 return 0;
47}
48
49//----------------------------------------------------------------------------//
50
52 vtkInformationVector **iVec,
53 vtkInformationVector *oVec) {
54
55 _outFld = vtkUnstructuredGrid::GetData(oVec, 0);
56 _outPts = vtkUnstructuredGrid::GetData(oVec, 1);
57
58 if(!preprocessFld(vtkDataSet::GetData(iVec[0])))
59 return 0;
60 if(!preprocessPts(vtkUnstructuredGrid::GetData(iVec[1]),
61 vtkUnstructuredGrid::GetData(iVec[2])))
62 return 0;
63 if(!process())
64 return 0;
65 if(!postprocess())
66 return 0;
67
68 return 1;
69}
70
71//----------------------------------------------------------------------------//
72
73bool ttkContourAroundPoint::preprocessFld(vtkDataSet *dataset) {
74 ttk::Triangulation *triangulation = ttkAlgorithm::GetTriangulation(dataset);
75 if(!triangulation)
76 return false;
77
78 auto scalars = GetInputArrayToProcess(0, dataset);
79 if(!scalars)
80 return false;
81
82 const double radius = ui_spherical ? -1. : 0.;
83
84 const auto errorCode = this->setInputField(
85 triangulation, ttkUtils::GetVoidPointer(scalars), ui_sizeFilter, radius);
86 if(errorCode < 0) {
87 printErr("super->setInputField failed with code "
88 + std::to_string(errorCode));
89 return false;
90 }
91
92 _triangTypeCode = triangulation->getType();
93 _scalarTypeCode = scalars->GetDataType();
94 _scalarsName = scalars->GetName();
95
96 std::ostringstream stream;
97 stream << "Scalar type: " << scalars->GetDataTypeAsString() << " (code "
98 << _scalarTypeCode << ")";
100
101 return true;
102}
103
104//----------------------------------------------------------------------------//
105
106bool ttkContourAroundPoint::preprocessPts(vtkUnstructuredGrid *nodes,
107 vtkUnstructuredGrid *arcs) {
108 // ---- Point data ---- //
109
110 auto points = nodes->GetPoints();
111 if(points->GetDataType() != VTK_FLOAT) {
112 printErr("The point coordinates must be of type float");
113 return false;
114 }
115 auto coords
116 = reinterpret_cast<float *>(ttkUtils::GetVoidPointer(points->GetData()));
117
118 auto pData = nodes->GetPointData();
119 auto scalarBuf = getBuffer<float>(pData, "Scalar", VTK_FLOAT, "float");
120 auto codeBuf = getBuffer<int>(pData, "CriticalType", VTK_INT, "int");
121 if(!scalarBuf || !codeBuf)
122 return false;
123
124 // ---- Cell data ---- //
125
126#ifndef NDEBUG // each arc should of course be defined by exactly two vertices
127 auto cells = arcs->GetCells();
128 const auto maxNvPerC = cells->GetMaxCellSize();
129 if(maxNvPerC != 2) {
130 printErr(
131 "The points must come in pairs but there is at least one cell with "
132 + std::to_string(maxNvPerC) + " points");
133 return false;
134 }
135 // NOTE Ideally check for minNvPerC != 2
136#endif
137
138 auto cData = arcs->GetCellData();
139 auto c2p = getBuffer<int>(cData, "upNodeId", VTK_INT, "int");
140 auto c2q = getBuffer<int>(cData, "downNodeId", VTK_INT, "int");
141 if(!c2p || !c2q)
142 return false;
143
144 // ---- Loop over pairs ---- //
145
146 static constexpr int minCode = 0;
147 static constexpr int maxCode = 3;
148 const double sadFac = ui_extension * 0.01; // saddle or mean of min and max
149 const double extFac = 1 - sadFac; // factor for the extreme point
150
151 _coords.resize(0);
152 _scalars.resize(0);
153 _isovals.resize(0);
154 _flags.resize(0);
155
156 auto addPoint = [this, coords, scalarBuf](int p, float isoval, int code) {
157 const auto point = &coords[p * 3];
158 _coords.push_back(point[0]);
159 _coords.push_back(point[1]);
160 _coords.push_back(point[2]);
161 _scalars.push_back(scalarBuf[p]);
162 _isovals.push_back(isoval);
163 _flags.push_back(code == minCode ? 0 : 1);
164 };
165
166 const vtkIdType nc = arcs->GetNumberOfCells();
167 for(vtkIdType c = 0; c < nc; ++c) {
168 const auto p = c2p[c];
169 const auto q = c2q[c];
170 const auto pCode = codeBuf[p];
171 const auto qCode = codeBuf[q];
172
173 const bool pIsSad = pCode != minCode && pCode != maxCode;
174 const bool qIsSad = qCode != minCode && qCode != maxCode;
175 if(pIsSad || qIsSad) {
176 if(pIsSad && qIsSad) // two saddles
177 continue;
178 // extremum and saddle
179 const auto ext = pIsSad ? q : p;
180 const auto sad = pIsSad ? p : q;
181 const float isoval = scalarBuf[ext] * extFac + scalarBuf[sad] * sadFac;
182 addPoint(ext, isoval, codeBuf[ext]);
183 } else { // min-max pair
184 printWrn("Arc " + std::to_string(c) + " joins a minimum and a maximum");
185 const auto pVal = scalarBuf[p];
186 const auto qVal = scalarBuf[q];
187 const auto cVal = (pVal + qVal) / 2;
188 addPoint(p, pVal * extFac + cVal * sadFac, pCode);
189 addPoint(q, qVal * extFac + cVal * sadFac, qCode);
190 }
191 }
192
193 auto np = _scalars.size();
194 const auto errorCode = this->setInputPoints(
195 _coords.data(), _scalars.data(), _isovals.data(), _flags.data(), np);
196 if(errorCode < 0) {
197 printErr("setInputPoints failed with code " + std::to_string(errorCode));
198 return false;
199 }
200 return true;
201}
202
203//----------------------------------------------------------------------------//
204
206 int errorCode = 0;
207 switch(_scalarTypeCode) {
208 vtkTemplateMacro((errorCode = this->execute<VTK_TT>()));
209 }
210 // ttkVtkTemplateMacro(
211 // _scalarTypeCode, _triangTypeCode,
212 // (errorCode = this->execute<VTK_TT, TTK_TT>())
213 // )
214 if(errorCode < 0) { // In TTK, negative is bad.
215 printErr("super->execute failed with code " + std::to_string(errorCode));
216 return false;
217 }
218 return true;
219}
220
221//----------------------------------------------------------------------------//
222
225 if(nc == 0) // very fine area filter
226 return true;
227
228 // ---- Cell data (output 0) ---- //
229
230 std::vector<int> ctypes(nc);
231
232 vtkIdType cinfoCounter = 0;
233 for(ttk::SimplexId c = 0; c < nc; ++c) {
234 const auto nvOfCell = _outContoursCinfos[cinfoCounter];
235 assert(nvOfCell >= 2 && nvOfCell <= 3); // ensured in super class
236 ctypes[c] = nvOfCell == 2 ? VTK_LINE : VTK_TRIANGLE;
237 cinfoCounter += nvOfCell + 1;
238 }
239
241 auto cinfoArr = vtkSmartPointer<vtkIdTypeArray>::New();
243 cinfoArr, _outContoursCinfos.data(), _outContoursCinfos.size(), 1);
244 cells->SetCells(nc, cinfoArr);
245 _outFld->SetCells(ctypes.data(), cells);
246
247 // ---- Point data (output 0) ---- //
248
249 if(vtkSmartPointer<vtkPoints>::New()->GetDataType() != VTK_FLOAT) {
250 printErr("The API has changed! We have expected the default "
251 "coordinate type to be float");
252 return false;
253 }
254
255 auto points = vtkSmartPointer<vtkPoints>::New();
256 auto coordArr = vtkSmartPointer<vtkFloatArray>::New();
257 coordArr->SetNumberOfComponents(3);
259 coordArr, _outContoursCoords.data(), _outContoursCoords.size(), 1);
260 points->SetData(coordArr);
261 _outFld->SetPoints(points);
262
263 auto scalarArr = vtkFloatArray::New();
265 scalarArr, _outContoursScalars.data(), _outContoursScalars.size(), 1);
266 scalarArr->SetName(_scalarsName);
267 _outFld->GetPointData()->AddArray(scalarArr);
268
269 auto flagArr = vtkIntArray::New();
271 flagArr, _outContoursFlags.data(), _outContoursFlags.size(), 1);
272 flagArr->SetName("isMax");
273 _outFld->GetPointData()->AddArray(flagArr);
274
275 // ---- Output 1 (added in a later revision of the algo) ---- //
276
279 coordArr->SetNumberOfComponents(3);
281 coordArr, _outCentroidsCoords.data(), _outCentroidsCoords.size(), 1);
282 points->SetData(coordArr);
283 _outPts->SetPoints(points);
284
285 scalarArr = vtkFloatArray::New();
287 scalarArr, _outCentroidsScalars.data(), _outCentroidsScalars.size(), 1);
288 scalarArr->SetName(_scalarsName);
289 _outPts->GetPointData()->AddArray(scalarArr);
290
291 flagArr = vtkIntArray::New();
293 flagArr, _outCentroidsFlags.data(), _outCentroidsFlags.size(), 1);
294 flagArr->SetName("isMax");
295 _outPts->GetPointData()->AddArray(flagArr);
296
297 return true;
298}
#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 wraps the contourAroundPoint processing package.
bool postprocess()
Assemble the output object from the results of the TTK module.
int RequestData(vtkInformation *request, vtkInformationVector **iVec, vtkInformationVector *oVec) override
bool preprocessPts(vtkUnstructuredGrid *nodes, vtkUnstructuredGrid *arcs)
int FillOutputPortInformation(int port, vtkInformation *info) override
bool preprocessFld(vtkDataSet *dataset)
int FillInputPortInformation(int port, vtkInformation *info) override
static void * GetVoidPointer(vtkDataArray *array, vtkIdType start=0)
Definition ttkUtils.cpp:226
static void SetVoidArray(vtkDataArray *array, void *data, vtkIdType size, int save)
Definition ttkUtils.cpp:280
int setInputField(triangulationType *triangulation, void *scalars, double sizeFilter, double radius=0.)
std::vector< float > _outCentroidsCoords
std::vector< float > _outContoursScalars
std::vector< float > _outContoursCoords
std::vector< int > _outCentroidsFlags
std::vector< float > _outCentroidsScalars
std::vector< LongSimplexId > _outContoursCinfos
int setInputPoints(float *coords, float *scalars, float *isovals, int *flags, std::size_t np)
std::vector< int > _outContoursFlags
int printWrn(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition Debug.h:159
int printErr(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition Debug.h:149
Triangulation is a class that provides time and memory efficient traversal methods on triangulations ...
Triangulation::Type getType() const
int SimplexId
Identifier type for simplices of any dimension.
Definition DataTypes.h:22
vtkStandardNewMacro(ttkContourAroundPoint)
printMsg(debug::output::BOLD+" | | | | | . \\ | | (__| | / __/| |_| / __/|__ _|"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)