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
118 vtkDataArray *GetOrderArray(vtkDataSet *const inputData,
119 const int scalarArrayIdx,
120 const int orderArrayIdx = 0,
121 const bool enforceOrderArrayIdx = false);
122
135 GetIdentifierArrayPtr(const bool &enforceArrayIndex,
136 const int &arrayIndex,
137 const std::string &arrayName,
138 vtkDataSet *const inputData,
139 std::vector<ttk::SimplexId> &spareStorage,
140 const int inputPort = 0,
141 const bool printErr = true);
142
164 ttk::Triangulation *GetTriangulation(vtkDataSet *dataSet);
165
171 static vtkInformationIntegerKey *SAME_DATA_TYPE_AS_INPUT_PORT();
172
181 int ProcessRequest(vtkInformation *request,
182 vtkInformationVector **inputVectors,
183 vtkInformationVector *outputVector) override;
184
188 vtkDataSet *GetOutput();
189 vtkDataSet *GetOutput(int);
190
196 void SetInputData(vtkDataSet *);
197 void SetInputData(int, vtkDataSet *);
198
204 void AddInputData(vtkDataSet *);
205 void AddInputData(int, vtkDataSet *);
206
218 template <typename inputType>
219 inline int checkEmptyMPIInput(inputType *input) {
220 if(!input) {
221#ifdef TTK_ENABLE_MPI
222 if(ttk::isRunningWithMPI()) {
223 return 1;
224 } else {
225#endif
226 return 0;
227#ifdef TTK_ENABLE_MPI
228 }
229#endif
230 }
231 return 2;
232 };
233
234protected:
236 ~ttkAlgorithm() override;
237 float CompactTriangulationCacheSize{0.2f};
238
239#ifdef TTK_ENABLE_MPI
248 int updateMPICommunicator(vtkDataSet *input);
249#endif
256 void MPIGhostPipelinePreconditioning(vtkDataSet *input);
257
263 void MPIPipelinePreconditioning(vtkDataSet *input,
264 std::vector<int> &neighbors,
265 ttk::Triangulation *triangulation = nullptr);
266
276 ttk::SimplexId simplexNumber,
277 unsigned char *ghost,
278 int *rankArray);
286 vtkDataSet *input,
287 std::unordered_map<ttk::SimplexId, ttk::SimplexId> &vertGtoL,
288 std::vector<int> &neighborRanks);
289
297 vtkDataSet *input);
298
306 virtual int RequestDataObject(vtkInformation *request,
307 vtkInformationVector **inputVectors,
308 vtkInformationVector *outputVector);
309
320 virtual int
321 RequestInformation(vtkInformation *ttkNotUsed(request),
322 vtkInformationVector **ttkNotUsed(inputVectors),
323 vtkInformationVector *ttkNotUsed(outputVector)) {
324 return 1;
325 }
326
333 virtual int
334 RequestUpdateTime(vtkInformation *ttkNotUsed(request),
335 vtkInformationVector **ttkNotUsed(inputVectors),
336 vtkInformationVector *ttkNotUsed(outputVector)) {
337 return 1;
338 }
339
347 vtkInformation *ttkNotUsed(request),
348 vtkInformationVector **ttkNotUsed(inputVectors),
349 vtkInformationVector *ttkNotUsed(outputVector)) {
350 return 1;
351 }
352
361 virtual int
362 RequestUpdateExtent(vtkInformation *ttkNotUsed(request),
363 vtkInformationVector **ttkNotUsed(inputVectors),
364 vtkInformationVector *ttkNotUsed(outputVector)) {
365 return 1;
366 }
367
375 virtual int
376 RequestDataNotGenerated(vtkInformation *ttkNotUsed(request),
377 vtkInformationVector **ttkNotUsed(inputVectors),
378 vtkInformationVector *ttkNotUsed(outputVector)) {
379 return 1;
380 }
381
390 virtual int RequestData(vtkInformation *ttkNotUsed(request),
391 vtkInformationVector **ttkNotUsed(inputVectors),
392 vtkInformationVector *ttkNotUsed(outputVector)) {
393 return 1;
394 }
395
405 vtkInformation *ttkNotUsed(info)) override {
406 return 0;
407 }
408
420 vtkInformation *ttkNotUsed(info)) override {
421 return 0;
422 }
423};
424
#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.
Definition: ttkAlgorithm.h:34
void SetCompactTriangulationCacheSize(float cacheSize)
Definition: ttkAlgorithm.h:83
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...
Definition: ttkAlgorithm.h:219
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))
Definition: ttkAlgorithm.h:334
~ttkAlgorithm() override
virtual int RequestDataNotGenerated(vtkInformation *ttkNotUsed(request), vtkInformationVector **ttkNotUsed(inputVectors), vtkInformationVector *ttkNotUsed(outputVector))
Definition: ttkAlgorithm.h:376
virtual int RequestData(vtkInformation *ttkNotUsed(request), vtkInformationVector **ttkNotUsed(inputVectors), vtkInformationVector *ttkNotUsed(outputVector))
Definition: ttkAlgorithm.h:390
void MPIGhostPipelinePreconditioning(vtkDataSet *input)
void SetThreadNumber(int threadNumber)
Definition: ttkAlgorithm.h:58
virtual int RequestInformation(vtkInformation *ttkNotUsed(request), vtkInformationVector **ttkNotUsed(inputVectors), vtkInformationVector *ttkNotUsed(outputVector))
Definition: ttkAlgorithm.h:321
virtual int RequestUpdateTimeDependentInformation(vtkInformation *ttkNotUsed(request), vtkInformationVector **ttkNotUsed(inputVectors), vtkInformationVector *ttkNotUsed(outputVector))
Definition: ttkAlgorithm.h:346
void SetDebugLevel(int debugLevel)
Definition: ttkAlgorithm.h:75
void SetUseAllCores(bool useAllCores)
Definition: ttkAlgorithm.h:66
virtual int RequestUpdateExtent(vtkInformation *ttkNotUsed(request), vtkInformationVector **ttkNotUsed(inputVectors), vtkInformationVector *ttkNotUsed(outputVector))
Definition: ttkAlgorithm.h:362
void MPITriangulationPreconditioning(ttk::Triangulation *triangulation, vtkDataSet *input)
int FillInputPortInformation(int ttkNotUsed(port), vtkInformation *ttkNotUsed(info)) override
Definition: ttkAlgorithm.h:404
void UpdateThreadNumber()
Definition: ttkAlgorithm.h:47
int GenerateGlobalIds(vtkDataSet *input, std::unordered_map< ttk::SimplexId, ttk::SimplexId > &vertGtoL, std::vector< int > &neighborRanks)
static ttkAlgorithm * New()
void MPIPipelinePreconditioning(vtkDataSet *input, std::vector< int > &neighbors, ttk::Triangulation *triangulation=nullptr)
int FillOutputPortInformation(int ttkNotUsed(port), vtkInformation *ttkNotUsed(info)) override
Definition: ttkAlgorithm.h:419
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 ...
Definition: Triangulation.h:48
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