DNDSR 0.2.1
Distributed Numeric Data Structure for CFV
Loading...
Searching...
No Matches
PeriodicInfo.hpp
Go to the documentation of this file.
1#pragma once
2
4#include "DNDS/ArrayPair.hpp"
5#include "Geometric.hpp"
8
9namespace DNDS::Geom
10{
12 {
13 uint8_t _v{0U};
14 [[nodiscard]] bool getP1() const { return _v & 0x01U; }
15 [[nodiscard]] bool getP2() const { return _v & 0x02U; }
16 [[nodiscard]] bool getP3() const { return _v & 0x04U; }
17 void setP1True() { _v |= 0x01U; }
18 void setP2True() { _v |= 0x02U; }
19 void setP3True() { _v |= 0x04U; }
21 DNDS_DEVICE_CALLABLE uint8_t operator^(const NodePeriodicBits &r) const
22 {
23 return uint8_t(_v ^ r._v);
24 }
26 {
27 return NodePeriodicBits{uint8_t(_v & r._v)};
28 }
29 DNDS_DEVICE_CALLABLE operator uint8_t() const
30 {
31 return uint8_t{_v};
32 }
33 DNDS_DEVICE_CALLABLE operator bool() const
34 {
35 return bool(_v);
36 }
38 {
39 return uint8_t(r) == uint8_t(*this);
40 }
41 static MPI_Datatype CommType() { return MPI_UINT8_T; }
42 static int CommMult() { return 1; }
43
44 friend std::ostream &operator<<(std::ostream &o, const NodePeriodicBits &b)
45 {
46 o << int(b._v);
47 return o;
48 }
49 };
50
51 static const NodePeriodicBits nodePB1{0x01U};
52 static const NodePeriodicBits nodePB2{0x02U};
53 static const NodePeriodicBits nodePB3{0x04U};
54
55}
56namespace DNDS
57{
58 // DNDS_DEVICE_STORAGE_BASE_DELETER_INST(Geom::NodePeriodicBits, extern)
59 // DNDS_DEVICE_STORAGE_INST(Geom::NodePeriodicBits, DeviceBackend::Host, extern)
60 // #ifdef DNDS_USE_CUDA
61 // DNDS_DEVICE_STORAGE_INST(Geom::NodePeriodicBits, DeviceBackend::CUDA, extern)
62 // #endif
63}
64namespace DNDS::Geom
65{
66
68 {
71
72 bool operator<(const NodeIndexPBI &r) const
73 {
74 if (i < r.i)
75 return true;
76 else if (i == r.i)
77 return uint8_t(pbi) < uint8_t(r.pbi);
78 else
79 return false;
80 }
81
82 bool operator==(const NodeIndexPBI &r) const
83 {
84 return r.i == i && r.pbi == pbi;
85 }
86
87 bool operator!=(const NodeIndexPBI &r) const { return !(*this == r); }
88 };
89
90 inline bool isCollaborativeNodePeriodicBits(const std::vector<NodePeriodicBits> &a, const std::vector<NodePeriodicBits> &b)
91 {
92 size_t n = a.size();
93 DNDS_assert(n == b.size());
94 if (n == 0)
95 return true;
96 auto v0 = a[0] ^ b[0];
97 for (size_t i = 1; i < n; i++)
98 if ((a.at(i) ^ b.at(i)) != v0)
99 return false;
100 return true;
101 }
102
103 class NodePeriodicBitsRow // instead of std::vector<NodePeriodicBits> for building on raw buffer as a "mapping" object
104 {
105 NodePeriodicBits *_p_indices;
106 rowsize _row_size;
107
108 public:
109 NodePeriodicBitsRow(NodePeriodicBits *ptr, rowsize siz) : _p_indices(ptr), _row_size(siz) {} // default actually
110
111 // --- Special members (cppcoreguidelines-special-member-functions) ---
112 // NodePeriodicBitsRow is a non-owning view (pointer + size).
113 // The copy-assignment operator deep-copies the pointed-to contents;
114 // all other special members are shallow (trivially copyable members).
119
121 {
122 DNDS_assert(j >= 0 && j < _row_size);
123 return _p_indices[j];
124 }
125
127 {
128 DNDS_assert(j >= 0 && j < _row_size);
129 return _p_indices[j];
130 }
131
132 operator std::vector<NodePeriodicBits>() const // copies to a new std::vector<index>
133 {
134 return {_p_indices, _p_indices + _row_size};
135 }
136
137 void operator=(const std::vector<NodePeriodicBits> &r)
138 {
139 DNDS_assert(_row_size == r.size());
140 std::copy(r.begin(), r.end(), _p_indices);
141 }
142
144 {
145 if (this == &r)
146 return;
147 DNDS_assert(_row_size == r.size());
148 std::copy(r.cbegin(), r.cend(), _p_indices);
149 }
150
152 {
154 ret.setP1True();
155 ret.setP2True();
156 ret.setP3True();
157 for (auto &v : *this)
158 ret = ret & v;
159 return ret;
160 }
161
162 NodePeriodicBits *begin() { return _p_indices; }
163 NodePeriodicBits *end() { return _p_indices + _row_size; } // past-end
164 [[nodiscard]] const NodePeriodicBits *cbegin() const { return _p_indices; }
165 [[nodiscard]] const NodePeriodicBits *cend() const { return _p_indices + _row_size; } // past-end
166 [[nodiscard]] rowsize size() const { return _row_size; }
167 };
168
169 template <rowsize _row_size = 1, rowsize _row_max = _row_size, rowsize _align = NoAlign>
170 class ArrayNodePeriodicBits : public ParArray<NodePeriodicBits, _row_size, _row_max, _align>
171 {
172 public:
174 using t_base::t_base;
175
177 {
178 DNDS_assert(i < this->Size()); //! disable past-end input
180 }
181
183 };
184 template <rowsize _row_size = 1, rowsize _row_max = _row_size, rowsize _align = NoAlign>
186
188 {
189 std::array<tGPointPortable, 4> rotation{};
190 std::array<tPointPortable, 4> translation{};
191 std::array<tPointPortable, 4> rotationCenter{};
194 {
195 for (auto &r : rotation)
196 r.map().setIdentity();
197 for (auto &r : rotationCenter)
198 r.map().setZero();
199 translation[0].map().setZero();
200 translation[1].map() = tPoint{1, 0, 0};
201 translation[2].map() = tPoint{0, 1, 0};
202 translation[3].map() = tPoint{0, 0, 1};
203
204 // translation[1] = tPoint{0, 0, 0};
205 // rotation[1] << 0, 1, 0,
206 // -1, 0, 0,
207 // 0, 0, 1;
208 }
209
210 DNDS_HOST void WriteSerializer(Serializer::SerializerBaseSSP serializerP, const std::string &name)
211 {
212 auto cwd = serializerP->GetCurrentPath();
213 serializerP->CreatePath(name);
214 serializerP->GoToPath(name);
215
216 for (int i = 1; i <= 3; i++)
217 {
218 serializerP->WriteRealVector("rotation" + std::to_string(i), Geom::JacobiToSTDVector(rotation.at(i).map()), Serializer::ArrayGlobalOffset_One);
219 serializerP->WriteRealVector("rotationCenter" + std::to_string(i), Geom::VectorToSTDVector(rotationCenter.at(i).map()), Serializer::ArrayGlobalOffset_One);
220 serializerP->WriteRealVector("translation" + std::to_string(i), Geom::VectorToSTDVector(translation.at(i).map()), Serializer::ArrayGlobalOffset_One);
221 }
222
223 serializerP->GoToPath(cwd);
224 }
225
226 DNDS_HOST void ReadSerializer(Serializer::SerializerBaseSSP serializerP, const std::string &name)
227 {
228 auto cwd = serializerP->GetCurrentPath();
229 // serializerP->CreatePath(name); // * no create
230 serializerP->GoToPath(name);
231
232 for (int i = 1; i <= 3; i++)
233 {
234 std::vector<real> rotRead, rotCRead, transRead;
235 Serializer::ArrayGlobalOffset offsetV = Serializer::ArrayGlobalOffset_One; // must explicitly set as One
236 serializerP->ReadRealVector("rotation" + std::to_string(i), rotRead, offsetV);
237 serializerP->ReadRealVector("rotationCenter" + std::to_string(i), rotCRead, offsetV);
238 serializerP->ReadRealVector("translation" + std::to_string(i), transRead, offsetV);
239 rotation.at(i).map() = Geom::STDVectorToJacobi(rotRead);
240 rotationCenter.at(i).map() = Geom::STDVectorToVector(rotCRead);
241 translation.at(i).map() = Geom::STDVectorToVector(transRead);
242 }
243
244 serializerP->GoToPath(cwd);
245 }
246
247 DNDS_DEVICE_CALLABLE [[nodiscard]] tPoint TransCoord(const tPoint &c, t_index id) const
248 {
250 t_index i{0};
251 if (FaceIDIsPeriodicDonor(id))
252 i = -id - 3;
253 else
254 i = -id;
255 return rotation.at(i).map() * (c - rotationCenter.at(i).map()) + rotationCenter.at(i).map() + translation.at(i).map();
256 }
257
258 DNDS_DEVICE_CALLABLE [[nodiscard]] tPoint TransCoordBack(const tPoint &c, t_index id) const
259 {
261 t_index i{0};
262 if (FaceIDIsPeriodicDonor(id))
263 i = -id - 3;
264 else
265 i = -id;
266 return rotation.at(i).map().transpose() * ((c - translation.at(i).map()) - rotationCenter.at(i).map()) + rotationCenter.at(i).map();
267 }
268
269 ///@todo //TODO: add support for cartesian tensor transformation
270
271 template <int dim, int nVec>
272 DNDS_DEVICE_CALLABLE Eigen::Matrix<real, dim, nVec> TransVector(const Eigen::Matrix<real, dim, nVec> &v, t_index id)
273 {
275 t_index i{0};
276 if (FaceIDIsPeriodicDonor(id))
277 i = -id - 3;
278 else
279 i = -id;
280 if constexpr (dim == 3)
281 return rotation.at(i).map() * v;
282 else
283 return rotation.at(i).map()({0, 1}, {0, 1}) * v;
284 }
285
286 template <int dim, int nVec>
287 DNDS_DEVICE_CALLABLE Eigen::Matrix<real, dim, nVec> TransVectorBack(const Eigen::Matrix<real, dim, nVec> &v, t_index id)
288 {
290 t_index i{0};
291 if (FaceIDIsPeriodicDonor(id))
292 i = -id - 3;
293 else
294 i = -id;
295 if constexpr (dim == 3)
296 return rotation.at(i).map().transpose() * v;
297 else
298 return rotation.at(i).map()({0, 1}, {0, 1}).transpose() * v;
299 }
300
301 template <int dim>
302 DNDS_DEVICE_CALLABLE Eigen::Matrix<real, dim, dim> TransMat(const Eigen::Matrix<real, dim, dim> &m, t_index id)
303 {
305 t_index i{0};
306 if (FaceIDIsPeriodicDonor(id))
307 i = -id - 3;
308 else
309 i = -id;
310 if constexpr (dim == 3)
311 return rotation.at(i).map() * m * rotation.at(i).map().transpose();
312 else
313 return rotation.at(i).map()({0, 1}, {0, 1}) * m * rotation.at(i).map()({0, 1}, {0, 1}).transpose();
314 }
315
316 template <int dim>
317 DNDS_DEVICE_CALLABLE Eigen::Matrix<real, dim, dim> TransMatBack(const Eigen::Matrix<real, dim, dim> &m, t_index id)
318 {
320 t_index i{0};
321 if (FaceIDIsPeriodicDonor(id))
322 i = -id - 3;
323 else
324 i = -id;
325 if constexpr (dim == 3)
326 return rotation.at(i).map().transpose() * m * rotation.at(i).map();
327 else
328 return rotation.at(i).map()({0, 1}, {0, 1}).transpose() * m * rotation.at(i).map()({0, 1}, {0, 1});
329 }
330
332 {
333 if (!bool(bits))
334 return c;
335 tPoint ret = c;
336 if (bits.getP3())
337 ret = this->TransCoord(ret, BC_ID_PERIODIC_3);
338 if (bits.getP2())
339 ret = this->TransCoord(ret, BC_ID_PERIODIC_2);
340 if (bits.getP1())
341 ret = this->TransCoord(ret, BC_ID_PERIODIC_1);
342 return ret;
343 }
344
345 template <int dim, int nVec>
346 DNDS_DEVICE_CALLABLE auto GetVectorByBits(const Eigen::Matrix<real, dim, nVec> &v, const NodePeriodicBits &bits)
347 {
348 if (!bool(bits))
349 return v;
350 Eigen::Matrix<real, dim, nVec> ret = v;
351 if (bits.getP3())
352 ret = this->TransVector<dim, nVec>(ret, BC_ID_PERIODIC_3);
353 if (bits.getP2())
354 ret = this->TransVector<dim, nVec>(ret, BC_ID_PERIODIC_2);
355 if (bits.getP1())
356 ret = this->TransVector<dim, nVec>(ret, BC_ID_PERIODIC_1);
357 return ret;
358 }
359
361 {
362 if (!bool(bits))
363 return c;
364 tPoint ret = c;
365 if (bits.getP1())
366 ret = this->TransCoordBack(ret, BC_ID_PERIODIC_1);
367 if (bits.getP2())
368 ret = this->TransCoordBack(ret, BC_ID_PERIODIC_2);
369 if (bits.getP3())
370 ret = this->TransCoordBack(ret, BC_ID_PERIODIC_3);
371 return ret;
372 }
373
374 template <int dim, int nVec>
375 DNDS_DEVICE_CALLABLE auto GetVectorBackByBits(const Eigen::Matrix<real, dim, nVec> &v, const NodePeriodicBits &bits)
376 {
377 if (!bool(bits))
378 return v;
379 Eigen::Matrix<real, dim, nVec> ret = v;
380 if (bits.getP1())
381 ret = this->TransVectorBack<dim, nVec>(ret, BC_ID_PERIODIC_1);
382 if (bits.getP2())
383 ret = this->TransVectorBack<dim, nVec>(ret, BC_ID_PERIODIC_2);
384 if (bits.getP3())
385 ret = this->TransVectorBack<dim, nVec>(ret, BC_ID_PERIODIC_3);
386 return ret;
387 }
388 };
389} // namespace DNDS::Geom
Father-son array pairs with device views and ghost communication.
ParArray (MPI-aware array) and ArrayTransformer (ghost/halo communication).
#define DNDS_DEVICE_TRIVIAL_COPY_DEFINE(T, T_Self)
Definition Defines.hpp:87
#define DNDS_DEVICE_CALLABLE
Definition Defines.hpp:76
#define DNDS_HOST
Definition Defines.hpp:78
#define DNDS_DEVICE_TRIVIAL_COPY_DEFINE_NO_EMPTY_CTOR(T, T_Self)
Definition Defines.hpp:95
#define DNDS_assert(expr)
Debug-only assertion (compiled out when DNDS_NDEBUG is defined). Prints the expression + file/line + ...
Definition Errors.hpp:112
Eigen::Matrix< real, 3, 3 > m
Base types and abstract interface for array serialization.
Non-owning device-callable view of an Array, specialised per DeviceBackend.
T * operator[](index iRow)
Return a raw pointer to the start of row iRow.
Definition Array.hpp:630
rowsize RowSize() const
Uniform row width for fixed layouts (no row index needed).
Definition Array.hpp:176
index Size() const
Number of rows currently stored. O(1).
Definition Array.hpp:171
NodePeriodicBits * rowPtr(index i)
NodePeriodicBitsRow operator[](index i)
NodePeriodicBitsRow(NodePeriodicBitsRow &&)=default
NodePeriodicBitsRow & operator=(NodePeriodicBitsRow &&)=default
void operator=(const NodePeriodicBitsRow &r)
const NodePeriodicBits * cbegin() const
NodePeriodicBitsRow(const NodePeriodicBitsRow &)=default
NodePeriodicBits & operator[](rowsize j)
void operator=(const std::vector< NodePeriodicBits > &r)
NodePeriodicBits operator[](rowsize j) const
const NodePeriodicBits * cend() const
NodePeriodicBitsRow(NodePeriodicBits *ptr, rowsize siz)
MPI-aware Array: adds a communicator, rank, and global index mapping.
Describes one rank's window into a globally-distributed dataset.
Eigen::VectorXd STDVectorToVector(const std::vector< real > &v)
Definition Geometric.hpp:94
tJacobi STDVectorToJacobi(const std::vector< real > &v)
Definition Geometric.hpp:80
int32_t t_index
Definition Geometric.hpp:6
Eigen::Vector3d tPoint
Definition Geometric.hpp:9
DNDS_DEVICE_CALLABLE bool FaceIDIsPeriodic(t_index id)
DNDS_DEVICE_CALLABLE bool FaceIDIsPeriodicDonor(t_index id)
bool isCollaborativeNodePeriodicBits(const std::vector< NodePeriodicBits > &a, const std::vector< NodePeriodicBits > &b)
std::vector< real > JacobiToSTDVector(const tJacobi &j)
Definition Geometric.hpp:72
std::vector< real > VectorToSTDVector(const Eigen::VectorXd &v)
Definition Geometric.hpp:89
ssp< SerializerBase > SerializerBaseSSP
the host side operators are provided as implemented
DNDS_CONSTANT const index UnInitIndex
Sentinel "not initialised" index value (= INT64_MIN).
Definition Defines.hpp:181
int32_t rowsize
Row-width / per-row element-count type (signed 32-bit).
Definition Defines.hpp:114
int64_t index
Global row / DOF index type (signed 64-bit; handles multi-billion-cell meshes).
Definition Defines.hpp:112
Convenience bundle of a father, son, and attached ArrayTransformer.
bool operator!=(const NodeIndexPBI &r) const
bool operator<(const NodeIndexPBI &r) const
bool operator==(const NodeIndexPBI &r) const
DNDS_DEVICE_CALLABLE bool operator==(const NodePeriodicBits &r) const
friend std::ostream & operator<<(std::ostream &o, const NodePeriodicBits &b)
DNDS_DEVICE_CALLABLE NodePeriodicBits operator&(const NodePeriodicBits &r) const
static MPI_Datatype CommType()
DNDS_HOST void WriteSerializer(Serializer::SerializerBaseSSP serializerP, const std::string &name)
DNDS_DEVICE_CALLABLE tPoint GetCoordBackByBits(const tPoint &c, const NodePeriodicBits &bits) const
std::array< tPointPortable, 4 > rotationCenter
DNDS_DEVICE_CALLABLE auto GetVectorByBits(const Eigen::Matrix< real, dim, nVec > &v, const NodePeriodicBits &bits)
DNDS_DEVICE_CALLABLE tPoint TransCoordBack(const tPoint &c, t_index id) const
DNDS_DEVICE_CALLABLE auto GetVectorBackByBits(const Eigen::Matrix< real, dim, nVec > &v, const NodePeriodicBits &bits)
std::array< tPointPortable, 4 > translation
std::array< tGPointPortable, 4 > rotation
DNDS_DEVICE_CALLABLE tPoint GetCoordByBits(const tPoint &c, const NodePeriodicBits &bits) const
DNDS_DEVICE_CALLABLE Eigen::Matrix< real, dim, nVec > TransVectorBack(const Eigen::Matrix< real, dim, nVec > &v, t_index id)
DNDS_DEVICE_CALLABLE Eigen::Matrix< real, dim, dim > TransMat(const Eigen::Matrix< real, dim, dim > &m, t_index id)
DNDS_DEVICE_CALLABLE Eigen::Matrix< real, dim, dim > TransMatBack(const Eigen::Matrix< real, dim, dim > &m, t_index id)
DNDS_HOST void ReadSerializer(Serializer::SerializerBaseSSP serializerP, const std::string &name)
DNDS_DEVICE_CALLABLE Eigen::Matrix< real, dim, nVec > TransVector(const Eigen::Matrix< real, dim, nVec > &v, t_index id)
DNDS_DEVICE_CALLABLE tPoint TransCoord(const tPoint &c, t_index id) const
Eigen::Matrix< real, 5, 1 > v
tVec r(NCells)
tVec b(NCells)
Eigen::Vector3d n(1.0, 0.0, 0.0)