TTK
Loading...
Searching...
No Matches
ttkAlgorithm.h
Go to the documentation of this file.
1
14
15#pragma once
16
17// VTK Module
18#include <ttkAlgorithmModule.h>
19
20// VTK Includes
21#include <vtkAlgorithm.h>
22class vtkDataSet;
23class vtkInformation;
24class vtkInformationIntegerKey;
25
26// Base Includes
27#include <Debug.h>
28
29namespace ttk {
30 class Triangulation;
31}
32
33class TTKALGORITHM_EXPORT ttkAlgorithm : public vtkAlgorithm,
34 virtual public ttk::Debug {
35private:
36 int ThreadNumber{1};
37 bool UseAllCores{true};
38
39public:
40 static ttkAlgorithm *New();
41 vtkTypeMacro(ttkAlgorithm, vtkAlgorithm);
42
48 // update ttk::BaseClass member
49 this->setThreadNumber(this->UseAllCores ? ttk::OsCall::getNumberOfCores()
50 : this->ThreadNumber);
51 this->Modified();
52 }
53
58 void SetThreadNumber(int threadNumber) {
59 this->ThreadNumber = threadNumber;
60 this->UpdateThreadNumber();
61 }
62
66 void SetUseAllCores(bool useAllCores) {
67 this->UseAllCores = useAllCores;
68 this->UpdateThreadNumber();
69 }
70
75 void SetDebugLevel(int debugLevel) {
76 this->setDebugLevel(debugLevel); // from ttk::Debug
77 this->Modified();
78 }
79
83 void SetCompactTriangulationCacheSize(float cacheSize) {
84 this->CompactTriangulationCacheSize = cacheSize;
85 this->Modified();
86 }
87
98 vtkDataArray *GetOptionalArray(const bool &enforceArrayIndex,
99 const int &arrayIndex,
100 const std::string &arrayName,
101 vtkDataSet *const inputData,
102 const int &inputPort = 0);
103
108 static std::string GetOrderArrayName(vtkDataArray *const array);
109
121 vtkDataArray *GetOrderArray(vtkDataSet *const inputData,
122 const int scalarArrayIdx,
123 ttk::Triangulation *triangulation,
124 const bool getGlobalOrder = false,
125 const int orderArrayIdx = 0,
126 const bool enforceOrderArrayIdx = false);
127
132 vtkDataArray *
133 checkForGlobalAndComputeOrderArray(vtkDataSet *const inputData,
134 vtkDataArray *scalarArray,
135 const int scalarArrayIdx,
136 const bool getGlobalOrder,
137 vtkDataArray *orderArray,
138 ttk::Triangulation *triangulation,
139 const bool enforceOrderArrayIdx);
140
148 vtkDataArray *ComputeOrderArray(vtkDataSet *const inputData,
149 vtkDataArray *scalarArray,
150 const int scalarArrayIdx,
151 const bool getGlobalOrder,
152 vtkDataArray *oldOrderArray,
153 ttk::Triangulation *triangulation);
166 GetIdentifierArrayPtr(const bool &enforceArrayIndex,
167 const int &arrayIndex,
168 const std::string &arrayName,
169 vtkDataSet *const inputData,
170 std::vector<ttk::SimplexId> &spareStorage,
171 const int inputPort = 0,
172 const bool printErr = true);
173
195 ttk::Triangulation *GetTriangulation(vtkDataSet *dataSet);
196
202 static vtkInformationIntegerKey *SAME_DATA_TYPE_AS_INPUT_PORT();
203
212 int ProcessRequest(vtkInformation *request,
213 vtkInformationVector **inputVectors,
214 vtkInformationVector *outputVector) override;
215
219 vtkDataSet *GetOutput();
220 vtkDataSet *GetOutput(int);
221
227 void SetInputData(vtkDataSet *);
228 void SetInputData(int, vtkDataSet *);
229
235 void AddInputData(vtkDataSet *);
236 void AddInputData(int, vtkDataSet *);
237
249 template <typename inputType>
250 inline int checkEmptyMPIInput(inputType *input) {
251 if(!input) {
252#ifdef TTK_ENABLE_MPI
253 if(ttk::isRunningWithMPI()) {
254 return 1;
255 } else {
256#endif
257 return 0;
258#ifdef TTK_ENABLE_MPI
259 }
260#endif
261 }
262 return 2;
263 };
264
265protected:
267 ~ttkAlgorithm() override;
268 float CompactTriangulationCacheSize{0.2f};
269
270#ifdef TTK_ENABLE_MPI
279 int updateMPICommunicator(vtkDataSet *input);
280#endif
287 void MPIGhostPipelinePreconditioning(vtkDataSet *input);
288
294 void MPIPipelinePreconditioning(vtkDataSet *input,
295 std::vector<int> &neighbors,
296 std::map<int, int> &neighToId,
297 ttk::Triangulation *triangulation = nullptr);
298
308 ttk::SimplexId simplexNumber,
309 unsigned char *ghost,
310 int *rankArray);
318 vtkDataSet *input,
319 std::unordered_map<ttk::SimplexId, ttk::SimplexId> &vertGtoL,
320 std::vector<int> &neighborRanks,
321 std::map<int, int> &neighborsToId);
322
330 vtkDataSet *input);
331
339 virtual int RequestDataObject(vtkInformation *request,
340 vtkInformationVector **inputVectors,
341 vtkInformationVector *outputVector);
342
353 virtual int
354 RequestInformation(vtkInformation *ttkNotUsed(request),
355 vtkInformationVector **ttkNotUsed(inputVectors),
356 vtkInformationVector *ttkNotUsed(outputVector)) {
357 return 1;
358 }
359
366 virtual int
367 RequestUpdateTime(vtkInformation *ttkNotUsed(request),
368 vtkInformationVector **ttkNotUsed(inputVectors),
369 vtkInformationVector *ttkNotUsed(outputVector)) {
370 return 1;
371 }
372
380 vtkInformation *ttkNotUsed(request),
381 vtkInformationVector **ttkNotUsed(inputVectors),
382 vtkInformationVector *ttkNotUsed(outputVector)) {
383 return 1;
384 }
385
394 virtual int
395 RequestUpdateExtent(vtkInformation *ttkNotUsed(request),
396 vtkInformationVector **ttkNotUsed(inputVectors),
397 vtkInformationVector *ttkNotUsed(outputVector)) {
398 return 1;
399 }
400
408 virtual int
409 RequestDataNotGenerated(vtkInformation *ttkNotUsed(request),
410 vtkInformationVector **ttkNotUsed(inputVectors),
411 vtkInformationVector *ttkNotUsed(outputVector)) {
412 return 1;
413 }
414
423 virtual int RequestData(vtkInformation *ttkNotUsed(request),
424 vtkInformationVector **ttkNotUsed(inputVectors),
425 vtkInformationVector *ttkNotUsed(outputVector)) {
426 return 1;
427 }
428
438 vtkInformation *ttkNotUsed(info)) override {
439 return 0;
440 }
441
453 vtkInformation *ttkNotUsed(info)) override {
454 return 0;
455 }
456};
457
#define ttkNotUsed(x)
Mark function/method parameters that are not used in the function body at all.
Definition BaseClass.h:47
Baseclass of all VTK filters that wrap ttk modules.
void SetCompactTriangulationCacheSize(float cacheSize)
static vtkInformationIntegerKey * SAME_DATA_TYPE_AS_INPUT_PORT()
int checkEmptyMPIInput(inputType *input)
This method tests whether the input is a nullptr. If the computation is being done on multiple proces...
bool checkGlobalIdValidity(ttk::LongSimplexId *globalIds, ttk::SimplexId simplexNumber, unsigned char *ghost, int *rankArray)
virtual int RequestUpdateTime(vtkInformation *ttkNotUsed(request), vtkInformationVector **ttkNotUsed(inputVectors), vtkInformationVector *ttkNotUsed(outputVector))
~ttkAlgorithm() override
virtual int RequestDataNotGenerated(vtkInformation *ttkNotUsed(request), vtkInformationVector **ttkNotUsed(inputVectors), vtkInformationVector *ttkNotUsed(outputVector))
virtual int RequestData(vtkInformation *ttkNotUsed(request), vtkInformationVector **ttkNotUsed(inputVectors), vtkInformationVector *ttkNotUsed(outputVector))
void MPIGhostPipelinePreconditioning(vtkDataSet *input)
void SetThreadNumber(int threadNumber)
virtual int RequestInformation(vtkInformation *ttkNotUsed(request), vtkInformationVector **ttkNotUsed(inputVectors), vtkInformationVector *ttkNotUsed(outputVector))
int GenerateGlobalIds(vtkDataSet *input, std::unordered_map< ttk::SimplexId, ttk::SimplexId > &vertGtoL, std::vector< int > &neighborRanks, std::map< int, int > &neighborsToId)
virtual int RequestUpdateTimeDependentInformation(vtkInformation *ttkNotUsed(request), vtkInformationVector **ttkNotUsed(inputVectors), vtkInformationVector *ttkNotUsed(outputVector))
void SetDebugLevel(int debugLevel)
void SetUseAllCores(bool useAllCores)
virtual int RequestUpdateExtent(vtkInformation *ttkNotUsed(request), vtkInformationVector **ttkNotUsed(inputVectors), vtkInformationVector *ttkNotUsed(outputVector))
void MPITriangulationPreconditioning(ttk::Triangulation *triangulation, vtkDataSet *input)
int FillInputPortInformation(int ttkNotUsed(port), vtkInformation *ttkNotUsed(info)) override
void UpdateThreadNumber()
static ttkAlgorithm * New()
void MPIPipelinePreconditioning(vtkDataSet *input, std::vector< int > &neighbors, std::map< int, int > &neighToId, ttk::Triangulation *triangulation=nullptr)
int FillOutputPortInformation(int ttkNotUsed(port), vtkInformation *ttkNotUsed(info)) override
virtual int setThreadNumber(const int threadNumber)
Definition BaseClass.h:80
Minimalist debugging class.
Definition Debug.h:88
virtual int setDebugLevel(const int &debugLevel)
Definition Debug.cpp:147
static int getNumberOfCores()
Definition Os.cpp:87
Triangulation is a class that provides time and memory efficient traversal methods on triangulations ...
The Topology ToolKit.
long long int LongSimplexId
Identifier type for simplices of any dimension.
Definition DataTypes.h:15
int SimplexId
Identifier type for simplices of any dimension.
Definition DataTypes.h:22