5#include <vtkCellData.h>
6#include <vtkIdTypeArray.h>
7#include <vtkImageData.h>
8#include <vtkInformation.h>
9#include <vtkIntArray.h>
11#include <vtkPointData.h>
12#include <vtkPolyData.h>
13#include <vtkSignedCharArray.h>
18 SetNumberOfInputPorts(1);
19 SetNumberOfOutputPorts(2);
23 vtkInformation *info) {
25 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(),
"vtkDataSet");
32 vtkInformation *info) {
33 if(port == 0 || port == 1) {
34 info->Set(vtkDataObject::DATA_TYPE_NAME(),
"vtkPolyData");
40template <
typename scalarType,
typename triangulationType>
41int ttkDiscreteGradient::fillCriticalPoints(
42 vtkPolyData *outputCriticalPoints,
43 vtkDataArray *
const inputScalars,
44 const triangulationType &triangulation) {
49 std::vector<std::array<float, 3>> critPoints_coords;
50 std::vector<char> critPoints_cellDimensions;
51 std::vector<SimplexId> critPoints_cellIds;
52 std::vector<char> critPoints_isOnBoundary;
53 std::vector<SimplexId> critPoints_PLVertexIdentifiers;
56 critPoints_cellIds, critPoints_isOnBoundary,
57 critPoints_PLVertexIdentifiers, triangulation);
58 const auto nPoints = critPoints_coords.size();
62 vtkNew<vtkPoints> points{};
63 points->SetNumberOfPoints(nPoints);
65 vtkNew<vtkSignedCharArray> cellDimensions{};
66 cellDimensions->SetNumberOfComponents(1);
67 cellDimensions->SetName(
"CellDimension");
68 cellDimensions->SetNumberOfTuples(nPoints);
70 vtkNew<ttkSimplexIdTypeArray> cellIds{};
71 cellIds->SetNumberOfComponents(1);
72 cellIds->SetName(
"CellId");
73 cellIds->SetNumberOfTuples(nPoints);
75 vtkSmartPointer<vtkDataArray>
const cellScalars{inputScalars->NewInstance()};
76#ifndef TTK_ENABLE_KAMIKAZE
78 this->
printErr(
"vtkDataArray allocation problem.");
82 cellScalars->SetNumberOfComponents(1);
83 cellScalars->SetName(inputScalars->GetName());
84 cellScalars->SetNumberOfTuples(nPoints);
86 vtkNew<vtkSignedCharArray> isOnBoundary{};
87 isOnBoundary->SetNumberOfComponents(1);
88 isOnBoundary->SetName(
"IsOnBoundary");
89 isOnBoundary->SetNumberOfTuples(nPoints);
91 vtkNew<ttkSimplexIdTypeArray> PLVertexIdentifiers{};
92 PLVertexIdentifiers->SetNumberOfComponents(1);
94 PLVertexIdentifiers->SetNumberOfTuples(nPoints);
96#ifdef TTK_ENABLE_OPENMP
97#pragma omp parallel for num_threads(this->threadNumber_)
99 for(
size_t i = 0; i < nPoints; ++i) {
100 points->SetPoint(i, critPoints_coords[i].data());
101 cellDimensions->SetTuple1(i, critPoints_cellDimensions[i]);
102 cellIds->SetTuple1(i, critPoints_cellIds[i]);
103 cellScalars->SetTuple1(i, scalars[critPoints_PLVertexIdentifiers[i]]);
104 isOnBoundary->SetTuple1(i, critPoints_isOnBoundary[i]);
106 if(ttk::hasInitializedMPI()) {
107 PLVertexIdentifiers->SetTuple1(
108 i, triangulation.getVertexGlobalId(critPoints_PLVertexIdentifiers[i]));
110 PLVertexIdentifiers->SetTuple1(i, critPoints_PLVertexIdentifiers[i]);
113 PLVertexIdentifiers->SetTuple1(i, critPoints_PLVertexIdentifiers[i]);
119 vtkPointData *pointData = outputCriticalPoints->GetPointData();
120#ifndef TTK_ENABLE_KAMIKAZE
122 this->
printErr(
"outputCriticalPoints has no point data");
127 pointData->SetScalars(cellDimensions);
128 pointData->AddArray(cellIds);
129 pointData->AddArray(cellScalars);
130 pointData->AddArray(isOnBoundary);
131 pointData->AddArray(PLVertexIdentifiers);
134 "Extracted critical points", 1.0, tm.getElapsedTime(), this->threadNumber_);
139template <
typename triangulationType>
140int ttkDiscreteGradient::fillGradientGlyphs(
141 vtkPolyData *
const outputGradientGlyphs,
142 const triangulationType &triangulation) {
146 std::vector<std::array<float, 3>> gradientGlyphs_points;
147 std::vector<char> gradientGlyphs_points_pairOrigins;
148 std::vector<char> gradientGlyphs_cells_pairTypes;
149 std::vector<SimplexId> gradientGlyphs_point_ids{};
150 std::vector<char> gradientGlyphs_point_dimensions{};
153 gradientGlyphs_points, gradientGlyphs_points_pairOrigins,
154 gradientGlyphs_cells_pairTypes, gradientGlyphs_point_ids,
155 gradientGlyphs_point_dimensions, triangulation);
157 const auto nPoints = gradientGlyphs_points.size();
159 vtkNew<vtkPoints> points{};
160 points->SetNumberOfPoints(nPoints);
161 vtkNew<vtkSignedCharArray> pairOrigins{};
162 pairOrigins->SetNumberOfComponents(1);
163 pairOrigins->SetName(
"PairOrigin");
164 pairOrigins->SetNumberOfTuples(nPoints);
166#ifdef TTK_ENABLE_OPENMP
167#pragma omp parallel for num_threads(this->threadNumber_)
169 for(
size_t i = 0; i < nPoints; ++i) {
170 points->SetPoint(i, gradientGlyphs_points[i].data());
171 pairOrigins->SetTuple1(i, gradientGlyphs_points_pairOrigins[i]);
173 outputGradientGlyphs->SetPoints(points);
175 const auto nCells = gradientGlyphs_cells_pairTypes.size();
177 vtkNew<vtkIdTypeArray> offsets{}, connectivity{};
178 offsets->SetNumberOfComponents(1);
179 offsets->SetNumberOfTuples(nCells + 1);
180 connectivity->SetNumberOfComponents(1);
181 connectivity->SetNumberOfTuples(2 * nCells);
182 vtkNew<vtkSignedCharArray> pairTypes{};
183 pairTypes->SetNumberOfComponents(1);
184 pairTypes->SetName(
"PairType");
185 pairTypes->SetNumberOfTuples(nCells);
186 vtkNew<ttkSimplexIdTypeArray> cellIds{};
187 cellIds->SetNumberOfComponents(1);
188 cellIds->SetName(
"CellId");
189 cellIds->SetNumberOfTuples(2 * nCells);
190 vtkNew<vtkSignedCharArray> cellDimensions{};
191 cellDimensions->SetNumberOfComponents(1);
192 cellDimensions->SetName(
"CellDimension");
193 cellDimensions->SetNumberOfTuples(2 * nCells);
195#ifdef TTK_ENABLE_OPENMP
196#pragma omp parallel for num_threads(this->threadNumber_)
198 for(
size_t i = 0; i < nCells; ++i) {
199 offsets->SetTuple1(i, 2 * i);
201 connectivity->SetTuple1(2 * i, 2 * i);
202 connectivity->SetTuple1(2 * i + 1, 2 * i + 1);
203 pairTypes->SetTuple1(i, gradientGlyphs_cells_pairTypes[i]);
204 cellIds->SetTuple1(2 * i + 0, gradientGlyphs_point_ids[2 * i + 0]);
205 cellIds->SetTuple1(2 * i + 1, gradientGlyphs_point_ids[2 * i + 1]);
206 cellDimensions->SetTuple1(
207 2 * i + 0, gradientGlyphs_point_dimensions[2 * i + 0]);
208 cellDimensions->SetTuple1(
209 2 * i + 1, gradientGlyphs_point_dimensions[2 * i + 1]);
211 offsets->SetTuple1(nCells, connectivity->GetNumberOfTuples());
213 vtkNew<vtkCellArray> cells{};
214 cells->SetData(offsets, connectivity);
215 outputGradientGlyphs->SetLines(cells);
217 vtkPointData *pointData = outputGradientGlyphs->GetPointData();
218 vtkCellData *cellData = outputGradientGlyphs->GetCellData();
220#ifndef TTK_ENABLE_KAMIKAZE
221 if(pointData ==
nullptr || cellData ==
nullptr) {
222 this->
printErr(
"In outputGradientGlyphs point or cell data");
227 pointData->AddArray(pairOrigins);
228 pointData->AddArray(cellIds);
229 pointData->AddArray(cellDimensions);
230 cellData->SetScalars(pairTypes);
233 "Computed gradient glyphs", 1.0, tm.getElapsedTime(), this->threadNumber_);
239 vtkInformationVector **inputVector,
240 vtkInformationVector *outputVector) {
242 auto input = vtkDataSet::GetData(inputVector[0]);
243 auto outputCriticalPoints = vtkPolyData::GetData(outputVector, 0);
244 auto outputGradientGlyphs = vtkPolyData::GetData(outputVector, 1);
252#ifndef TTK_ENABLE_KAMIKAZE
254 this->
printErr(
"Input pointer is NULL.");
257 if(!outputCriticalPoints or !outputGradientGlyphs) {
258 this->
printErr(
"Output pointer is NULL.");
261 if(!input->GetNumberOfPoints()) {
262 this->
printErr(
"Input has no point.");
269 const auto inputScalars = this->GetInputArrayToProcess(0, input);
271 input, 0, triangulation,
false, 1, this->ForceInputOffsetScalarField);
273 if(inputScalars ==
nullptr || inputOffsets ==
nullptr) {
274 this->
printErr(
"Input scalar arrays are NULL");
280#ifndef TTK_ENABLE_KAMIKAZE
281 if(inputOffsets->GetDataType() != VTK_INT
282 and inputOffsets->GetDataType() != VTK_ID_TYPE) {
283 this->
printErr(
"input offset field type not supported.");
294 const auto imageDataInput = vtkImageData::SafeDownCast(input);
299 if(Backend == 1 && !imageDataInput) {
301 this->
printWrn(
"The stochastic gradient (IEEE TVCG 2012) can only");
302 this->
printWrn(
"be used on vtkImageData (.vti).");
303 this->
printWrn(
"Defaulting to homotopic expansion (IEE PAMI 2011)");
305 if(Backend == 1 && imageDataInput) {
308 this->
setSeed(StochasticGradientSeed);
311 inputScalars->GetDataType(), triangulation->getType(),
312 (ret = this->buildGradient<VTK_TT, TTK_TT>(
313 *
static_cast<TTK_TT *
>(triangulation->getData()),
true,
nullptr)));
316 this->
printErr(
"DiscreteGradient.buildGradient() error code: "
317 + std::to_string(ret));
323 (fillCriticalPoints<VTK_TT, TTK_TT>(
324 outputCriticalPoints, inputScalars,
325 *
static_cast<TTK_TT *
>(triangulation->getData()))));
328 if(ComputeGradientGlyphs) {
330 (fillGradientGlyphs<TTK_TT>(
331 outputGradientGlyphs,
332 *
static_cast<TTK_TT *
>(triangulation->getData()))));
#define ttkTemplateMacro(triangulationType, call)
#define ttkNotUsed(x)
Mark function/method parameters that are not used in the function body at all.
int checkEmptyMPIInput(inputType *input)
This method tests whether the input is a nullptr. If the computation is being done on multiple proces...
ttk::Triangulation * GetTriangulation(vtkDataSet *dataSet)
vtkDataArray * GetOrderArray(vtkDataSet *const inputData, const int scalarArrayIdx, ttk::Triangulation *triangulation, const bool getGlobalOrder=false, const int orderArrayIdx=0, const bool enforceOrderArrayIdx=false)
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)
static int CellVertexFromPoints(vtkDataSet *const dataSet, vtkPoints *const points)
int printWrn(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
int printErr(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
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
void setSeed(const unsigned int &newSeed)
void setBackend(const BACKEND newBackend)
int SimplexId
Identifier type for simplices of any dimension.
const char VertexScalarFieldName[]
default name for vertex scalar field
vtkStandardNewMacro(ttkDiscreteGradient)
#define ttkVtkTemplateMacro(dataType, triangulationType, call)
printMsg(debug::output::BOLD+" | | | | | . \\ | | (__| | / __/| |_| / __/| (_) |"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)