TTK
Loading...
Searching...
No Matches
ttkPeriodicGhostsGeneration.h
Go to the documentation of this file.
1
30
31#pragma once
32
33// VTK Module
34#include <ttkPeriodicGhostsGenerationModule.h>
35
36// VTK Includes
38#include <Triangulation.h>
39#include <ttkAlgorithm.h>
40#include <vtkCellData.h>
41#include <vtkCharArray.h>
42#include <vtkDataSet.h>
43#include <vtkImageData.h>
44#include <vtkInformation.h>
45#include <vtkInformationVector.h>
46#include <vtkPointData.h>
47#include <vtkStreamingDemandDrivenPipeline.h>
48#include <vtkStructuredPoints.h>
49
50#ifdef TTK_ENABLE_MPI
51#include <cstring>
52#include <vtkCommunicator.h>
53#include <vtkExtentTranslator.h>
54#include <vtkExtractVOI.h>
55#endif
56
57namespace ttk {
58
59 namespace periodicGhosts {
61 unsigned char isBound{0};
62 double x{0};
63 double y{0};
64 double z{0};
65 };
66 } // namespace periodicGhosts
67} // namespace ttk
68
69class TTKPERIODICGHOSTSGENERATION_EXPORT ttkPeriodicGhostsGeneration
70 : public ttkAlgorithm {
71
72private:
73#ifdef TTK_ENABLE_MPI
74 std::array<int, 6> outExtent_; // Global extent size before computation
75 std::array<int, 6> inExtent_; // Global extent size after computation
76 std::array<double, 6> boundsWithoutGhosts_; // Local Bounds without ghosts
77 std::array<double, 6> globalBounds_;
78 std::array<double, 3> origin_;
79 std::array<double, 3> spacing_;
80 bool isOutputExtentComputed_{false};
81 std::array<ttk::periodicGhosts::partialGlobalBound, 6> localGlobalBounds_;
82 std::vector<int> neighbors_;
83 std::vector<std::array<ttk::SimplexId, 6>> neighborVertexBBoxes_;
84 std::array<unsigned char, 6> isBoundaryPeriodic_{};
85#endif
86public:
89
91 ~ttkPeriodicGhostsGeneration() override = default;
92
93#ifdef TTK_ENABLE_MPI
94
95 int RequestUpdateExtent(vtkInformation *ttkNotUsed(request),
96 vtkInformationVector **inputVector,
97 vtkInformationVector *outputVector) override;
98
99 int RequestInformation(vtkInformation *request,
100 vtkInformationVector **inputVectors,
101 vtkInformationVector *outputVector) override;
108 int MPIPeriodicGhostPipelinePreconditioning(vtkImageData *imageIn,
109 vtkImageData *imageOut);
110
111 inline std::vector<int> getNeighbors() {
112 std::sort(neighbors_.begin(), neighbors_.end());
113 return neighbors_;
114 }
115
116 inline std::map<int, int> getNeighborsToId() {
117 std::map<int, int> neighborsToId{};
118 int neighborNumber = neighbors_.size();
119 for(int i = 0; i < neighborNumber; i++) {
120 neighborsToId[neighbors_[i]] = i;
121 }
122 return neighborsToId;
123 }
124
125 inline std::array<unsigned char, 6> getIsBoundaryPeriodic() {
126 return isBoundaryPeriodic_;
127 }
128
129#endif
130
131protected:
132 int FillInputPortInformation(int port, vtkInformation *info) override;
133
134 int FillOutputPortInformation(int port, vtkInformation *info) override;
135
136#ifdef TTK_ENABLE_MPI
137
144 int ComputeOutputExtent();
145
165 template <int matchesSize, int metaDataSize>
166 int MarshalAndSendRecv(
167 vtkImageData *imageIn,
168 std::vector<std::vector<vtkSmartPointer<vtkCharArray>>>
169 &charArrayBoundaries,
170 std::vector<std::vector<std::array<ttk::SimplexId, metaDataSize>>>
171 &charArrayBoundariesMetaData,
172 std::vector<std::array<ttk::SimplexId, matchesSize>> &matches,
173 std::vector<vtkSmartPointer<vtkCharArray>> &charArrayBoundariesReceived,
174 std::vector<std::array<ttk::SimplexId, metaDataSize>>
175 &charArrayBoundariesMetaDataReceived,
176 int dim);
177
189 int MergeImageAndSlice(vtkImageData *image,
190 vtkImageData *slice,
191 vtkImageData *mergedImage,
192 int direction);
193
207 template <typename boundaryType>
208 int UnMarshalAndMerge(
209 std::vector<boundaryType> &metaDataReceived,
210 std::vector<vtkSmartPointer<vtkCharArray>> &boundariesReceived,
211 boundaryType direction,
212 int mergeDirection,
213 vtkImageData *mergedImage);
214
227 template <typename boundaryType>
228 int UnMarshalAndCopy(
229 std::vector<boundaryType> &metaDataReceived,
230 std::vector<vtkSmartPointer<vtkCharArray>> &boundariesReceived,
231 boundaryType direction,
232 vtkImageData *mergedImage);
233
249 int MergeDataArrays(vtkDataArray *imageArray,
250 vtkDataArray *sliceArray,
251 vtkSmartPointer<vtkDataArray> &currentArray,
252 int direction,
253 int dims[3],
254 unsigned char ghostValue,
255 ttk::SimplexId numberOfSimplices,
256 ttk::SimplexId numberOfTuples);
257
258#endif
259
260 int RequestData(vtkInformation *request,
261 vtkInformationVector **inputVector,
262 vtkInformationVector *outputVector) override;
263};
#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.
virtual int RequestData(vtkInformation *ttkNotUsed(request), vtkInformationVector **ttkNotUsed(inputVectors), vtkInformationVector *ttkNotUsed(outputVector))
virtual int RequestInformation(vtkInformation *ttkNotUsed(request), vtkInformationVector **ttkNotUsed(inputVectors), vtkInformationVector *ttkNotUsed(outputVector))
virtual int RequestUpdateExtent(vtkInformation *ttkNotUsed(request), vtkInformationVector **ttkNotUsed(inputVectors), vtkInformationVector *ttkNotUsed(outputVector))
int FillInputPortInformation(int ttkNotUsed(port), vtkInformation *ttkNotUsed(info)) override
int FillOutputPortInformation(int ttkNotUsed(port), vtkInformation *ttkNotUsed(info)) override
TTK VTK-filter that generates an outside ghost layer for periodic implicit grids when used in a distr...
~ttkPeriodicGhostsGeneration() override=default
static ttkPeriodicGhostsGeneration * New()
The Topology ToolKit.
int SimplexId
Identifier type for simplices of any dimension.
Definition DataTypes.h:22