DNDSR 0.2.1
Distributed Numeric Data Structure for CFV
Loading...
Searching...
No Matches
Mesh_ReadSerializeDistributed.cpp
Go to the documentation of this file.
1#include "Mesh.hpp"
4#include "Geom/Metis.hpp"
5
6#include <array>
7#include <unordered_set>
8
9namespace DNDS::Geom
10{
11 /**
12 * Helper: read a single ArrayPair's father using EvenSplit offset.
13 * Creates father+son, reads father from H5 with even-split.
14 * The son is allocated but empty (size 0).
15 */
16 template <class TArray>
17 static void EvenSplitReadFather(
18 ssp<TArray> &father, ssp<TArray> &son,
19 const MPIInfo &mpi,
21 const std::string &name)
22 {
23 father = std::make_shared<TArray>(mpi);
24 son = std::make_shared<TArray>(mpi);
25 auto offset = Serializer::ArrayGlobalOffset_EvenSplit;
26 father->ReadSerializer(serializerP, name, offset);
27 }
28
29 /// Overload for arrays that need CommType/CommMult (ElemInfo, NodePeriodicBits).
30 template <class TArray>
31 static void EvenSplitReadFather(
32 ssp<TArray> &father, ssp<TArray> &son,
33 const MPIInfo &mpi,
35 const std::string &name,
36 MPI_Datatype commType, int commMult)
37 {
38 father = std::make_shared<TArray>(commType, commMult, mpi);
39 son = std::make_shared<TArray>(commType, commMult, mpi);
40 auto offset = Serializer::ArrayGlobalOffset_EvenSplit;
41 father->ReadSerializer(serializerP, name, offset);
42 }
43
44 // =================================================================
45 // Top-level orchestrator
46 // =================================================================
47
51 const std::string &name,
53 {
55 "ReadSerializeAndDistribute requires a collective (H5) serializer");
56
57 if (mpi.rank == 0)
58 log() << "UnstructuredMesh === ReadSerializeAndDistribute: begin" << std::endl;
59
60 auto cwd = serializerP->GetCurrentPath();
61 serializerP->GoToPath(name);
62
63 // Step 0+1: read metadata and all primary arrays with even-split
64 ReadDistributed_EvenSplitRead(serializerP);
65 serializerP->GoToPath(cwd);
66
67 // Step 2+3: build adjacencies and filter to facial cell2cell
68 auto cell2cellFacial = ReadDistributed_BuildFacialCell2Cell();
69
70 // Step 4: ParMetis repartition
71 auto cellPartition = ReadDistributed_PartitionParMetis(
73 cell2cellFacial.reset(); // free before redistribution
74
75 // Step 5: derive entity partitions and redistribute
76 auto partitions = ReadDistributed_DeriveEntityPartitions(std::move(cellPartition));
77 ReadDistributed_Redistribute(partitions);
78
79 // Step 6: log final stats
80 coords.father->createGlobalMapping();
81 cell2node.father->createGlobalMapping();
82 bnd2node.father->createGlobalMapping();
83
84 {
85 index nCellG = this->NumCellGlobal();
86 index nNodeG = this->NumNodeGlobal();
87 index nBndG = this->NumBndGlobal();
88 if (mpi.rank == 0)
89 {
90 log() << "UnstructuredMesh === ReadSerializeAndDistribute done, "
91 << fmt::format("nCellGlobal [{}], nNodeGlobal [{}], nBndGlobal [{}]",
93 << std::endl;
94 }
95 MPISerialDo(mpi, [&]()
96 { log() << fmt::format(" Rank {}: nCell {}, nNode {}, nBnd {}",
97 mpi.rank, this->NumCell(), this->NumNode(), this->NumBnd())
98 << std::endl; });
99 }
100 }
101
102 // =================================================================
103 // Step 0+1: Even-split read
104 // =================================================================
105
109 {
110 // Read scalar metadata
111 {
112 std::string meshRead;
113 index dimRead{0}, sizeRead{0};
114 int isPeriodicRead = 0;
115 serializerP->ReadString("mesh", meshRead);
116 serializerP->ReadIndex("dim", dimRead);
117 serializerP->ReadIndex("MPISize", sizeRead);
118 serializerP->ReadInt("isPeriodic", isPeriodicRead);
120 DNDS_assert(meshRead == "UnstructuredMesh");
122 }
123
124 if (mpi.rank == 0)
125 log() << "UnstructuredMesh === ReadSerializeAndDistribute: even-split read" << std::endl;
126
127 auto innerCwd = serializerP->GetCurrentPath();
128
129 auto readFatherAt = [&](auto &father, auto &son, const std::string &group)
130 {
131 serializerP->GoToPath(innerCwd);
132 serializerP->GoToPath(group);
133 EvenSplitReadFather(father, son, mpi, serializerP, "father");
134 };
135
136 auto readFatherAtTyped = [&](auto &father, auto &son, const std::string &group,
138 {
139 serializerP->GoToPath(innerCwd);
140 serializerP->GoToPath(group);
141 EvenSplitReadFather(father, son, mpi, serializerP, "father", commType, commMult);
142 };
143
144 readFatherAt(coords.father, coords.son, "coords");
145 readFatherAt(cell2node.father, cell2node.son, "cell2node");
148 readFatherAt(bnd2node.father, bnd2node.son, "bnd2node");
151
152 if (isPeriodic)
153 {
158 serializerP->GoToPath(innerCwd);
159 periodicInfo.ReadSerializer(serializerP, "periodicInfo");
160 }
161
165
167 cell2node.idx.markGlobal();
168 bnd2node.idx.markGlobal();
169
170 if (mpi.rank == 0)
171 log() << "UnstructuredMesh === ReadSerializeAndDistribute: even-split read done, "
172 << "nCellLocal=" << cell2node.father->Size()
173 << " nNodeLocal=" << coords.father->Size()
174 << " nBndLocal=" << bnd2node.father->Size() << std::endl;
175 }
176
177 // =================================================================
178 // Step 2+3: Build facial cell2cell
179 // =================================================================
180
181 ssp<tAdj::element_type> UnstructuredMesh::
183 {
184 if (mpi.rank == 0)
185 log() << "UnstructuredMesh === ReadSerializeAndDistribute: building distributed cell2cell" << std::endl;
186
187 // Step 2: build node-neighbor cell2cell and bnd2cell via DSL.
190
191 if (mpi.rank == 0)
192 log() << "UnstructuredMesh === ReadSerializeAndDistribute: building facial cell2cell" << std::endl;
193
194 // Ghost-pull cell2node and cellElemInfo for all neighbor cells.
195 {
196 cell2cell.TransAttach();
197 cell2cell.trans.createFatherGlobalMapping();
198
199 MeshConnectivity dagTmp;
201
203 auto cellResult = dagTmp.evaluateGhostTree(
205
206 auto it = cellResult.ghostIndices.find(EntityKind::Cell);
207 std::vector<index> ghostCells;
208 if (it != cellResult.ghostIndices.end())
209 ghostCells = std::move(it->second);
210
211 cell2node.TransAttach();
212 cell2node.trans.createFatherGlobalMapping();
214 cellElemInfo.trans.createFatherGlobalMapping();
215
216 cell2cell.trans.createGhostMapping(ghostCells);
217 cell2node.trans.BorrowGGIndexing(cell2cell.trans);
218 cellElemInfo.trans.BorrowGGIndexing(cell2cell.trans);
219
220 cell2cell.trans.createMPITypes();
221 cell2node.trans.createMPITypes();
222 cellElemInfo.trans.createMPITypes();
223
224 cell2cell.trans.pullOnce();
225 cell2node.trans.pullOnce();
226 cellElemInfo.trans.pullOnce();
227 }
228
229 // Step 3: filter cell2cell to face-sharing neighbors (O1 vertex intersection >= dim).
230 auto cell2cellFacial = make_ssp<tAdj::element_type>(ObjName{"cell2cellFacialTmp"}, mpi);
231 cell2cellFacial->Resize(cell2cell.father->Size(), 6);
232
233 for (index iCell = 0; iCell < cell2cell.father->Size(); iCell++)
234 {
235 auto elemI = Elem::Element{cellElemInfo(iCell, 0).getElemType()};
236 rowsize nVertI = elemI.GetNumVertices();
237 std::vector<index> vertI(cell2node[iCell].begin(),
238 cell2node[iCell].begin() + nVertI);
239 std::sort(vertI.begin(), vertI.end());
240
241 std::vector<index> facialNeighbors;
242 facialNeighbors.reserve(6);
243 for (rowsize ic2c = 0; ic2c < cell2cell.father->RowSize(iCell); ic2c++)
244 {
245 index iCellOther = (*cell2cell.father)(iCell, ic2c);
246
247 auto [found, rank, localAppend] =
248 cell2cell.trans.pLGhostMapping->search_indexAppend(iCellOther);
249 DNDS_assert_info(found, "cell2cell neighbor not found in ghost mapping");
250
251 auto elemJ = Elem::Element{cellElemInfo(localAppend, 0).getElemType()};
252 rowsize nVertJ = elemJ.GetNumVertices();
253 std::vector<index> vertJ(cell2node[localAppend].begin(),
254 cell2node[localAppend].begin() + nVertJ);
255 std::sort(vertJ.begin(), vertJ.end());
256
257 std::vector<index> intersect;
258 intersect.reserve(9);
259 std::set_intersection(vertI.begin(), vertI.end(),
260 vertJ.begin(), vertJ.end(),
261 std::back_inserter(intersect));
262 if (static_cast<int>(intersect.size()) >= dim)
263 facialNeighbors.push_back(iCellOther);
264 }
265 cell2cellFacial->ResizeRow(iCell, facialNeighbors.size());
268 }
269 cell2cellFacial->Compress();
270
271 return cell2cellFacial;
272 }
273
274 // =================================================================
275 // Step 4: ParMetis repartition
276 // =================================================================
277
278 std::vector<MPI_int> UnstructuredMesh::
280 const ssp<tAdj::element_type> &cell2cellFacial,
281 const PartitionOptions &partitionOptions)
282 {
283 if (mpi.rank == 0)
284 log() << "UnstructuredMesh === ReadSerializeAndDistribute: ParMetis repartition" << std::endl;
285
286 cell2cellFacial->createGlobalMapping();
287 idx_t nPart = mpi.size;
288
289 std::vector<idx_t> vtxdist(mpi.size + 1);
290 for (MPI_int r = 0; r <= mpi.size; r++)
291 vtxdist[r] = METIS::indexToIdx(cell2cellFacial->pLGlobalMapping->ROffsets().at(r));
292
293 std::vector<idx_t> xadj(cell2cellFacial->Size() + 1);
294 for (index i = 0; i <= cell2cellFacial->Size(); i++)
295 xadj[i] = METIS::indexToIdx(cell2cellFacial->rowPtr(i) - cell2cellFacial->rowPtr(0));
296
297 std::vector<idx_t> adjncy(xadj.back());
298 for (index i = 0; i < xadj.back(); i++)
299 adjncy[i] = METIS::indexToIdx(cell2cellFacial->data()[i]);
300
301 if (adjncy.empty())
302 adjncy.resize(1, -1); // cope with zero-sized data
303
304 std::vector<MPI_int> cellPartition(cell2cellFacial->Size());
305
306 if (nPart > 1)
307 {
308 for (int i = 0; i < mpi.size; i++)
309 DNDS_assert_info(vtxdist[i + 1] - vtxdist[i] > 0,
310 "ParMetis requires > 0 cells on each proc");
311
312 idx_t nCon{1};
313 idx_t wgtflag{0}, numflag{0};
314 std::vector<real_t> tpWeights(static_cast<size_t>(nPart) * nCon, 1.0 / nPart);
315 std::array<real_t, 1> ubVec{1.05};
316 std::array<idx_t, 3> optsC{1, 0, static_cast<idx_t>(partitionOptions.metisSeed)};
317 idx_t objval = 0;
318 std::vector<idx_t> partOut(cell2cellFacial->Size());
319 if (partOut.empty())
320 partOut.resize(1, 0);
321
323 vtxdist.data(), xadj.data(), adjncy.data(),
324 nullptr, nullptr, &wgtflag, &numflag,
325 &nCon, &nPart, tpWeights.data(), ubVec.data(), optsC.data(),
326 &objval, partOut.data(), &mpi.comm);
328 fmt::format("ParMETIS_V3_PartKway returned {}", ret));
329
330 for (index i = 0; i < cell2cellFacial->Size(); i++)
331 cellPartition[i] = static_cast<MPI_int>(partOut[i]);
332 }
333 else
334 {
335 std::fill(cellPartition.begin(), cellPartition.end(), 0);
336 }
337
338 // Print partition stats
339 {
340 std::vector<index> localPartCnt(nPart, 0);
341 for (auto p : cellPartition)
342 localPartCnt.at(p)++;
343 std::vector<index> globalPartCnt(nPart, 0);
345 if (mpi.rank == 0)
346 {
347 auto [minIt, maxIt] = std::minmax_element(globalPartCnt.begin(), globalPartCnt.end());
348 index nCellGlobal = 0;
349 for (auto c : globalPartCnt)
350 nCellGlobal += c;
351 log() << "UnstructuredMesh === ReadSerializeAndDistribute: partition done, "
352 << fmt::format("nCellGlobal [{}], ave [{}], min [{}], max [{}], ratio [{:.4f}]",
353 nCellGlobal, real(nCellGlobal) / nPart, *minIt, *maxIt, real(*minIt) / *maxIt)
354 << std::endl;
355 }
356 }
357
358 return cellPartition;
359 }
360
361 // =================================================================
362 // Step 5a: Derive node and bnd partitions from cell partition
363 // =================================================================
364
365 UnstructuredMesh::EntityPartitions UnstructuredMesh::
367 std::vector<MPI_int> cellPartition)
368 {
369 EntityPartitions result;
370 result.cellPartition = std::move(cellPartition);
371
372 // Node partition: each node goes to the min cell partition that claims it.
373 result.nodePartition.assign(coords.father->Size(), static_cast<MPI_int>(INT32_MAX));
374 for (index iCell = 0; iCell < cell2node.father->Size(); iCell++)
375 {
376 for (rowsize ic2n = 0; ic2n < cell2node.father->RowSize(iCell); ic2n++)
377 {
378 index iNode = (*cell2node.father)(iCell, ic2n);
379 MPI_int rank = UnInitMPIInt;
380 index val = UnInitIndex;
381 bool found = coords.father->pLGlobalMapping->search(iNode, rank, val);
383 if (rank == mpi.rank)
384 result.nodePartition[val] = std::min(result.nodePartition[val],
385 result.cellPartition.at(iCell));
386 }
387 }
388 for (auto &p : result.nodePartition)
390 p = 0;
391
392 // Bnd partition: bnd goes to same rank as its owner cell (bnd2cell(iBnd, 0)).
393 // Build a distributed lookup for the owner cell's target partition.
394 auto cellPartArr = make_ssp<tAdj1::element_type>(ObjName{"cellPartArr"}, mpi);
395 auto cellPartArrGhost = make_ssp<tAdj1::element_type>(ObjName{"cellPartArrGhost"}, mpi);
396 cellPartArr->Resize(cell2node.father->Size());
397 for (index i = 0; i < cellPartArr->Size(); i++)
398 (*cellPartArr)(i, 0) = static_cast<index>(result.cellPartition.at(i));
399 cellPartArr->createGlobalMapping();
400
401 std::vector<index> bndOwnerGhostQuery;
402 for (index iBnd = 0; iBnd < bnd2node.father->Size(); iBnd++)
403 {
407 bool found = cellPartArr->pLGlobalMapping->search(iOwnerCell, ownerRank, ownerVal);
409 if (ownerRank != mpi.rank)
411 }
414 cellPartTrans.createGhostMapping(bndOwnerGhostQuery);
415 cellPartTrans.createMPITypes();
416 cellPartTrans.pullOnce();
417
418 result.bndPartition.resize(bnd2node.father->Size());
419 for (index iBnd = 0; iBnd < bnd2node.father->Size(); iBnd++)
420 {
423 fmt::format("bnd {} has no owner cell after RecoverCell2CellAndBnd2Cell", iBnd));
426 bool foundOwner = cellPartArr->pLGlobalMapping->search(iOwnerCell, ownerRank, ownerVal);
429 if (ownerRank == mpi.rank)
430 targetPart = (*cellPartArr)(ownerVal, 0);
431 else
432 {
433 auto [found, gRank, gVal] = cellPartTrans.pLGhostMapping->search_indexAppend(iOwnerCell);
435 if (gRank == -1)
436 targetPart = (*cellPartArr)(gVal, 0);
437 else
438 targetPart = (*cellPartArrGhost)(gVal - cellPartArr->Size(), 0);
439 }
440 result.bndPartition[iBnd] = static_cast<MPI_int>(targetPart);
441 }
442
443 return result;
444 }
445
446 // =================================================================
447 // Step 5b-5d: Redistribute arrays to new partition
448 // =================================================================
449
452 const EntityPartitions &partitions)
453 {
454 if (mpi.rank == 0)
455 log() << "UnstructuredMesh === ReadSerializeAndDistribute: redistributing" << std::endl;
456
457 // Free temporary adjacencies no longer needed (same as legacy).
458 cell2cell.reset();
459 node2cell.reset();
460 node2bnd.reset();
461 bnd2cell.reset();
462
463 // Use ReorderEntities with all three partitions as explicit maps.
464 // (No follow computation needed — partitions are pre-computed by
465 // ReadDistributed_DeriveEntityPartitions.)
466 ReorderInput input;
467 input.explicitMaps.push_back(EntityReorderMap{EntityKind::Cell, partitions.cellPartition});
468 input.explicitMaps.push_back(EntityReorderMap{EntityKind::Node, partitions.nodePartition});
469 input.explicitMaps.push_back(EntityReorderMap{EntityKind::Bnd, partitions.bndPartition});
470 // No follows (all explicit). No destroyKinds (no faces at this stage).
471
472 this->ReorderEntities(input);
473
474 if (mpi.rank == 0)
475 log() << "UnstructuredMesh === ReadSerializeAndDistribute: redistribute done" << std::endl;
476 }
477
480 const EntityPartitions &partitions)
481 {
482 if (mpi.rank == 0)
483 log() << "UnstructuredMesh === ReadSerializeAndDistribute: redistributing" << std::endl;
484
485 // Free temporary adjacencies no longer needed.
486 cell2cell.reset();
487 node2cell.reset();
488 node2bnd.reset();
489 bnd2cell.reset();
490
491 // Compute push indices and serial-to-global reordering.
492 std::vector<index> cell_push, cell_pushStart;
493 std::vector<index> node_push, node_pushStart;
494 std::vector<index> bnd_push, bnd_pushStart;
498
499 std::vector<index> cell_Serial2Global, node_Serial2Global;
502
503 // Convert adjacency indices to new global numbering.
506
507 // Transfer all arrays to new partition.
508 auto coordsOld = coords.father;
509 auto cell2nodeOld = cell2node.father;
511 auto bnd2nodeOld = bnd2node.father;
516
517 coords.InitPair("coords", mpi);
518 cell2node.InitPair("cell2node", mpi);
520 bnd2node.InitPair("bnd2node", mpi);
522 cell2cellOrig.InitPair("cell2cellOrig", mpi);
523 node2nodeOrig.InitPair("node2nodeOrig", mpi);
524 bnd2bndOrig.InitPair("bnd2bndOrig", mpi);
525
528
532
536
537 if (isPeriodic)
538 {
545 }
546
548 cell2node.idx.markGlobal();
549 bnd2node.idx.markGlobal();
550 }
551
552} // namespace DNDS::Geom
#define DNDS_assert_info(expr, info)
Debug-only assertion with an extra std::string info message.
Definition Errors.hpp:117
#define DNDS_assert(expr)
Debug-only assertion (compiled out when DNDS_NDEBUG is defined). Prints the expression + file/line + ...
Definition Errors.hpp:112
#define DNDS_check_throw_info(expr, info)
Same as DNDS_check_throw but attaches a user-supplied info message to the thrown std::runtime_error.
Definition Errors.hpp:100
Layered DAG of mesh adjacency relations with composable DSL operations.
Internal helpers for mesh partition redistribution.
constexpr AdjKind Cell2Cell
void Partition2LocalIdx(const std::vector< TPartitionIdx > &partition, std::vector< DNDS::index > &localPush, std::vector< DNDS::index > &localPushStart, const DNDS::MPIInfo &mpi)
void Partition2Serial2Global(const std::vector< TPartitionIdx > &partition, std::vector< DNDS::index > &serial2Global, const DNDS::MPIInfo &mpi, DNDS::MPI_int nPart)
void TransferDataSerial2Global(TArr &arraySerial, TArr &arrayDist, const std::vector< DNDS::index > &pushIndex, const std::vector< DNDS::index > &pushIndexStart, const DNDS::MPIInfo &mpi)
void ConvertAdjSerial2Global(TAdj &arraySerialAdj, const std::vector< DNDS::index > &partitionJSerial2Global, const DNDS::MPIInfo &mpi)
MPI_int Allreduce(const void *sendbuf, void *recvbuf, MPI_int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
Wrapper over MPI_Allreduce.
Definition MPI.cpp:202
ssp< SerializerBase > SerializerBaseSSP
const MPI_Datatype DNDS_MPI_INDEX
MPI datatype matching index (= MPI_INT64_T).
Definition MPI.hpp:106
DNDS_CONSTANT const index UnInitIndex
Sentinel "not initialised" index value (= INT64_MIN).
Definition Defines.hpp:181
void MPISerialDo(const MPIInfo &mpi, F f)
Execute f on each rank serially, in rank order.
Definition MPI.hpp:749
int32_t rowsize
Row-width / per-row element-count type (signed 32-bit).
Definition Defines.hpp:114
constexpr MPI_int UnInitMPIInt
Sentinel "not initialised" MPI_int value (= -1, invalid rank).
Definition MPI.hpp:66
int64_t index
Global row / DOF index type (signed 64-bit; handles multi-billion-cell meshes).
Definition Defines.hpp:112
double real
Canonical floating-point scalar used throughout DNDSR (double precision).
Definition Defines.hpp:110
std::ostream & log()
Return the current DNDSR log stream (either std::cout or the installed file).
Definition Defines.cpp:50
int MPI_int
MPI counterpart type for MPI_int (= C int). Used for counts and ranks in MPI calls.
Definition MPI.hpp:54
void TransAttach()
Bind the transformer to the current father / son pointers.
ssp< TArray > father
Owned-side array (must be resized before ghost setup).
void InitPair(const std::string &name, Args &&...args)
Allocate both father and son arrays, forwarding all args to TArray constructor.
ssp< TArray > son
Ghost-side array (sized automatically by createMPITypes / BorrowAndPull).
TTrans trans
Ghost-communication engine bound to father and son.
ArrayTransformer< typename TArray::value_type, TArray::rs, TArray::rm, TArray::al > Type
static CompiledGhostTree compile(const GhostSpec &spec)
static MPI_Datatype CommType()
static MPI_Datatype CommType()
DNDS_HOST void ReadSerializer(Serializer::SerializerBaseSSP serializerP, const std::string &name)
UnstructuredMeshDeviceView< B > t_deviceView
Definition Mesh.hpp:1036
tPbiPair cell2nodePbi
periodic only, after reader
Definition Mesh.hpp:68
void fillRegistry(MeshConnectivity &dag) const
Populate a MeshConnectivity registry from this mesh's currently-built adjacencies.
Definition Mesh.cpp:1743
void ReorderEntities(const ReorderInput &input)
Reorder entities using the general framework.
tCoordPair coords
reader
Definition Mesh.hpp:57
void RecoverNode2CellAndNode2Bnd()
only requires father part of cell2node, bnd2node and coords generates node2cell and node2bnd (father ...
Definition Mesh.cpp:379
AdjPairTracked< tAdjPair > cell2cell
Definition Mesh.hpp:61
MeshAdjState adjPrimaryState
Definition Mesh.hpp:44
AdjPairTracked< tAdjPair > cell2node
Definition Mesh.hpp:58
void RecoverCell2CellAndBnd2Cell()
needs to use RecoverNode2CellAndNode2Bnd before doing this. Requires node2cell.father and builds a ve...
Definition Mesh.cpp:422
tElemInfoArrayPair cellElemInfo
Definition Mesh.hpp:62
AdjPairTracked< tAdjPair > node2bnd
Definition Mesh.hpp:87
void ReadSerializeAndDistribute(Serializer::SerializerBaseSSP serializerP, const std::string &name, const PartitionOptions &partitionOptions)
Reads mesh from an H5 serializer using even-split distribution, then repartitions via ParMetis for lo...
AdjPairTracked< tAdjPair > bnd2node
Definition Mesh.hpp:59
AdjPairTracked< tAdjPair > node2cell
inverse relations
Definition Mesh.hpp:86
tElemInfoArrayPair bndElemInfo
Definition Mesh.hpp:63
AdjPairTracked< tAdj2Pair > bnd2cell
Definition Mesh.hpp:60
int size
Number of ranks in comm (-1 until initialised).
Definition MPI.hpp:237
int rank
This rank's 0-based index within comm (-1 until initialised).
Definition MPI.hpp:235
MPI_Comm comm
The underlying MPI communicator handle.
Definition MPI.hpp:233
if(DEFINED ENV{DNDS_TEST_NP_LIST}) set(DNDS_TEST_NP_LIST "$ENV
tVec r(NCells)
auto result
std::vector< GhostCell > ghostCells
std::vector< MPI_int > nodePartition(nNode, mpi.rank)
std::vector< MPI_int > cellPartition(nCellBefore)
ReorderInput input
const tPoint const tPoint const tPoint & p