DNDSR 0.1.0.dev1+gcd065ad
Distributed Numeric Data Structure for CFV
Loading...
Searching...
No Matches
test_MeshIndexConversion.cpp
Go to the documentation of this file.
1/**
2 * @file test_MeshIndexConversion.cpp
3 * @brief Parameterized doctest-based unit tests for UnstructuredMesh index
4 * conversion and adjacency Global2Local / Local2Global methods.
5 *
6 * Tests run against 4 mesh configurations (same as test_MeshPipeline):
7 * [0] UniformSquare_10 -- 2D, 100 quad cells, non-periodic
8 * [1] IV10_10 -- 2D, 100 quad cells, periodic (isentropic vortex)
9 * [2] NACA0012_H2 -- 2D, 20816 unstructured quad cells, non-periodic
10 * [3] IV10U_10 -- 2D, 322 unstructured tri cells, periodic
11 *
12 * Each mesh is partitioned with Metis, ghost layers built, then exercises:
13 * 1. Per-entity index conversions (Global2Local, Local2Global, _NoSon)
14 * 2. Adjacency state-machine round-trips (Primary, ForBnd)
15 * 3. Round-trip stability under repeated conversions
16 * 4. UnInitIndex pass-through
17 * 5. Negative encoding for not-found globals
18 * 6. Local index range validity after AdjGlobal2LocalPrimary
19 *
20 * All meshes are pre-built in main() and destroyed together after tests.
21 * This avoids an OpenMPI resource exhaustion observed when building/destroying
22 * 20+ meshes in a single process at np=8.
23 *
24 * Run under mpirun with 1, 2, 4, and 8 ranks.
25 */
26
27#define DOCTEST_CONFIG_IMPLEMENT
28#include "doctest.h"
29
30#include "Geom/Mesh.hpp"
31#include <string>
32#include <vector>
33
34using namespace DNDS;
35using namespace DNDS::Geom;
36
37// NOTE: DNDS::index, DNDS::real, DNDS::rowsize clash with POSIX symbols
38// from <strings.h> (pulled in by doctest/MPI). Every declaration of these
39// types must use the fully-qualified DNDS:: prefix. See AGENTS.md.
40
41// ---------------------------------------------------------------------------
42// Mesh configuration (shared with test_MeshPipeline)
43// ---------------------------------------------------------------------------
44struct MeshConfig
45{
46 const char *name; ///< Human-readable label
47 const char *file; ///< Filename relative to data/mesh/
48 int dim; ///< Spatial dimension
49 bool periodic; ///< Has periodic boundaries
50 tPoint translation1; ///< Periodic translation axis 1
51 tPoint translation2; ///< Periodic translation axis 2
52 tPoint translation3; ///< Periodic translation axis 3
53};
54
55static const MeshConfig g_configs[] = {
56 {"UniformSquare_10", "UniformSquare_10.cgns", 2, false,
57 {0, 0, 0}, {0, 0, 0}, {0, 0, 0}},
58 {"IV10_10", "IV10_10.cgns", 2, true,
59 {10, 0, 0}, {0, 10, 0}, {0, 0, 10}},
60 {"NACA0012_H2", "NACA0012_H2.cgns", 2, false,
61 {0, 0, 0}, {0, 0, 0}, {0, 0, 0}},
62 {"IV10U_10", "IV10U_10.cgns", 2, true,
63 {10, 0, 0}, {0, 10, 0}, {0, 0, 10}},
64 {"Ball2", "Ball2.cgns", 3, false,
65 {0, 0, 0}, {0, 0, 0}, {0, 0, 0}}, // 3D mixed-element mesh (Ball2 is np=8 only)
66};
67static constexpr int N_CONFIGS = sizeof(g_configs) / sizeof(g_configs[0]);
68
69// ---------------------------------------------------------------------------
70// Globals
71// ---------------------------------------------------------------------------
72static MPIInfo g_mpi;
73
74/// One mesh per config for read-only tests.
75static ssp<UnstructuredMesh> g_meshes[N_CONFIGS];
76
77// ---------------------------------------------------------------------------
78// Helpers
79// ---------------------------------------------------------------------------
80
81static std::string meshPath(const std::string &name)
82{
83 std::string f(__FILE__);
84 for (int i = 0; i < 4; i++)
85 {
86 auto pos = f.rfind('/');
87 if (pos == std::string::npos)
88 pos = f.rfind('\\');
89 if (pos != std::string::npos)
90 f = f.substr(0, pos);
91 }
92 return f + "/data/mesh/" + name;
93}
94
95static ssp<UnstructuredMesh> buildMesh(int meshId, const MeshConfig &cfg)
96{
97 if (g_mpi.rank == 0)
98 std::cout << "[buildMesh " << meshId << " " << cfg.name << "] START" << std::endl;
99
100 auto mesh = std::make_shared<UnstructuredMesh>(g_mpi, cfg.dim);
101 UnstructuredMeshSerialRW reader(mesh, 0);
102
103 if (cfg.periodic)
104 {
105 tPoint zero{0, 0, 0};
106 mesh->SetPeriodicGeometry(
107 cfg.translation1, zero, zero,
108 cfg.translation2, zero, zero,
109 cfg.translation3, zero, zero);
110 }
111
112 if (g_mpi.rank == 0)
113 std::cout << "[buildMesh " << meshId << " " << cfg.name << "] ReadFromCGNSSerial..." << std::endl;
114 reader.ReadFromCGNSSerial(meshPath(cfg.file));
115
116 if (g_mpi.rank == 0)
117 std::cout << "[buildMesh " << meshId << " " << cfg.name << "] Deduplicate1to1Periodic..." << std::endl;
118 reader.Deduplicate1to1Periodic(1e-8);
119
120 if (g_mpi.rank == 0)
121 std::cout << "[buildMesh " << meshId << " " << cfg.name << "] BuildCell2Cell..." << std::endl;
122 reader.BuildCell2Cell();
123
125 pOpt.metisType = "KWAY";
126 pOpt.metisUfactor = 30;
127 pOpt.metisSeed = 42;
128 pOpt.metisNcuts = 1;
129
130 if (g_mpi.rank == 0)
131 std::cout << "[buildMesh " << meshId << " " << cfg.name << "] MeshPartitionCell2Cell..." << std::endl;
132 reader.MeshPartitionCell2Cell(pOpt);
133
134 if (g_mpi.rank == 0)
135 std::cout << "[buildMesh " << meshId << " " << cfg.name << "] PartitionReorderToMeshCell2Cell..." << std::endl;
136 reader.PartitionReorderToMeshCell2Cell();
137
138 if (g_mpi.rank == 0)
139 std::cout << "[buildMesh " << meshId << " " << cfg.name << "] RecoverNode2CellAndNode2Bnd..." << std::endl;
140 mesh->RecoverNode2CellAndNode2Bnd();
141
142 if (g_mpi.rank == 0)
143 std::cout << "[buildMesh " << meshId << " " << cfg.name << "] RecoverCell2CellAndBnd2Cell..." << std::endl;
144 mesh->RecoverCell2CellAndBnd2Cell();
145
146 if (g_mpi.rank == 0)
147 std::cout << "[buildMesh " << meshId << " " << cfg.name << "] BuildGhostPrimary..." << std::endl;
148 mesh->BuildGhostPrimary();
149
150 if (g_mpi.rank == 0)
151 std::cout << "[buildMesh " << meshId << " " << cfg.name << "] AdjGlobal2LocalPrimary..." << std::endl;
152 mesh->AdjGlobal2LocalPrimary();
153
154 if (g_mpi.rank == 0)
155 std::cout << "[buildMesh " << meshId << " " << cfg.name << "] DONE" << std::endl;
156
157 return mesh;
158}
159
160static std::vector<std::vector<DNDS::index>> snapshotAdj(
161 const tAdjPair &adj, DNDS::index nRows)
162{
163 std::vector<std::vector<DNDS::index>> out(nRows);
164 for (DNDS::index i = 0; i < nRows; i++)
165 {
166 out[i].resize(adj.RowSize(i));
167 for (DNDS::rowsize j = 0; j < adj.RowSize(i); j++)
168 out[i][j] = adj(i, j);
169 }
170 return out;
171}
172
173// ---------------------------------------------------------------------------
174// Parameterized test macro
175// ---------------------------------------------------------------------------
176
177#define FOR_EACH_MESH_CONFIG(body) \
178 for (int _ci = 0; _ci < N_CONFIGS; _ci++) \
179 { \
180 /* Ball2 (index 4) is only tested with np=8 */ \
181 if (_ci == 4 && g_mpi.size != 8) continue; \
182 CAPTURE(_ci); \
183 const auto &cfg = g_configs[_ci]; \
184 auto m = g_meshes[_ci]; \
185 SUBCASE(cfg.name) { body } \
186 }
187
188// ---------------------------------------------------------------------------
189int main(int argc, char **argv)
190{
191 MPI_Init(&argc, &argv);
192 g_mpi.setWorld();
193
194 // Pre-build all mesh configs for read-only tests.
195 for (int i = 0; i < N_CONFIGS; i++)
196 {
197 // Skip Ball2 (index 4) for np < 8 to avoid timeout
198 if (i == 4 && g_mpi.size != 8)
199 continue;
200 g_meshes[i] = buildMesh(i, g_configs[i]);
201 }
202
203 doctest::Context ctx;
204 ctx.applyCommandLine(argc, argv);
205 int res = ctx.run();
206
207 // Destroy all meshes together before MPI_Finalize.
208 for (auto &m : g_meshes)
209 m.reset();
210 MPI_Finalize();
211 return res;
212}
213
214// ===========================================================================
215// Mesh sanity checks (parameterized over all configs)
216// ===========================================================================
217TEST_CASE("Mesh setup: every rank has cells and nodes")
218{
220 CHECK(m->NumCell() >= 1);
221 CHECK(m->NumNode() >= 1);
222 CHECK(m->adjPrimaryState == Adj_PointToLocal);
223 })
224}
225
226TEST_CASE("Mesh setup: multi-rank produces ghosts")
227{
228 if (g_mpi.size < 2)
229 return;
231 CHECK(m->NumCellGhost() > 0);
232 CHECK(m->NumNodeGhost() > 0);
233 })
234}
235
236TEST_CASE("Mesh setup: periodic state matches configuration")
237{
239 // cfg.periodic indicates if the mesh was configured as periodic
240 // After Deduplicate1to1Periodic, isPeriodic should match
241 CHECK(m->isPeriodic == cfg.periodic);
242 })
243}
244
245TEST_CASE("Ball2: mixed 3D element types")
246{
247 // This test only runs for Ball2 (index 4) with np=8
248 if (g_mpi.size != 8)
249 return;
250
251 auto m = g_meshes[4]; // Ball2 is index 4
252 DNDS_assert(m != nullptr);
253
254 // Count element types locally
255 std::map<Elem::ElemType, int> typeCounts;
256 for (DNDS::index iC = 0; iC < m->NumCell(); iC++)
257 {
258 auto elem = m->GetCellElement(iC);
259 typeCounts[elem.type]++;
260 }
261
262 // Gather all unique element types from all ranks
263 std::vector<int> localTypes;
264 for (const auto& [type, count] : typeCounts)
265 localTypes.push_back(static_cast<int>(type));
266
267 int nLocalTypes = localTypes.size();
268 std::vector<int> allRankSizes(g_mpi.size);
269 MPI_Gather(&nLocalTypes, 1, MPI_INT, allRankSizes.data(), 1, MPI_INT, 0, g_mpi.comm);
270
271 std::vector<int> allTypes;
272 std::vector<int> displs;
273 int totalTypes = 0;
274
275 if (g_mpi.rank == 0)
276 {
277 displs.resize(g_mpi.size);
278 for (int i = 0; i < g_mpi.size; i++)
279 {
280 displs[i] = totalTypes;
281 totalTypes += allRankSizes[i];
282 }
283 allTypes.resize(totalTypes);
284 }
285
286 MPI_Gatherv(localTypes.data(), nLocalTypes, MPI_INT,
287 allTypes.data(), allRankSizes.data(), displs.data(), MPI_INT, 0, g_mpi.comm);
288
289 // Build set of unique types (on rank 0)
290 std::set<Elem::ElemType> uniqueTypes;
291 if (g_mpi.rank == 0)
292 {
293 for (int typeInt : allTypes)
294 uniqueTypes.insert(static_cast<Elem::ElemType>(typeInt));
295 }
296
297 // Broadcast unique types to all ranks
298 int nUnique = uniqueTypes.size();
299 MPI_Bcast(&nUnique, 1, MPI_INT, 0, g_mpi.comm);
300 std::vector<int> uniqueTypesVec;
301 if (g_mpi.rank == 0)
302 {
303 for (auto type : uniqueTypes)
304 uniqueTypesVec.push_back(static_cast<int>(type));
305 }
306 uniqueTypesVec.resize(nUnique);
307 MPI_Bcast(uniqueTypesVec.data(), nUnique, MPI_INT, 0, g_mpi.comm);
308
309 // All ranks iterate over the same types in the same order
310 std::map<Elem::ElemType, int> globalCounts;
311 for (int typeInt : uniqueTypesVec)
312 {
313 Elem::ElemType type = static_cast<Elem::ElemType>(typeInt);
314 int localCount = typeCounts.count(type) ? typeCounts[type] : 0;
315 int globalCount;
316 MPI_Reduce(&localCount, &globalCount, 1, MPI_INT, MPI_SUM, 0, g_mpi.comm);
317 if (g_mpi.rank == 0)
318 globalCounts[type] = globalCount;
319 }
320
321 if (g_mpi.rank == 0)
322 {
323 std::cout << "\nBall2 element type counts:" << std::endl;
324 for (const auto& [type, count] : globalCounts)
325 {
326 std::cout << " Type " << static_cast<int>(type) << ": " << count << std::endl;
327 }
328
329 // Verify expected total
330 int total = 0;
331 for (const auto& [type, count] : globalCounts)
332 total += count;
333 CHECK(total == 958994);
334 }
335}
336
337// ===========================================================================
338// Per-entity index conversions -- Node (parameterized)
339// ===========================================================================
340TEST_CASE("NodeIndex: Global2Local round-trip on father nodes")
341{
343 for (DNDS::index iNode = 0; iNode < m->NumNode(); iNode++)
344 {
345 DNDS::index g = m->NodeIndexLocal2Global(iNode);
346 CHECK(g >= 0);
347 CHECK(m->NodeIndexGlobal2Local(g) == iNode);
348 }
349 })
350}
351
352TEST_CASE("NodeIndex: Global2Local round-trip on ghost nodes")
353{
354 if (g_mpi.size < 2)
355 return;
357 for (DNDS::index i = m->NumNode(); i < m->NumNode() + m->NumNodeGhost(); i++)
358 {
359 DNDS::index g = m->NodeIndexLocal2Global(i);
360 CHECK(g >= 0);
361 CHECK(m->NodeIndexGlobal2Local(g) == i);
362 }
363 })
364}
365
366TEST_CASE("NodeIndex: _NoSon round-trip on father nodes")
367{
369 for (DNDS::index iNode = 0; iNode < m->NumNode(); iNode++)
370 {
371 DNDS::index g = m->NodeIndexLocal2Global_NoSon(iNode);
372 CHECK(g >= 0);
373 CHECK(m->NodeIndexGlobal2Local_NoSon(g) == iNode);
374 }
375 })
376}
377
378TEST_CASE("NodeIndex: _NoSon returns negative for non-local global")
379{
380 if (g_mpi.size < 2)
381 return;
383 MPI_int other = (g_mpi.rank + 1) % g_mpi.size;
384 DNDS::index g = (*m->coords.father->pLGlobalMapping)(other, 0);
385 CHECK(m->NodeIndexGlobal2Local_NoSon(g) < 0);
386 })
387}
388
389TEST_CASE("NodeIndex: Global2Local returns negative for unknown global")
390{
392 DNDS::index bogus = m->coords.father->globalSize() + 999;
393 DNDS::index result = m->NodeIndexGlobal2Local(bogus);
394 CHECK(result < 0);
395 CHECK(result == -1 - bogus);
396 })
397}
398
399// ===========================================================================
400// Per-entity index conversions -- Cell (parameterized)
401// ===========================================================================
402TEST_CASE("CellIndex: Global2Local round-trip on father cells")
403{
405 for (DNDS::index iCell = 0; iCell < m->NumCell(); iCell++)
406 {
407 DNDS::index g = m->CellIndexLocal2Global(iCell);
408 CHECK(g >= 0);
409 CHECK(m->CellIndexGlobal2Local(g) == iCell);
410 }
411 })
412}
413
414TEST_CASE("CellIndex: Global2Local round-trip on ghost cells")
415{
416 if (g_mpi.size < 2)
417 return;
419 CHECK(m->NumCellGhost() > 0);
420 for (DNDS::index i = m->NumCell(); i < m->NumCell() + m->NumCellGhost(); i++)
421 {
422 DNDS::index g = m->CellIndexLocal2Global(i);
423 CHECK(g >= 0);
424 CHECK(m->CellIndexGlobal2Local(g) == i);
425 }
426 })
427}
428
429TEST_CASE("CellIndex: _NoSon round-trip on father cells")
430{
432 for (DNDS::index iCell = 0; iCell < m->NumCell(); iCell++)
433 {
434 DNDS::index g = m->CellIndexLocal2Global_NoSon(iCell);
435 CHECK(g >= 0);
436 CHECK(m->CellIndexGlobal2Local_NoSon(g) == iCell);
437 }
438 })
439}
440
441TEST_CASE("CellIndex: _NoSon returns negative for non-local global")
442{
443 if (g_mpi.size < 2)
444 return;
446 MPI_int other = (g_mpi.rank + 1) % g_mpi.size;
447 DNDS::index g = (*m->cell2node.father->pLGlobalMapping)(other, 0);
448 CHECK(m->CellIndexGlobal2Local_NoSon(g) < 0);
449 })
450}
451
452// ===========================================================================
453// Per-entity index conversions -- Bnd (parameterized)
454// ===========================================================================
455TEST_CASE("BndIndex: Global2Local round-trip on father bnds")
456{
458 for (DNDS::index iBnd = 0; iBnd < m->NumBnd(); iBnd++)
459 {
460 DNDS::index g = m->BndIndexLocal2Global(iBnd);
461 CHECK(g >= 0);
462 CHECK(m->BndIndexGlobal2Local(g) == iBnd);
463 }
464 })
465}
466
467TEST_CASE("BndIndex: _NoSon round-trip on father bnds")
468{
470 for (DNDS::index iBnd = 0; iBnd < m->NumBnd(); iBnd++)
471 {
472 DNDS::index g = m->BndIndexLocal2Global_NoSon(iBnd);
473 CHECK(g >= 0);
474 CHECK(m->BndIndexGlobal2Local_NoSon(g) == iBnd);
475 }
476 })
477}
478
479// ===========================================================================
480// UnInitIndex pass-through and special encodings (parameterized)
481// ===========================================================================
482TEST_CASE("UnInitIndex pass-through for all 12 conversion methods")
483{
485 CHECK(m->NodeIndexGlobal2Local(UnInitIndex) == UnInitIndex);
486 CHECK(m->CellIndexGlobal2Local(UnInitIndex) == UnInitIndex);
487 CHECK(m->BndIndexGlobal2Local(UnInitIndex) == UnInitIndex);
488
489 CHECK(m->NodeIndexLocal2Global(UnInitIndex) == UnInitIndex);
490 CHECK(m->CellIndexLocal2Global(UnInitIndex) == UnInitIndex);
491 CHECK(m->BndIndexLocal2Global(UnInitIndex) == UnInitIndex);
492
493 CHECK(m->NodeIndexLocal2Global_NoSon(UnInitIndex) == UnInitIndex);
494 CHECK(m->CellIndexLocal2Global_NoSon(UnInitIndex) == UnInitIndex);
495 CHECK(m->BndIndexLocal2Global_NoSon(UnInitIndex) == UnInitIndex);
496
497 CHECK(m->NodeIndexGlobal2Local_NoSon(UnInitIndex) == UnInitIndex);
498 CHECK(m->CellIndexGlobal2Local_NoSon(UnInitIndex) == UnInitIndex);
499 CHECK(m->BndIndexGlobal2Local_NoSon(UnInitIndex) == UnInitIndex);
500 })
501}
502
503TEST_CASE("Local2Global: negative local index decodes via -1-x encoding")
504{
506 DNDS::index fakeGlobal = 42;
507 DNDS::index neg = -1 - fakeGlobal;
508
509 CHECK(m->NodeIndexLocal2Global(neg) == fakeGlobal);
510 CHECK(m->CellIndexLocal2Global(neg) == fakeGlobal);
511 CHECK(m->BndIndexLocal2Global(neg) == fakeGlobal);
512 })
513}
514
515// ===========================================================================
516// Local index range validity (parameterized)
517// ===========================================================================
518TEST_CASE("AdjPrimary: cell2node local indices are in valid range")
519{
521 REQUIRE(m->adjPrimaryState == Adj_PointToLocal);
522 DNDS::index totalNodes = m->NumNode() + m->NumNodeGhost();
523 for (DNDS::index iC = 0; iC < m->NumCell(); iC++)
524 for (DNDS::rowsize j = 0; j < m->cell2node.RowSize(iC); j++)
525 {
526 DNDS::index iNode = m->cell2node(iC, j);
527 CHECK(iNode >= 0);
528 CHECK(iNode < totalNodes);
529 }
530 })
531}
532
533TEST_CASE("AdjPrimary: cell2cell local entries are valid or not-found")
534{
536 REQUIRE(m->adjPrimaryState == Adj_PointToLocal);
537 DNDS::index totalCells = m->NumCell() + m->NumCellGhost();
538 for (DNDS::index iC = 0; iC < m->NumCell(); iC++)
539 for (DNDS::rowsize j = 0; j < m->cell2cell.RowSize(iC); j++)
540 {
541 DNDS::index n = m->cell2cell(iC, j);
542 if (n >= 0)
543 CHECK(n < totalCells);
544 }
545 })
546}
547
548TEST_CASE("AdjPrimary: bnd2cell owner cell is a local father cell")
549{
551 REQUIRE(m->adjPrimaryState == Adj_PointToLocal);
552 for (DNDS::index ib = 0; ib < m->NumBnd(); ib++)
553 {
554 DNDS::index owner = m->bnd2cell(ib, 0);
555 CHECK(owner >= 0);
556 CHECK(owner < m->NumCell());
557 }
558 })
559}
560
561TEST_CASE("AdjPrimary: bnd2node local indices are in valid range")
562{
564 REQUIRE(m->adjPrimaryState == Adj_PointToLocal);
565 DNDS::index totalNodes = m->NumNode() + m->NumNodeGhost();
566 for (DNDS::index ib = 0; ib < m->NumBnd(); ib++)
567 for (DNDS::rowsize j = 0; j < m->bnd2node.RowSize(ib); j++)
568 {
569 DNDS::index iNode = m->bnd2node(ib, j);
570 CHECK(iNode >= 0);
571 CHECK(iNode < totalNodes);
572 }
573 })
574}
575
576// ===========================================================================
577// Adjacency state machine -- parameterized over all configs
578// ===========================================================================
579TEST_CASE("AdjPrimary: Local2Global then Global2Local is identity")
580{
582 REQUIRE(m->adjPrimaryState == Adj_PointToLocal);
583
584 auto c2nSnap = snapshotAdj(m->cell2node, m->NumCell());
585 auto c2cSnap = snapshotAdj(m->cell2cell, m->NumCell());
586
587 m->AdjLocal2GlobalPrimary();
588 CHECK(m->adjPrimaryState == Adj_PointToGlobal);
589
590 for (DNDS::index iC = 0; iC < m->NumCell(); iC++)
591 for (DNDS::rowsize j = 0; j < m->cell2node.RowSize(iC); j++)
592 {
593 DNDS::index gNode = m->cell2node(iC, j);
594 CHECK(gNode >= 0);
595 CHECK(gNode < m->coords.father->globalSize());
596 }
597
598 m->AdjGlobal2LocalPrimary();
599 CHECK(m->adjPrimaryState == Adj_PointToLocal);
600
601 for (DNDS::index iC = 0; iC < m->NumCell(); iC++)
602 {
603 REQUIRE(static_cast<DNDS::rowsize>(c2nSnap[iC].size()) == m->cell2node.RowSize(iC));
604 for (DNDS::rowsize j = 0; j < m->cell2node.RowSize(iC); j++)
605 CHECK(m->cell2node(iC, j) == c2nSnap[iC][j]);
606 REQUIRE(static_cast<DNDS::rowsize>(c2cSnap[iC].size()) == m->cell2cell.RowSize(iC));
607 for (DNDS::rowsize j = 0; j < m->cell2cell.RowSize(iC); j++)
608 CHECK(m->cell2cell(iC, j) == c2cSnap[iC][j]);
609 }
610 })
611}
612
613TEST_CASE("AdjPrimary: three consecutive round-trips are stable")
614{
616 auto snap = snapshotAdj(m->cell2node, m->NumCell());
617
618 for (int trip = 0; trip < 3; trip++)
619 {
620 m->AdjLocal2GlobalPrimary();
621 CHECK(m->adjPrimaryState == Adj_PointToGlobal);
622 m->AdjGlobal2LocalPrimary();
623 CHECK(m->adjPrimaryState == Adj_PointToLocal);
624 }
625
626 for (DNDS::index iC = 0; iC < m->NumCell(); iC++)
627 for (DNDS::rowsize j = 0; j < m->cell2node.RowSize(iC); j++)
628 CHECK(m->cell2node(iC, j) == snap[iC][j]);
629 })
630}
631
632TEST_CASE("AdjPrimaryForBnd: round-trip on cell2node only")
633{
635 REQUIRE(m->adjPrimaryState == Adj_PointToLocal);
636
637 auto snap = snapshotAdj(m->cell2node, m->NumCell());
638
639 m->AdjLocal2GlobalPrimaryForBnd();
640 CHECK(m->adjPrimaryState == Adj_PointToGlobal);
641
642 m->AdjGlobal2LocalPrimaryForBnd();
643 CHECK(m->adjPrimaryState == Adj_PointToLocal);
644
645 for (DNDS::index iC = 0; iC < m->NumCell(); iC++)
646 for (DNDS::rowsize j = 0; j < m->cell2node.RowSize(iC); j++)
647 CHECK(m->cell2node(iC, j) == snap[iC][j]);
648 })
649}
#define DNDS_assert(expr)
Debug-only assertion (compiled out when DNDS_NDEBUG is defined). Prints the expression + file/line + ...
Definition Errors.hpp:108
int main()
Definition dummy.cpp:3
Eigen::Vector3d tPoint
Definition Geometric.hpp:9
the host side operators are provided as implemented
DNDS_CONSTANT const index UnInitIndex
Sentinel "not initialised" index value (= INT64_MIN).
Definition Defines.hpp:176
int32_t rowsize
Row-width / per-row element-count type (signed 32-bit).
Definition Defines.hpp:109
int64_t index
Global row / DOF index type (signed 64-bit; handles multi-billion-cell meshes).
Definition Defines.hpp:107
std::shared_ptr< T > ssp
Shortened alias for std::shared_ptr used pervasively in DNDSR.
Definition Defines.hpp:138
int MPI_int
MPI counterpart type for MPI_int (= C int). Used for counts and ranks in MPI calls.
Definition MPI.hpp:43
Convenience bundle of a father, son, and attached ArrayTransformer.
Lightweight bundle of an MPI communicator and the calling rank's coordinates.
Definition MPI.hpp:215
int size
Number of ranks in comm (-1 until initialised).
Definition MPI.hpp:221
int rank
This rank's 0-based index within comm (-1 until initialised).
Definition MPI.hpp:219
MPI_Comm comm
The underlying MPI communicator handle.
Definition MPI.hpp:217
void setWorld()
Initialise the object to MPI_COMM_WORLD. Requires MPI_Init to have run.
Definition MPI.hpp:242
tPoint translation1
Periodic translation axis 1.
bool periodic
Has periodic boundaries.
int dim
Spatial dimension.
tPoint translation3
Periodic translation axis 3.
const char * file
Filename relative to data/mesh/.
tPoint translation2
Periodic translation axis 2.
const char * name
Human-readable label.
auto adj
CHECK(result.size()==3)
auto result
#define FOR_EACH_MESH_CONFIG(body)
auto res
Definition test_ODE.cpp:196
TEST_CASE("3D: VFV P2 HQM error < P1 on sinCos3D")
Eigen::Vector3d n(1.0, 0.0, 0.0)