TTK
Loading...
Searching...
No Matches
ttkDiscreteGradient.cpp
Go to the documentation of this file.
2#include <ttkMacros.h>
3#include <ttkUtils.h>
4
5#include <vtkCellData.h>
6#include <vtkIdTypeArray.h>
7#include <vtkInformation.h>
8#include <vtkIntArray.h>
9#include <vtkLine.h>
10#include <vtkPointData.h>
11#include <vtkPolyData.h>
12#include <vtkSignedCharArray.h>
13
15
17 SetNumberOfInputPorts(1);
18 SetNumberOfOutputPorts(2);
19}
20
22 vtkInformation *info) {
23 if(port == 0) {
24 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkDataSet");
25 return 1;
26 }
27 return 0;
28}
29
31 vtkInformation *info) {
32 if(port == 0 || port == 1) {
33 info->Set(vtkDataObject::DATA_TYPE_NAME(), "vtkPolyData");
34 return 1;
35 }
36 return 0;
37}
38
39template <typename scalarType, typename triangulationType>
40int ttkDiscreteGradient::fillCriticalPoints(
41 vtkPolyData *outputCriticalPoints,
42 vtkDataArray *const inputScalars,
43 const triangulationType &triangulation) {
44
45 ttk::Timer tm{};
46
47 // critical points
48 std::vector<std::array<float, 3>> critPoints_coords;
49 std::vector<char> critPoints_cellDimensions;
50 std::vector<SimplexId> critPoints_cellIds;
51 std::vector<char> critPoints_isOnBoundary;
52 std::vector<SimplexId> critPoints_PLVertexIdentifiers;
53
54 this->setCriticalPoints(critPoints_coords, critPoints_cellDimensions,
55 critPoints_cellIds, critPoints_isOnBoundary,
56 critPoints_PLVertexIdentifiers, triangulation);
57 const auto nPoints = critPoints_coords.size();
58 const auto scalars
59 = static_cast<scalarType *>(ttkUtils::GetVoidPointer(inputScalars));
60
61 vtkNew<vtkPoints> points{};
62 points->SetNumberOfPoints(nPoints);
63
64 vtkNew<vtkSignedCharArray> cellDimensions{};
65 cellDimensions->SetNumberOfComponents(1);
66 cellDimensions->SetName("CellDimension");
67 cellDimensions->SetNumberOfTuples(nPoints);
68
69 vtkNew<ttkSimplexIdTypeArray> cellIds{};
70 cellIds->SetNumberOfComponents(1);
71 cellIds->SetName("CellId");
72 cellIds->SetNumberOfTuples(nPoints);
73
74 vtkSmartPointer<vtkDataArray> cellScalars{inputScalars->NewInstance()};
75#ifndef TTK_ENABLE_KAMIKAZE
76 if(!cellScalars) {
77 this->printErr("vtkDataArray allocation problem.");
78 return -1;
79 }
80#endif
81 cellScalars->SetNumberOfComponents(1);
82 cellScalars->SetName(inputScalars->GetName());
83 cellScalars->SetNumberOfTuples(nPoints);
84
85 vtkNew<vtkSignedCharArray> isOnBoundary{};
86 isOnBoundary->SetNumberOfComponents(1);
87 isOnBoundary->SetName("IsOnBoundary");
88 isOnBoundary->SetNumberOfTuples(nPoints);
89
90 vtkNew<ttkSimplexIdTypeArray> PLVertexIdentifiers{};
91 PLVertexIdentifiers->SetNumberOfComponents(1);
92 PLVertexIdentifiers->SetName(ttk::VertexScalarFieldName);
93 PLVertexIdentifiers->SetNumberOfTuples(nPoints);
94
95#ifdef TTK_ENABLE_OPENMP
96#pragma omp parallel for num_threads(this->threadNumber_)
97#endif // TTK_ENABLE_OPENMP
98 for(size_t i = 0; i < nPoints; ++i) {
99 points->SetPoint(i, critPoints_coords[i].data());
100 cellDimensions->SetTuple1(i, critPoints_cellDimensions[i]);
101 cellIds->SetTuple1(i, critPoints_cellIds[i]);
102 cellScalars->SetTuple1(i, scalars[critPoints_PLVertexIdentifiers[i]]);
103 isOnBoundary->SetTuple1(i, critPoints_isOnBoundary[i]);
104#ifdef TTK_ENABLE_MPI
105 if(ttk::hasInitializedMPI()) {
106 PLVertexIdentifiers->SetTuple1(
107 i, triangulation.getVertexGlobalId(critPoints_PLVertexIdentifiers[i]));
108 } else {
109 PLVertexIdentifiers->SetTuple1(i, critPoints_PLVertexIdentifiers[i]);
110 }
111#else
112 PLVertexIdentifiers->SetTuple1(i, critPoints_PLVertexIdentifiers[i]);
113#endif
114 }
115
116 ttkUtils::CellVertexFromPoints(outputCriticalPoints, points);
117
118 vtkPointData *pointData = outputCriticalPoints->GetPointData();
119#ifndef TTK_ENABLE_KAMIKAZE
120 if(!pointData) {
121 this->printErr("outputCriticalPoints has no point data");
122 return -1;
123 }
124#endif
125
126 pointData->SetScalars(cellDimensions);
127 pointData->AddArray(cellIds);
128 pointData->AddArray(cellScalars);
129 pointData->AddArray(isOnBoundary);
130 pointData->AddArray(PLVertexIdentifiers);
131
132 this->printMsg(
133 "Extracted critical points", 1.0, tm.getElapsedTime(), this->threadNumber_);
134
135 return 0;
136}
137
138template <typename triangulationType>
139int ttkDiscreteGradient::fillGradientGlyphs(
140 vtkPolyData *const outputGradientGlyphs,
141 const triangulationType &triangulation) {
142
143 ttk::Timer tm{};
144
145 std::vector<std::array<float, 3>> gradientGlyphs_points;
146 std::vector<char> gradientGlyphs_points_pairOrigins;
147 std::vector<char> gradientGlyphs_cells_pairTypes;
148 std::vector<SimplexId> gradientGlyphs_point_ids{};
149 std::vector<char> gradientGlyphs_point_dimensions{};
150
151 this->setGradientGlyphs(
152 gradientGlyphs_points, gradientGlyphs_points_pairOrigins,
153 gradientGlyphs_cells_pairTypes, gradientGlyphs_point_ids,
154 gradientGlyphs_point_dimensions, triangulation);
155
156 const auto nPoints = gradientGlyphs_points.size();
157
158 vtkNew<vtkPoints> points{};
159 points->SetNumberOfPoints(nPoints);
160 vtkNew<vtkSignedCharArray> pairOrigins{};
161 pairOrigins->SetNumberOfComponents(1);
162 pairOrigins->SetName("PairOrigin");
163 pairOrigins->SetNumberOfTuples(nPoints);
164
165#ifdef TTK_ENABLE_OPENMP
166#pragma omp parallel for num_threads(this->threadNumber_)
167#endif // TTK_ENABLE_OPENMP
168 for(size_t i = 0; i < nPoints; ++i) {
169 points->SetPoint(i, gradientGlyphs_points[i].data());
170 pairOrigins->SetTuple1(i, gradientGlyphs_points_pairOrigins[i]);
171 }
172 outputGradientGlyphs->SetPoints(points);
173
174 const auto nCells = gradientGlyphs_cells_pairTypes.size();
175
176 vtkNew<vtkIdTypeArray> offsets{}, connectivity{};
177 offsets->SetNumberOfComponents(1);
178 offsets->SetNumberOfTuples(nCells + 1);
179 connectivity->SetNumberOfComponents(1);
180 connectivity->SetNumberOfTuples(2 * nCells);
181 vtkNew<vtkSignedCharArray> pairTypes{};
182 pairTypes->SetNumberOfComponents(1);
183 pairTypes->SetName("PairType");
184 pairTypes->SetNumberOfTuples(nCells);
185 vtkNew<ttkSimplexIdTypeArray> cellIds{};
186 cellIds->SetNumberOfComponents(1);
187 cellIds->SetName("CellId");
188 cellIds->SetNumberOfTuples(2 * nCells);
189 vtkNew<vtkSignedCharArray> cellDimensions{};
190 cellDimensions->SetNumberOfComponents(1);
191 cellDimensions->SetName("CellDimension");
192 cellDimensions->SetNumberOfTuples(2 * nCells);
193
194#ifdef TTK_ENABLE_OPENMP
195#pragma omp parallel for num_threads(this->threadNumber_)
196#endif // TTK_ENABLE_OPENMP
197 for(size_t i = 0; i < nCells; ++i) {
198 offsets->SetTuple1(i, 2 * i);
199 // each glyph/line has unique points
200 connectivity->SetTuple1(2 * i, 2 * i);
201 connectivity->SetTuple1(2 * i + 1, 2 * i + 1);
202 pairTypes->SetTuple1(i, gradientGlyphs_cells_pairTypes[i]);
203 cellIds->SetTuple1(2 * i + 0, gradientGlyphs_point_ids[2 * i + 0]);
204 cellIds->SetTuple1(2 * i + 1, gradientGlyphs_point_ids[2 * i + 1]);
205 cellDimensions->SetTuple1(
206 2 * i + 0, gradientGlyphs_point_dimensions[2 * i + 0]);
207 cellDimensions->SetTuple1(
208 2 * i + 1, gradientGlyphs_point_dimensions[2 * i + 1]);
209 }
210 offsets->SetTuple1(nCells, connectivity->GetNumberOfTuples());
211
212 vtkNew<vtkCellArray> cells{};
213 cells->SetData(offsets, connectivity);
214 outputGradientGlyphs->SetLines(cells);
215
216 vtkPointData *pointData = outputGradientGlyphs->GetPointData();
217 vtkCellData *cellData = outputGradientGlyphs->GetCellData();
218
219#ifndef TTK_ENABLE_KAMIKAZE
220 if(pointData == nullptr || cellData == nullptr) {
221 this->printErr("In outputGradientGlyphs point or cell data");
222 return -1;
223 }
224#endif
225
226 pointData->AddArray(pairOrigins);
227 pointData->AddArray(cellIds);
228 pointData->AddArray(cellDimensions);
229 cellData->SetScalars(pairTypes);
230
231 this->printMsg(
232 "Computed gradient glyphs", 1.0, tm.getElapsedTime(), this->threadNumber_);
233
234 return 0;
235}
236
237int ttkDiscreteGradient::RequestData(vtkInformation *ttkNotUsed(request),
238 vtkInformationVector **inputVector,
239 vtkInformationVector *outputVector) {
240
241 auto input = vtkDataSet::GetData(inputVector[0]);
242 auto outputCriticalPoints = vtkPolyData::GetData(outputVector, 0);
243 auto outputGradientGlyphs = vtkPolyData::GetData(outputVector, 1);
244
245 auto triangulation = ttkAlgorithm::GetTriangulation(input);
246
247 int keepGoing = checkEmptyMPIInput<ttk::Triangulation>(triangulation);
248 if(keepGoing < 2) {
249 return keepGoing;
250 }
251#ifndef TTK_ENABLE_KAMIKAZE
252 if(!input) {
253 this->printErr("Input pointer is NULL.");
254 return -1;
255 }
256 if(!outputCriticalPoints or !outputGradientGlyphs) {
257 this->printErr("Output pointer is NULL.");
258 return -1;
259 }
260 if(!input->GetNumberOfPoints()) {
261 this->printErr("Input has no point.");
262 return -1;
263 }
264#endif
265
266 this->preconditionTriangulation(triangulation);
267
268 const auto inputScalars = this->GetInputArrayToProcess(0, input);
269 auto inputOffsets = ttkAlgorithm::GetOrderArray(
270 input, 0, 1, this->ForceInputOffsetScalarField);
271
272 if(inputScalars == nullptr || inputOffsets == nullptr) {
273 this->printErr("Input scalar arrays are NULL");
274 return 0;
275 }
276
277 int ret{};
278
279#ifndef TTK_ENABLE_KAMIKAZE
280 if(inputOffsets->GetDataType() != VTK_INT
281 and inputOffsets->GetDataType() != VTK_ID_TYPE) {
282 this->printErr("input offset field type not supported.");
283 return -1;
284 }
285#endif
286
287 // baseCode processing
289 ttkUtils::GetVoidPointer(inputScalars), inputScalars->GetMTime());
290 this->setInputOffsets(
291 static_cast<SimplexId *>(ttkUtils::GetVoidPointer(inputOffsets)));
292#ifdef TTK_ENABLE_MPI_TIME
293 ttk::Timer t_mpi;
294 ttk::startMPITimer(t_mpi, ttk::MPIrank_, ttk::MPIsize_);
295#endif
296 ttkTemplateMacro(triangulation->getType(),
297 (ret = this->buildGradient<TTK_TT>(
298 *static_cast<TTK_TT *>(triangulation->getData()), true)));
299#ifdef TTK_ENABLE_MPI_TIME
300 double elapsedTime = ttk::endMPITimer(t_mpi, ttk::MPIrank_, ttk::MPIsize_);
301 if(ttk::MPIrank_ == 0) {
302 printMsg("Computation performed using " + std::to_string(ttk::MPIsize_)
303 + " MPI processes lasted :" + std::to_string(elapsedTime));
304 }
305#endif
306 if(ret != 0) {
307 this->printErr("DiscreteGradient.buildGradient() error code: "
308 + std::to_string(ret));
309 return 0;
310 }
311
312 // critical points
313 ttkVtkTemplateMacro(inputScalars->GetDataType(), triangulation->getType(),
314 (fillCriticalPoints<VTK_TT, TTK_TT>(
315 outputCriticalPoints, inputScalars,
316 *static_cast<TTK_TT *>(triangulation->getData()))));
317
318 // gradient glyphs
319 if(ComputeGradientGlyphs) {
320 ttkTemplateMacro(triangulation->getType(),
321 (fillGradientGlyphs<TTK_TT>(
322 outputGradientGlyphs,
323 *static_cast<TTK_TT *>(triangulation->getData()))));
324 }
325
326 return 1;
327}
#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
vtkDataArray * GetOrderArray(vtkDataSet *const inputData, const int scalarArrayIdx, const int orderArrayIdx=0, const bool enforceOrderArrayIdx=false)
ttk::Triangulation * GetTriangulation(vtkDataSet *dataSet)
TTK VTK-filter that wraps the discreteGradient processing package.
int FillInputPortInformation(int port, vtkInformation *info) override
int FillOutputPortInformation(int port, vtkInformation *info) override
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
static void * GetVoidPointer(vtkDataArray *array, vtkIdType start=0)
Definition: ttkUtils.cpp:225
static int CellVertexFromPoints(vtkDataSet *const dataSet, vtkPoints *const points)
Definition: ttkUtils.cpp:327
int printMsg(const std::string &msg, const debug::Priority &priority=debug::Priority::INFO, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cout) const
Definition: Debug.h:118
int printErr(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition: Debug.h:149
void setInputOffsets(const SimplexId *const data)
int setGradientGlyphs(std::vector< std::array< float, 3 > > &points, std::vector< char > &points_pairOrigins, std::vector< char > &cells_pairTypes, std::vector< SimplexId > &cellsIds, std::vector< char > &cellsDimensions, const triangulationType &triangulation) const
void setInputScalarField(const void *const data, const size_t mTime)
void preconditionTriangulation(AbstractTriangulation *const data)
int setCriticalPoints(const std::array< std::vector< SimplexId >, 4 > &criticalCellsByDim, std::vector< std::array< float, 3 > > &points, std::vector< char > &cellDimensions, std::vector< SimplexId > &cellIds, std::vector< char > &isOnBoundary, std::vector< SimplexId > &PLVertexIdentifiers, const triangulationType &triangulation) const
COMMON_EXPORTS int MPIsize_
Definition: BaseClass.cpp:10
COMMON_EXPORTS int MPIrank_
Definition: BaseClass.cpp:9
const char VertexScalarFieldName[]
default name for vertex scalar field
Definition: DataTypes.h:35
vtkStandardNewMacro(ttkDiscreteGradient)
#define ttkVtkTemplateMacro(dataType, triangulationType, call)
Definition: ttkMacros.h:69