Loading...
Searching...
No Matches
ParaTree.hpp
1/*---------------------------------------------------------------------------*\
2 *
3 * bitpit
4 *
5 * Copyright (C) 2015-2021 OPTIMAD engineering Srl
6 *
7 * -------------------------------------------------------------------------
8 * License
9 * This file is part of bitpit.
10 *
11 * bitpit is free software: you can redistribute it and/or modify it
12 * under the terms of the GNU Lesser General Public License v3 (LGPL)
13 * as published by the Free Software Foundation.
14 *
15 * bitpit is distributed in the hope that it will be useful, but WITHOUT
16 * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
17 * FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
18 * License for more details.
19 *
20 * You should have received a copy of the GNU Lesser General Public License
21 * along with bitpit. If not, see <http://www.gnu.org/licenses/>.
22 *
23 \*---------------------------------------------------------------------------*/
24
25#ifndef __BITPIT_PABLO_PARA_TREE_HPP__
26#define __BITPIT_PABLO_PARA_TREE_HPP__
27
28// =================================================================================== //
29// INCLUDES //
30// =================================================================================== //
31
32// Octant header should be included before bitpit_communications header. If
33// this is not the case, octants will be streamed in the communication buffer
34// using the standard bitpit operator instead of the custom operator defined
35// for the octant. GCC will work also if the Octant header is included after
36// bitpit_communications header, but that's a bug in GCC.
37//
38// https://gcc.gnu.org/bugzilla/show_bug.cgi?id=51577
39// https://bugs.llvm.org/show_bug.cgi?id=47967
40//
41// For a dependent name used in a template definition, the lookup is postponed
42// until the template arguments are known, at which time ADL examines function
43// declarations with external linkage (until C++11) that are visible from the
44// template definition context as well as in the template instantiation context
45// (in other words, ADL will only look in namespaces associated with the
46// argument types), while unqualified lookup only examines function declarations
47// with external linkage (until C++11) that are visible from the template
48// definition context (in other words, adding a new function declaration after
49// template definition does not make it visible except via ADL). The stream
50// operator template is defined in the global namespace, whereas the octant
51// class is defined in the bitpit namespace. Hence, custom stream operator
52// defined for the octant class cannot be selected by ADL (ADL will only look
53// in the namaspace of its arguments, i.e., bitpit, whereas the custom stream
54// operator is defined in the global namespace). The only chance to select
55// the custom operator is through unqualified lookup, and for that to happen,
56// the stream operator defined for the octant class should be defined before
57// the template.
58
59#include "Octant.hpp"
60
61#if BITPIT_ENABLE_MPI==1
62#include <mpi.h>
63#include "DataLBInterface.hpp"
64#include "DataCommInterface.hpp"
65#include "bitpit_communications.hpp"
66#endif
67#include "tree_constants.hpp"
68#include "LocalTree.hpp"
69#include "Map.hpp"
70#include <map>
71#include <unordered_map>
72#include <unordered_set>
73#include <set>
74#include <bitset>
75#include <algorithm>
76
77#include "bitpit_common.hpp"
78
79namespace bitpit {
80
81 // =================================================================================== //
82 // TYPEDEFS //
83 // =================================================================================== //
84 typedef std::vector<bool> bvector;
85 typedef std::vector<int> ivector;
86 typedef std::bitset<72> octantID;
87 typedef std::vector<Octant*> ptroctvector;
88 typedef ptroctvector::iterator octantIterator;
89
90 // =================================================================================== //
91 // CLASS DEFINITION //
92 // =================================================================================== //
93
113 class ParaTree{
114
115 // =================================================================================== //
116 // MEMBERS //
117 // =================================================================================== //
118 public:
119 BITPIT_PUBLIC_API static const std::string DEFAULT_LOG_FILE;
121 typedef std::unordered_map<int, std::array<uint32_t, 2>> ExchangeRanges;
122
123 enum Operation {
124 OP_NONE,
125 OP_INIT,
126 OP_PRE_ADAPT,
127 OP_ADAPT_MAPPED,
128 OP_ADAPT_UNMAPPED,
129 OP_LOADBALANCE_FIRST,
130 OP_LOADBALANCE
131 };
132
135 ACTION_UNDEFINED = -1,
136 ACTION_NONE,
137 ACTION_SEND,
138 ACTION_DELETE,
139 ACTION_RECEIVE
140 };
141
142 ExchangeAction sendAction;
143 ExchangeRanges sendRanges;
144
145 ExchangeAction recvAction;
146 ExchangeRanges recvRanges;
147
149 LoadBalanceRanges(bool serial, const ExchangeRanges &_sendRanges, const ExchangeRanges &_recvRanges);
150
151 void clear();
152 };
153
154 private:
155 typedef std::unordered_map<int, std::array<uint64_t, 2>> PartitionIntersections;
156
157 struct AccretionData {
158 int targetRank;
159 std::unordered_map<uint64_t, int> internalSeeds;
160 std::unordered_map<uint64_t, int> foreignSeeds;
161 std::unordered_map<uint64_t, int> population;
162 };
163
164 //undistributed members
165 std::vector<uint64_t> m_partitionFirstDesc;
166 std::vector<uint64_t> m_partitionLastDesc;
167 std::vector<uint64_t> m_partitionRangeGlobalIdx;
168 std::vector<uint64_t> m_partitionRangeGlobalIdx0;
169 uint64_t m_globalNumOctants;
170 int m_nproc;
171 int8_t m_maxDepth;
172 const TreeConstants *m_treeConstants;
173 std::size_t m_nofGhostLayers;
175 //distributed members
176 int m_rank;
177 LocalTree m_octree;
178 std::map<int,u32vector> m_bordersPerProc;
179 ptroctvector m_internals;
180 ptroctvector m_pborders;
182 //distributed adpapting memebrs
183 u32vector m_mapIdx;
188 //elements sent during last loadbalance operation
189 LoadBalanceRanges m_loadBalanceRanges;
191 //auxiliary members
192 int m_errorFlag;
193 bool m_serial;
194 double m_tol;
196 //map members
197 Map m_trans;
198 uint8_t m_dim;
200 //boundary conditions members
201 bvector m_periodic;
203 //info member
204 uint64_t m_status;
206 Operation m_lastOp;
208 //log member
209 Logger* m_log;
211 //communicator
212#if BITPIT_ENABLE_MPI==1
213 //TODO Duplicate communicator
214 MPI_Comm m_comm;
215#endif
216
217 // =================================================================================== //
218 // CONSTRUCTORS AND OPERATORS //
219 // =================================================================================== //
220 public:
221#if BITPIT_ENABLE_MPI==1
222 ParaTree(const std::string &logfile = DEFAULT_LOG_FILE, MPI_Comm comm = MPI_COMM_WORLD);
223 ParaTree(uint8_t dim, const std::string &logfile = DEFAULT_LOG_FILE, MPI_Comm comm = MPI_COMM_WORLD);
224 ParaTree(std::istream &stream, const std::string &logfile = DEFAULT_LOG_FILE, MPI_Comm comm = MPI_COMM_WORLD);
225#else
226 ParaTree(const std::string &logfile = DEFAULT_LOG_FILE);
227 ParaTree(uint8_t dim, const std::string &logfile = DEFAULT_LOG_FILE);
228 ParaTree(std::istream &stream, const std::string &logfile = DEFAULT_LOG_FILE);
229#endif
230 virtual ~ParaTree();
231
232 ParaTree(const ParaTree & other);
233
234 private:
235 ParaTree & operator=(const ParaTree & other) = delete;
236
237 // =================================================================================== //
238 // METHODS //
239 // =================================================================================== //
240 public:
241 virtual void reset();
242
243 virtual int getDumpVersion() const;
244 virtual void dump(std::ostream &stream, bool full = true);
245 virtual void restore(std::istream &stream);
246
247 void printHeader();
248
249 // =================================================================================== //
250 // BASIC GET/SET METHODS //
251 // =================================================================================== //
252 uint8_t getDim() const;
253 uint64_t getGlobalNumOctants() const;
254 bool getSerial() const;
255 bool getParallel() const;
256 int getRank() const;
257 int getNproc() const;
258 Logger& getLog();
259 Operation getLastOperation() const;
260#if BITPIT_ENABLE_MPI==1
261 void setComm(MPI_Comm communicator);
262 void replaceComm(MPI_Comm communicator, MPI_Comm *previousCommunicator = nullptr);
263 void freeComm();
264 MPI_Comm getComm() const;
265 bool isCommSet() const;
266#endif
267 const std::vector<uint64_t> &getPartitionRangeGlobalIdx() const;
268 const std::vector<uint64_t> &getPartitionFirstDesc() const;
269 const std::vector<uint64_t> &getPartitionLastDesc() const;
270 darray3 getOrigin() const;
271 double getX0() const;
272 double getY0() const;
273 double getZ0() const;
274 double getL() const;
275 int getMaxLevel() const;
276 uint32_t getMaxLength() const;
277 uint8_t getNnodes() const;
278 uint8_t getNfaces() const;
279 uint8_t getNedges() const;
280 uint8_t getNchildren() const;
281 uint8_t getNnodesperface() const;
282 void getNormals(int8_t normals[6][3]) const;
283 void getOppface(uint8_t oppface[6]) const;
284 void getFacenode(uint8_t facenode[6][4]) const;
285 void getNodeface(uint8_t nodeface[8][3]) const;
286 void getEdgeface(uint8_t edgeface[12][2]) const;
287 void getNodecoeffs(int8_t nodecoeffs[8][3]) const;
288 void getEdgecoeffs(int8_t edgecoeffs[12][3]) const;
289 const int8_t (*getNormals() const)[3];
290 const uint8_t *getOppface() const;
291 const uint8_t (*getFacenode() const)[4];
292 const uint8_t (*getNodeface() const)[3];
293 const uint8_t (*getEdgeface() const)[2];
294 const int8_t (*getNodecoeffs() const)[3];
295 const int8_t (*getEdgecoeffs() const)[3];
296 bvector getPeriodic() const;
297 double getTol() const;
298 bool getPeriodic(uint8_t i) const;
299 void setPeriodic(uint8_t i);
300 void setTol(double tol = 1.0e-14);
301
302 // =================================================================================== //
303 // INDEX BASED METHODS //
304 // =================================================================================== //
305 darray3 getCoordinates(uint32_t idx) const;
306 double getX(uint32_t idx) const;
307 double getY(uint32_t idx) const;
308 double getZ(uint32_t idx) const;
309 double getSize(uint32_t idx) const;
310 double getArea(uint32_t idx) const;
311 double getVolume(uint32_t idx) const;
312 void getCenter(uint32_t idx, darray3& centerCoords) const;
313 darray3 getCenter(uint32_t idx) const;
314 darray3 getFaceCenter(uint32_t idx, uint8_t face) const;
315 void getFaceCenter(uint32_t idx, uint8_t face, darray3& centerCoords) const;
316 darray3 getNode(uint32_t idx, uint8_t node) const;
317 void getNode(uint32_t idx, uint8_t node, darray3& nodeCoords) const;
318 void getNodes(uint32_t idx, darr3vector & nodesCoords) const;
319 darr3vector getNodes(uint32_t idx) const;
320 void getNormal(uint32_t idx, uint8_t face, darray3 & normal) const;
321 darray3 getNormal(uint32_t idx, uint8_t face) const;
322 int8_t getMarker(uint32_t idx) const;
323 uint8_t getLevel(uint32_t idx) const;
324 uint64_t getMorton(uint32_t idx) const;
325 uint64_t computeNodePersistentKey(uint32_t idx, uint8_t node) const;
326 bool getBalance(uint32_t idx) const;
327 bool getBound(uint32_t idx, uint8_t face) const;
328 bool getBound(uint32_t idx) const;
329 bool getPbound(uint32_t idx, uint8_t face) const;
330 bool getPbound(uint32_t idx) const;
331 bool getIsNewR(uint32_t idx) const;
332 bool getIsNewC(uint32_t idx) const;
333 uint64_t getGlobalIdx(uint32_t idx) const;
334 uint64_t getGhostGlobalIdx(uint32_t idx) const;
335 uint32_t getLocalIdx(uint64_t gidx) const;
336 uint32_t getLocalIdx(uint64_t gidx,int rank) const;
337 void getLocalIdx(uint64_t gidx,uint32_t & lidx,int & rank) const;
338 uint32_t getGhostLocalIdx(uint64_t gidx) const;
339 octantID getPersistentIdx(uint32_t idx) const;
340 int8_t getPreMarker(uint32_t idx);
341 void setMarker(uint32_t idx, int8_t marker);
342 void setBalance(uint32_t idx, bool balance);
343
344 // =================================================================================== //
345 // POINTER BASED METHODS //
346 // =================================================================================== //
347 darray3 getCoordinates(const Octant* oct) const;
348 double getX(const Octant* oct) const;
349 double getY(const Octant* oct) const;
350 double getZ(const Octant* oct) const;
351 double getSize(const Octant* oct) const;
352 double getArea(const Octant* oct) const;
353 double getVolume(const Octant* oct) const;
354 void getCenter(const Octant* oct, darray3& centerCoords) const;
355 darray3 getCenter(const Octant* oct) const;
356 darray3 getFaceCenter(const Octant* oct, uint8_t face) const;
357 void getFaceCenter(const Octant* oct, uint8_t face, darray3& centerCoords) const;
358 darray3 getNode(const Octant* oct, uint8_t node) const;
359 void getNode(const Octant* oct, uint8_t node, darray3& nodeCoords) const;
360 void getNodes(const Octant* oct, darr3vector & nodesCoords) const;
361 darr3vector getNodes(const Octant* oct) const;
362 void getNormal(const Octant* oct, uint8_t face, darray3 & normal) const;
363 darray3 getNormal(const Octant* oct, uint8_t face) const;
364 int8_t getMarker(const Octant* oct) const;
365 uint8_t getLevel(const Octant* oct) const;
366 uint64_t getMorton(const Octant* oct) const;
367 uint64_t getLastDescMorton(const Octant* oct) const;
368 uint64_t computeNodePersistentKey(const Octant* oct, uint8_t node) const;
369 bool getBalance(const Octant* oct) const;
370 bool getBound(const Octant* oct, uint8_t face) const;
371 bool getBound(const Octant* oct) const;
372 bool getPbound(const Octant* oct, uint8_t face) const;
373 bool getPbound(const Octant* oct) const;
374 bool getIsNewR(const Octant* oct) const;
375 bool getIsNewC(const Octant* oct) const;
376 uint32_t getIdx(const Octant* oct) const;
377 uint64_t getGlobalIdx(const Octant* oct) const;
378 octantID getPersistentIdx(const Octant* oct) const;
379 int8_t getPreMarker(Octant* oct);
380 void setMarker(Octant* oct, int8_t marker);
381 void setBalance(Octant* oct, bool balance);
382
383 // =================================================================================== //
384 // LOCAL TREE GET/SET METHODS //
385 // =================================================================================== //
386 uint64_t getStatus() const;
387 uint32_t getNumOctants() const;
388 uint32_t getNumGhosts() const;
389 uint32_t getNumNodes() const;
390 uint8_t getLocalMaxDepth() const;
391 double getLocalMaxSize() const;
392 double getLocalMinSize() const;
393 uint8_t getBalanceCodimension() const;
394 uint64_t getFirstDescMorton() const;
395 uint64_t getLastDescMorton() const;
396 uint64_t getLastDescMorton(uint32_t idx) const;
397 octantIterator getInternalOctantsBegin();
398 octantIterator getInternalOctantsEnd();
399 octantIterator getPboundOctantsBegin();
400 octantIterator getPboundOctantsEnd();
401 void setBalanceCodimension(uint8_t b21codim);
402
403 // =================================================================================== //
404 // INTERSECTION GET/SET METHODS //
405 // =================================================================================== //
406 uint32_t getNumIntersections() const;
407 Intersection* getIntersection(uint32_t idx);
408 uint8_t getLevel(const Intersection* inter) const;
409 bool getFiner(const Intersection* inter) const;
410 bool getBound(const Intersection* inter) const;
411 bool getIsGhost(const Intersection* inter) const;
412 bool getPbound(const Intersection* inter) const;
413 uint8_t getFace(const Intersection* inter) const;
414 u32vector getOwners(const Intersection* inter) const;
415 uint32_t getIn(const Intersection* inter) const;
416 uint32_t getOut(const Intersection* inter) const;
417 bool getOutIsGhost(const Intersection* inter) const;
418 double getSize(const Intersection* inter) const;
419 double getArea(const Intersection* inter) const;
420 darray3 getCenter(const Intersection* inter) const;
421 darr3vector getNodes(const Intersection* inter) const;
422 darray3 getNormal(const Intersection* inter) const;
423
424 // =================================================================================== //
425 // OTHER GET/SET METHODS //
426 // =================================================================================== //
427 Octant* getOctant(uint32_t idx);
428 const Octant* getOctant(uint32_t idx) const;
429 Octant* getGhostOctant(uint32_t idx);
430 const Octant* getGhostOctant(uint32_t idx) const;
431 bool getIsGhost(const Octant* oct) const;
432 int getGhostLayer(const Octant* oct) const;
433 const LoadBalanceRanges & getLoadBalanceRanges() const;
434 std::size_t getNofGhostLayers() const;
435 void setNofGhostLayers(std::size_t nofGhostLayers);
436 const std::map<int, std::vector<uint32_t>> & getBordersPerProc() const;
437
438 // =================================================================================== //
439 // PRIVATE GET/SET METHODS //
440 // =================================================================================== //
441 private:
442 void setDim(uint8_t dim);
443#if BITPIT_ENABLE_MPI==1
444 void updateGlobalFirstDescMorton();
445 void updateGlobalLasttDescMorton();
446#endif
447
448 // =================================================================================== //
449 // OTHER METHODS //
450 // =================================================================================== //
451
452#if BITPIT_ENABLE_MPI==1
453 void _initializeCommunicator(MPI_Comm comm);
454#else
455 void _initializeSerialCommunicator();
456#endif
457 void _initializePartitions();
458 void _initialize(uint8_t dim, const std::string &logfile);
459
460#if BITPIT_ENABLE_MPI==1
461 void initialize(const std::string &logfile, MPI_Comm comm);
462 void initialize(uint8_t dim, const std::string &logfile, MPI_Comm comm);
463#else
464 void initialize(const std::string &logfile);
465 void initialize(uint8_t dim, const std::string &logfile);
466#endif
467 void initializeLogger(const std::string &logfile);
468 void reinitialize(uint8_t dim, const std::string &logfile);
469
470 void reset(bool createRoot);
471
472 uint64_t evalPointAnchorMorton(const double * point) const;
473
474 // =================================================================================== //
475 // OTHER OCTANT BASED METHODS //
476 // =================================================================================== //
477
478 void findNeighbours(const Octant* oct, uint8_t face, uint8_t codim, u32vector & neighbours, bvector & isghost, bool onlyinternals, bool append) const;
479 public:
480 void findNeighbours(uint32_t idx, uint8_t face, uint8_t codim, u32vector & neighbours, bvector & isghost) const;
481 void findNeighbours(uint32_t idx, uint8_t face, uint8_t codim, u32vector & neighbours) const;
482 void findNeighbours(const Octant* oct, uint8_t face, uint8_t codim, u32vector & neighbours, bvector & isghost) const ;
483 void findGhostNeighbours(uint32_t idx, uint8_t face, uint8_t codim, u32vector & neighbours) const;
484 void findGhostNeighbours(uint32_t idx, uint8_t face, uint8_t codim, u32vector & neighbours, bvector & isghost) const;
485 void findGhostNeighbours(const Octant* oct, uint8_t face, uint8_t codim, u32vector & neighbours, bvector & isghost) const;
486 void findAllNodeNeighbours(uint32_t idx, uint32_t node, u32vector & neighbours, bvector & isghost);
487 void findAllNodeNeighbours(const Octant* oct, uint32_t node, u32vector & neighbours, bvector & isghost) const;
488 void findAllCodimensionNeighbours(uint32_t idx, u32vector & neighbours, bvector & isghost);
489 void findAllCodimensionNeighbours(const Octant* oct, u32vector & neighbours, bvector & isghost);
490 void findGhostAllCodimensionNeighbours(uint32_t idx, u32vector & neighbours, bvector & isghost);
491 void findGhostAllCodimensionNeighbours(Octant* oct, u32vector & neighbours, bvector & isghost);
492 Octant* getPointOwner(const dvector &point);
493 Octant* getPointOwner(const dvector &point, bool & isghost);
494 Octant* getPointOwner(const darray3 &point);
495 Octant* getPointOwner(const darray3 &point, bool & isghost);
496 uint32_t getPointOwnerIdx(const double * point) const;
497 uint32_t getPointOwnerIdx(const double * point, bool & isghost) const;
498 uint32_t getPointOwnerIdx(const dvector &point) const;
499 uint32_t getPointOwnerIdx(const dvector &point, bool & isghost) const;
500 uint32_t getPointOwnerIdx(const darray3 &point) const;
501 uint32_t getPointOwnerIdx(const darray3 &point, bool & isghost) const;
502 void getMapping(uint32_t idx, u32vector & mapper, bvector & isghost) const;
503 void getMapping(uint32_t idx, u32vector & mapper, bvector & isghost, ivector & rank) const;
504 void getPreMapping(u32vector & idx, std::vector<int8_t> & markers, std::vector<bool> & isghost);
505 bool isNodeOnOctant(const Octant* nodeOctant, uint8_t nodeIndex, const Octant* octant) const;
506 bool isEdgeOnOctant(const Octant* edgeOctant, uint8_t edgeIndex, const Octant* octant) const;
507 bool isFaceOnOctant(const Octant* faceOctant, uint8_t faceIndex, const Octant* octant) const;
508 int getPointOwnerRank(const darray3 &point);
509 uint8_t getFamilySplittingNode(const Octant*) const;
510 void expectedOctantAdapt(const Octant* octant, int8_t marker, octvector* result) const;
511
512 // =================================================================================== //
513 // OTHER PARATREE BASED METHODS //
514 // =================================================================================== //
515 int8_t getMaxDepth() const;
516 int findOwner(uint64_t morton) const;
517 int getOwnerRank(uint64_t globalIdx) const;
518 void settleMarkers();
519 void preadapt();
520 bool checkToAdapt();
521 bool adapt(bool mapper_flag = false);
522 bool adaptGlobalRefine(bool mapper_flag = false);
523 bool adaptGlobalCoarse(bool mapper_flag = false);
524 void computeConnectivity();
525 void clearConnectivity(bool release = true);
526 void updateConnectivity();
527 const u32vector2D & getConnectivity() const;
528 const u32vector & getConnectivity(uint32_t idx) const;
529 const u32vector & getConnectivity(Octant* oct) const;
530 const u32arr3vector & getNodes() const;
531 const u32array3 & getNodeLogicalCoordinates(uint32_t node) const;
532 darray3 getNodeCoordinates(uint32_t node) const;
533 const u32vector2D & getGhostConnectivity() const;
534 const u32vector & getGhostConnectivity(uint32_t idx) const;
535 const u32vector & getGhostConnectivity(const Octant* oct) const;
536 bool check21Balance();
537#if BITPIT_ENABLE_MPI==1
538 void loadBalance(const dvector* weight = NULL);
539 void loadBalance(uint8_t & level, const dvector* weight = NULL);
540
541 LoadBalanceRanges evalLoadBalanceRanges(dvector *weights);
542 LoadBalanceRanges evalLoadBalanceRanges(uint8_t level, dvector *weights);
543 private:
544 LoadBalanceRanges evalLoadBalanceRanges(const uint32_t *updatedPartition);
545
546 ExchangeRanges evalLoadBalanceSendRanges(const uint32_t *updatedPartition);
547 ExchangeRanges evalLoadBalanceRecvRanges(const uint32_t *updatedPartition);
548
549 PartitionIntersections evalPartitionIntersections(const uint32_t *schema_A, int rank_A, const uint32_t *schema_B);
550#endif
551 public:
552 double levelToSize(uint8_t level) const;
553
554 // =================================================================================== //
555 // OTHER INTERSECTION BASED METHODS //
556 // =================================================================================== //
558
559 // =================================================================================== //
560 // OTHER PRIVATE METHODS //
561 // =================================================================================== //
562 private:
563 Octant& extractOctant(uint32_t idx);
564 bool private_adapt_mapidx(bool mapflag);
565 void updateAdapt();
566#if BITPIT_ENABLE_MPI==1
567 void computePartition(uint32_t *partition);
568 void computePartition(const dvector *weight, uint32_t *partition);
569 void computePartition(uint8_t level_, const dvector *weight, uint32_t *partition);
570 void updateLoadBalance();
571 void setPboundGhosts();
572 void buildGhostOctants(const std::map<int, u32vector> &bordersPerProc, const std::vector<AccretionData> &accretions);
573 void initializeGhostHaloAccretions(std::vector<AccretionData> *accretions);
574 void growGhostHaloAccretions(std::unordered_map<uint32_t, std::vector<uint64_t>> *cachedOneRings, std::vector<AccretionData> *accretions);
575 void exchangeGhostHaloAccretions(DataCommunicator *dataCommunicator, std::vector<AccretionData> *accretions);
576
577 void computeGhostHalo();
578 bool commMarker();
579#endif
580 void updateAfterCoarse();
581 void balance21(bool verbose, bool balanceNewOctants);
582 void createPartitionInfo();
583
584 bool isInternal(uint64_t gidx) const;
585
586 // =================================================================================== //
587 // TESTING OUTPUT METHODS //
588 // =================================================================================== //
589 public:
590 void write(const std::string &filename);
591 void writeTest(const std::string &filename, dvector data);
592
593 // =================================================================================== //
594 // TEMPLATE METHODS //
595 // =================================================================================== //
596#if BITPIT_ENABLE_MPI==1
597
600 template<class Impl>
601 void
603 //BUILD SEND BUFFERS
604 DataCommunicator communicator(m_comm);
605 size_t fixedDataSize = userData.fixedSize();
606 std::map<int,u32vector >::iterator bitend = m_bordersPerProc.end();
607 std::map<int,u32vector >::iterator bitbegin = m_bordersPerProc.begin();
608 for(std::map<int,u32vector >::iterator bit = bitbegin; bit != bitend; ++bit){
609 int key = bit->first;
610 const u32vector & pborders = bit->second;
611 size_t buffSize = 0;
612 size_t nofPbordersPerProc = pborders.size();
613 if(fixedDataSize != 0){
614 buffSize = fixedDataSize*nofPbordersPerProc;
615 }
616 else{
617 for(size_t i = 0; i < nofPbordersPerProc; ++i){
618 buffSize += userData.size(pborders[i]);
619 }
620 }
621 //enlarge buffer to store number of pborders from this proc
622 buffSize += sizeof(size_t);
623 //build buffer for this proc
624 communicator.setSend(key,buffSize);
625 SendBuffer & sendBuffer = communicator.getSendBuffer(key);
626 //store number of pborders from this proc at the begining
627 sendBuffer << nofPbordersPerProc;
628
629 //WRITE SEND BUFFERS
630 for(size_t j = 0; j < nofPbordersPerProc; ++j){
631 userData.gather(sendBuffer,pborders[j]);
632 }
633 }
634
635 communicator.discoverRecvs();
636 communicator.startAllRecvs();
637
638 communicator.startAllSends();
639
640 //READ RECEIVE BUFFERS
641 int ghostOffset = 0;
642 std::vector<int> recvRanks = communicator.getRecvRanks();
643 std::sort(recvRanks.begin(),recvRanks.end());
644 for(int rank : recvRanks){
645 communicator.waitRecv(rank);
646 RecvBuffer & recvBuffer = communicator.getRecvBuffer(rank);
647 size_t nofGhostFromThisProc = 0;
648 recvBuffer >> nofGhostFromThisProc;
649 for(size_t k = 0; k < nofGhostFromThisProc; ++k){
650 userData.scatter(recvBuffer, k+ghostOffset);
651 }
652 ghostOffset += nofGhostFromThisProc;
653 }
654 communicator.waitAllSends();
655 }
656
664 template<class Impl>
665 void
666 loadBalance(DataLBInterface<Impl> & userData, dvector* weight = NULL){
667 //Write info on log
668 (*m_log) << "---------------------------------------------" << std::endl;
669 (*m_log) << " LOAD BALANCE " << std::endl;
670
671 m_lastOp = OP_LOADBALANCE;
672 if (m_nproc>1){
673
674 std::vector<uint32_t> partition(m_nproc);
675 if (weight == NULL)
676 computePartition(partition.data());
677 else
678 computePartition(weight, partition.data());
679
680 weight = NULL;
681
682 privateLoadBalance(partition.data(), &userData);
683
684 //Write info of final partition on log
685 (*m_log) << " " << std::endl;
686 (*m_log) << " Final Parallel partition : " << std::endl;
687 (*m_log) << " Octants for proc "+ std::to_string(static_cast<unsigned long long>(0))+" : " + std::to_string(static_cast<unsigned long long>(m_partitionRangeGlobalIdx[0]+1)) << std::endl;
688 for(int ii=1; ii<m_nproc; ii++){
689 (*m_log) << " Octants for proc "+ std::to_string(static_cast<unsigned long long>(ii))+" : " + std::to_string(static_cast<unsigned long long>(m_partitionRangeGlobalIdx[ii]-m_partitionRangeGlobalIdx[ii-1])) << std::endl;
690 }
691 (*m_log) << " " << std::endl;
692 (*m_log) << "---------------------------------------------" << std::endl;
693
694 }
695 else{
696 m_loadBalanceRanges.clear();
697
698 (*m_log) << " " << std::endl;
699 (*m_log) << " Serial partition : " << std::endl;
700 (*m_log) << " Octants for proc "+ std::to_string(static_cast<unsigned long long>(0))+" : " + std::to_string(static_cast<unsigned long long>(m_partitionRangeGlobalIdx[0]+1)) << std::endl;
701 (*m_log) << " " << std::endl;
702 (*m_log) << "---------------------------------------------" << std::endl;
703 }
704
705 }
706
715 template<class Impl>
716 void
717 loadBalance(DataLBInterface<Impl> & userData, uint8_t & level, dvector* weight = NULL){
718
719 //Write info on log
720 (*m_log) << "---------------------------------------------" << std::endl;
721 (*m_log) << " LOAD BALANCE " << std::endl;
722
723 m_lastOp = OP_LOADBALANCE;
724 if (m_nproc>1){
725
726 std::vector<uint32_t> partition(m_nproc);
727 computePartition(level, weight, partition.data());
728
729 privateLoadBalance(partition.data(), &userData);
730
731 //Write info of final partition on log
732 (*m_log) << " " << std::endl;
733 (*m_log) << " Final Parallel partition : " << std::endl;
734 (*m_log) << " Octants for proc "+ std::to_string(static_cast<unsigned long long>(0))+" : " + std::to_string(static_cast<unsigned long long>(m_partitionRangeGlobalIdx[0]+1)) << std::endl;
735 for(int ii=1; ii<m_nproc; ii++){
736 (*m_log) << " Octants for proc "+ std::to_string(static_cast<unsigned long long>(ii))+" : " + std::to_string(static_cast<unsigned long long>(m_partitionRangeGlobalIdx[ii]-m_partitionRangeGlobalIdx[ii-1])) << std::endl;
737 }
738 (*m_log) << " " << std::endl;
739 (*m_log) << "---------------------------------------------" << std::endl;
740
741 }
742 else{
743 m_loadBalanceRanges.clear();
744
745 (*m_log) << " " << std::endl;
746 (*m_log) << " Serial partition : " << std::endl;
747 (*m_log) << " Octants for proc "+ std::to_string(static_cast<unsigned long long>(0))+" : " + std::to_string(static_cast<unsigned long long>(m_partitionRangeGlobalIdx[0]+1)) << std::endl;
748 (*m_log) << " " << std::endl;
749 (*m_log) << "---------------------------------------------" << std::endl;
750 }
751
752 }
753
763 template<class Impl>
764 void
765 privateLoadBalance(const uint32_t *partition, DataLBInterface<Impl> *userData = nullptr){
766
767 (*m_log) << " " << std::endl;
768 if (m_serial) {
769 (*m_log) << " Initial Serial distribution : " << std::endl;
770 } else {
771 (*m_log) << " Initial Parallel partition : " << std::endl;
772 }
773 (*m_log) << " Octants for proc "+ std::to_string(static_cast<unsigned long long>(0))+" : " + std::to_string(static_cast<unsigned long long>(m_partitionRangeGlobalIdx[0]+1)) << std::endl;
774 for(int ii=1; ii<m_nproc; ii++){
775 (*m_log) << " Octants for proc "+ std::to_string(static_cast<unsigned long long>(ii))+" : " + std::to_string(static_cast<unsigned long long>(m_partitionRangeGlobalIdx[ii]-m_partitionRangeGlobalIdx[ii-1])) << std::endl;
776 }
777
778 // Compute load balance ranges
779 std::unordered_map<int, std::array<uint32_t, 2>> sendRanges = evalLoadBalanceSendRanges(partition);
780 std::unordered_map<int, std::array<uint32_t, 2>> recvRanges = evalLoadBalanceRecvRanges(partition);
781
782 m_loadBalanceRanges = LoadBalanceRanges(m_serial, sendRanges, recvRanges);
783
784 // Compute information about the new partitioning
785 assert(m_nproc > 0);
786
787 uint32_t newSizeOctants = partition[m_rank];
788
789 uint64_t newLastOctantGlobalIdx;
790 uint64_t newFirstOctantGlobalIdx;
791 if (newSizeOctants > 0) {
792 newLastOctantGlobalIdx = partition[0];
793 for (int p = 1; p < m_rank + 1; ++p) {
794 newLastOctantGlobalIdx += partition[p];
795 }
796 --newLastOctantGlobalIdx;
797
798 newFirstOctantGlobalIdx = newLastOctantGlobalIdx - newSizeOctants + 1;
799 } else {
800 newLastOctantGlobalIdx = 0;
801 newFirstOctantGlobalIdx = 1;
802 }
803
804 // Partition internal octants
805 if (m_serial) {
806 m_lastOp = OP_LOADBALANCE_FIRST;
807
808 if (newFirstOctantGlobalIdx != 0) {
809 for (uint32_t i = 0; i < newSizeOctants; ++i){
810 m_octree.m_octants[i] = m_octree.m_octants[newFirstOctantGlobalIdx + i];
811 if (userData) {
812 userData->move(newFirstOctantGlobalIdx + i, i);
813 }
814 }
815 }
816
817 m_octree.m_octants.resize(newSizeOctants);
818 m_octree.m_octants.shrink_to_fit();
819
820 if (userData) {
821 userData->resize(newSizeOctants);
822 userData->shrink();
823 }
824 } else {
825 // Compute information about the current partitioning
826 uint64_t lastOctantGlobalIdx = m_partitionRangeGlobalIdx[m_rank];
827 uint64_t firstOctantGlobalIdx = lastOctantGlobalIdx - m_octree.m_octants.size() + 1;
828
829 // Initialize communications
830 DataCommunicator lbCommunicator(m_comm);
831
832 if (!userData || userData->fixedSize()) {
833 for (const auto &entry : recvRanges) {
834 int rank = entry.first;
835 uint32_t beginRecvIdx = entry.second[0];
836 uint32_t endRecvIdx = entry.second[1];
837
838 uint32_t nOctantsToReceive = endRecvIdx - beginRecvIdx;
839 std::size_t buffSize = nOctantsToReceive * Octant::getBinarySize();
840 if (userData) {
841 buffSize += nOctantsToReceive * userData->fixedSize();
842 }
843
844 lbCommunicator.setRecv(rank, buffSize);
845 lbCommunicator.startRecv(rank);
846 }
847 }
848
849 for (const auto &entry : sendRanges) {
850 int rank = entry.first;
851 uint32_t beginSendIdx = entry.second[0];
852 uint32_t endSendIdx = entry.second[1];
853
854 uint32_t nOctantsToSend = endSendIdx - beginSendIdx;
855 std::size_t buffSize = nOctantsToSend * Octant::getBinarySize();
856 if (userData) {
857 if (userData->fixedSize()) {
858 buffSize += nOctantsToSend * userData->fixedSize();
859 } else {
860 for (uint32_t i = beginSendIdx; i < endSendIdx; ++i) {
861 buffSize += userData->size(i);
862 }
863 }
864 }
865 lbCommunicator.setSend(rank, buffSize);
866
867 SendBuffer &sendBuffer = lbCommunicator.getSendBuffer(rank);
868 for (uint32_t i = beginSendIdx; i < endSendIdx; ++i) {
869 sendBuffer << m_octree.m_octants[i];
870 userData->gather(sendBuffer, i);
871 }
872
873 lbCommunicator.startSend(rank);
874 }
875
876 if (userData && !userData->fixedSize()) {
877 lbCommunicator.discoverRecvs();
878 lbCommunicator.startAllRecvs();
879 }
880
881 // Move resident octants into place
882 if (newSizeOctants > m_octree.m_octants.size()) {
883 m_octree.m_octants.reserve(newSizeOctants);
884 m_octree.m_octants.resize(newSizeOctants);
885
886 if (userData) {
887 userData->resize(newSizeOctants);
888 }
889 }
890
891 bool hasResidentOctants;
892 if (newSizeOctants > 0) {
893 hasResidentOctants = true;
894 if (newFirstOctantGlobalIdx > lastOctantGlobalIdx) {
895 hasResidentOctants = false;
896 } else if (newLastOctantGlobalIdx < firstOctantGlobalIdx) {
897 hasResidentOctants = false;
898 }
899 } else {
900 hasResidentOctants = false;
901 }
902
903 if (hasResidentOctants) {
904 uint64_t firstResidentGlobalOctantIdx = std::max(firstOctantGlobalIdx, newFirstOctantGlobalIdx);
905 uint64_t lastResidentGlobalOctantIdx = std::min(lastOctantGlobalIdx, newLastOctantGlobalIdx);
906
907 uint32_t nofResidents = static_cast<uint32_t>(lastResidentGlobalOctantIdx - firstResidentGlobalOctantIdx + 1);
908 uint32_t newFirstResidentOffsetIdx = static_cast<uint32_t>(firstResidentGlobalOctantIdx - newFirstOctantGlobalIdx);
909 uint32_t firstResidentOffsetIdx = static_cast<uint32_t>(firstResidentGlobalOctantIdx - firstOctantGlobalIdx);
910
911 // If residents are moved closer to the head, we need to move
912 // them from the first to the last. Otherwise, if resident are
913 // moved closer to the tail, we need to move them in reversed
914 // order, i.e., from the last to the first (otherwise octants
915 // will be overwritten during the relocaiton).
916 if (newFirstResidentOffsetIdx != firstResidentOffsetIdx) {
917 for (uint32_t i = 0; i < nofResidents; ++i) {
918 int residentIdx;
919 if (newFirstResidentOffsetIdx < firstResidentOffsetIdx) {
920 residentIdx = i;
921 } else {
922 residentIdx = nofResidents - i - 1;
923 }
924
925 m_octree.m_octants[newFirstResidentOffsetIdx + residentIdx] = m_octree.m_octants[firstResidentOffsetIdx + residentIdx];
926 if (userData) {
927 userData->move(firstResidentOffsetIdx + residentIdx, newFirstResidentOffsetIdx + residentIdx);
928 }
929 }
930 }
931 }
932
933 if (newSizeOctants < m_octree.m_octants.size()) {
934 m_octree.m_octants.resize(newSizeOctants);
935 m_octree.m_octants.shrink_to_fit();
936
937 if (userData) {
938 userData->resize(newSizeOctants);
939 userData->shrink();
940 }
941 }
942
943 // Read buffers and build new octants
944 int nCompletedRecvs = 0;
945 while (nCompletedRecvs < lbCommunicator.getRecvCount()) {
946 int senderRank = lbCommunicator.waitAnyRecv();
947 RecvBuffer &recvBuffer = lbCommunicator.getRecvBuffer(senderRank);
948
949 const std::array<uint32_t, 2> &recvRange = recvRanges.at(senderRank);
950 uint32_t beginRecvIdx = recvRange[0];
951 uint32_t endRecvIdx = recvRange[1];
952 assert(userData || ((endRecvIdx - beginRecvIdx) == (recvBuffer.getSize() / (Octant::getBinarySize()))));
953 assert(!userData || !userData->fixedSize() || ((endRecvIdx - beginRecvIdx) == (recvBuffer.getSize() / (Octant::getBinarySize() + userData->fixedSize()))));
954
955 for (uint32_t i = beginRecvIdx; i < endRecvIdx; ++i) {
956 recvBuffer >> m_octree.m_octants[i];
957 if (userData) {
958 userData->scatter(recvBuffer, i);
959 }
960 }
961
962 ++nCompletedRecvs;
963 }
964
965 lbCommunicator.waitAllSends();
966 }
967
968 // Update load balance information
969 updateLoadBalance();
970
971 // Update ghosts
972 m_octree.m_ghosts.clear();
973
974 computeGhostHalo();
975
976 if (userData) {
977 userData->resizeGhost(m_octree.getNumGhosts());
978 }
979 }
980#endif
981
982 // =============================================================================== //
983
984 };
985
986}
987
988#endif /* __BITPIT_PARA_TREE_HPP__ */
Base class for data communications.
void scatter(Buffer &buff, const uint32_t e)
void gather(Buffer &buff, const uint32_t e)
size_t size(const uint32_t e) const
The DataCommunicator class provides the infrastructure needed to exchange data among processes.
void setRecv(int rank, long length=0)
void startRecv(int srcRank)
SendBuffer & getSendBuffer(int rank)
const std::vector< int > & getRecvRanks() const
int waitAnyRecv(const std::vector< int > &blackList=std::vector< int >())
void setSend(int rank, long length=0)
RecvBuffer & getRecvBuffer(int rank)
void startSend(int dstRank)
Base class for data communications.
static unsigned int getBinarySize()
Definition Octant.cpp:658
Para Tree is the user interface class.
Definition ParaTree.hpp:113
uint8_t getNnodesperface() const
Octant * getPointOwner(const dvector &point)
void write(const std::string &filename)
uint32_t getLocalIdx(uint64_t gidx) const
double levelToSize(uint8_t level) const
double getTol() const
uint32_t getGhostLocalIdx(uint64_t gidx) const
uint8_t getDim() const
Definition ParaTree.cpp:780
void setMarker(uint32_t idx, int8_t marker)
bool isFaceOnOctant(const Octant *faceOctant, uint8_t faceIndex, const Octant *octant) const
const u32array3 & getNodeLogicalCoordinates(uint32_t node) const
uint32_t getNumIntersections() const
int8_t getMaxDepth() const
uint32_t getNumOctants() const
uint8_t getNchildren() const
void getPreMapping(u32vector &idx, std::vector< int8_t > &markers, std::vector< bool > &isghost)
const uint8_t(* getFacenode() const)[4]
bool getPbound(uint32_t idx, uint8_t face) const
void expectedOctantAdapt(const Octant *octant, int8_t marker, octvector *result) const
bool getParallel() const
Definition ParaTree.cpp:804
const u32arr3vector & getNodes() const
bool adaptGlobalCoarse(bool mapper_flag=false)
uint8_t getFamilySplittingNode(const Octant *) const
bool adapt(bool mapper_flag=false)
const int8_t(* getNormals() const)[3]
Logger & getLog()
Definition ParaTree.cpp:828
void findAllNodeNeighbours(uint32_t idx, uint32_t node, u32vector &neighbours, bvector &isghost)
void findGhostAllCodimensionNeighbours(uint32_t idx, u32vector &neighbours, bvector &isghost)
Octant * getOctant(uint32_t idx)
void getMapping(uint32_t idx, u32vector &mapper, bvector &isghost) const
void loadBalance(DataLBInterface< Impl > &userData, dvector *weight=NULL)
Definition ParaTree.hpp:666
int getMaxLevel() const
octantID getPersistentIdx(uint32_t idx) const
virtual void restore(std::istream &stream)
Definition ParaTree.cpp:588
bool getBalance(uint32_t idx) const
bool getSerial() const
Definition ParaTree.cpp:796
double getLocalMinSize() const
void loadBalance(const dvector *weight=NULL)
void setNofGhostLayers(std::size_t nofGhostLayers)
uint64_t getMorton(uint32_t idx) const
void computeConnectivity()
ParaTree(const std::string &logfile=DEFAULT_LOG_FILE, MPI_Comm comm=MPI_COMM_WORLD)
Definition ParaTree.cpp:136
darray3 getFaceCenter(uint32_t idx, uint8_t face) const
uint8_t getNnodes() const
uint64_t getLastDescMorton() const
int getPointOwnerRank(const darray3 &point)
bool getBound(uint32_t idx, uint8_t face) const
int findOwner(uint64_t morton) const
virtual ~ParaTree()
Definition ParaTree.cpp:211
int getRank() const
Definition ParaTree.cpp:812
static BITPIT_PUBLIC_API const std::string DEFAULT_LOG_FILE
Definition ParaTree.hpp:119
uint32_t getIn(const Intersection *inter) const
const uint8_t(* getEdgeface() const)[2]
Operation getLastOperation() const
Definition ParaTree.cpp:836
const std::vector< uint64_t > & getPartitionLastDesc() const
Definition ParaTree.cpp:975
uint32_t getMaxLength() const
bool adaptGlobalRefine(bool mapper_flag=false)
octantIterator getInternalOctantsBegin()
double getSize(uint32_t idx) const
bool isCommSet() const
Definition ParaTree.cpp:936
const std::vector< uint64_t > & getPartitionFirstDesc() const
Definition ParaTree.cpp:965
LoadBalanceRanges evalLoadBalanceRanges(dvector *weights)
double getY(uint32_t idx) const
const int8_t(* getEdgecoeffs() const)[3]
const std::vector< uint64_t > & getPartitionRangeGlobalIdx() const
Definition ParaTree.cpp:955
virtual void reset()
Definition ParaTree.cpp:448
bvector getPeriodic() const
double getZ0() const
uint32_t getPointOwnerIdx(const double *point) const
void setBalanceCodimension(uint8_t b21codim)
darray3 getNode(uint32_t idx, uint8_t node) const
void writeTest(const std::string &filename, dvector data)
bool isNodeOnOctant(const Octant *nodeOctant, uint8_t nodeIndex, const Octant *octant) const
double getZ(uint32_t idx) const
const std::map< int, std::vector< uint32_t > > & getBordersPerProc() const
double getY0() const
Definition ParaTree.cpp:999
const LoadBalanceRanges & getLoadBalanceRanges() const
MPI_Comm getComm() const
Definition ParaTree.cpp:944
void setTol(double tol=1.0e-14)
bool getFiner(const Intersection *inter) const
uint64_t computeNodePersistentKey(uint32_t idx, uint8_t node) const
double getArea(uint32_t idx) const
double getL() const
bool getIsGhost(const Intersection *inter) const
uint32_t getIdx(const Octant *oct) const
bool isEdgeOnOctant(const Octant *edgeOctant, uint8_t edgeIndex, const Octant *octant) const
std::size_t getNofGhostLayers() const
void findAllCodimensionNeighbours(uint32_t idx, u32vector &neighbours, bvector &isghost)
uint64_t getFirstDescMorton() const
void findGhostNeighbours(uint32_t idx, uint8_t face, uint8_t codim, u32vector &neighbours) const
darray3 getCoordinates(uint32_t idx) const
octantIterator getPboundOctantsBegin()
void loadBalance(DataLBInterface< Impl > &userData, uint8_t &level, dvector *weight=NULL)
Definition ParaTree.hpp:717
uint8_t getFace(const Intersection *inter) const
int8_t getMarker(uint32_t idx) const
uint32_t getNumGhosts() const
const u32vector2D & getGhostConnectivity() const
uint8_t getNedges() const
const uint8_t(* getNodeface() const)[3]
int getOwnerRank(uint64_t globalIdx) const
octantIterator getInternalOctantsEnd()
uint64_t getGlobalNumOctants() const
Definition ParaTree.cpp:788
uint32_t getNumNodes() const
uint32_t getOut(const Intersection *inter) const
void updateConnectivity()
int getNproc() const
Definition ParaTree.cpp:820
uint8_t getNfaces() const
bool getOutIsGhost(const Intersection *inter) const
void communicate(DataCommInterface< Impl > &userData)
Definition ParaTree.hpp:602
uint8_t getLocalMaxDepth() const
int getGhostLayer(const Octant *oct) const
double getVolume(uint32_t idx) const
virtual int getDumpVersion() const
Definition ParaTree.cpp:488
const int8_t(* getNodecoeffs() const)[3]
uint8_t getLevel(uint32_t idx) const
const uint8_t * getOppface() const
void clearConnectivity(bool release=true)
double getX(uint32_t idx) const
void computeIntersections()
void replaceComm(MPI_Comm communicator, MPI_Comm *previousCommunicator=nullptr)
Definition ParaTree.cpp:878
bool getIsNewC(uint32_t idx) const
uint64_t getGlobalIdx(uint32_t idx) const
void getCenter(uint32_t idx, darray3 &centerCoords) const
darray3 getNodeCoordinates(uint32_t node) const
u32vector getOwners(const Intersection *inter) const
void setPeriodic(uint8_t i)
const u32vector2D & getConnectivity() const
void setBalance(uint32_t idx, bool balance)
int8_t getPreMarker(uint32_t idx)
void privateLoadBalance(const uint32_t *partition, DataLBInterface< Impl > *userData=nullptr)
Definition ParaTree.hpp:765
virtual void dump(std::ostream &stream, bool full=true)
Definition ParaTree.cpp:502
uint64_t getGhostGlobalIdx(uint32_t idx) const
void setComm(MPI_Comm communicator)
Definition ParaTree.cpp:848
uint8_t getBalanceCodimension() const
void getNormal(uint32_t idx, uint8_t face, darray3 &normal) const
octantIterator getPboundOctantsEnd()
double getX0() const
Definition ParaTree.cpp:991
uint64_t getStatus() const
Intersection * getIntersection(uint32_t idx)
Octant * getGhostOctant(uint32_t idx)
double getLocalMaxSize() const
darray3 getOrigin() const
Definition ParaTree.cpp:983
bool getIsNewR(uint32_t idx) const
Buffer to be used for receive communications.
Buffer to be used for send communications.
std::unordered_map< int, std::array< uint32_t, 2 > > ExchangeRanges
Definition ParaTree.hpp:121
--- layout: doxygen_footer ---