DNDSR 0.1.0.dev1+gcd065ad
Distributed Numeric Data Structure for CFV
Loading...
Searching...
No Matches
test_Reconstruction.cpp
Go to the documentation of this file.
1/**
2 * @file test_Reconstruction.cpp
3 * @brief Phase-0 regression tests for CFV Variational Reconstruction.
4 *
5 * Parameterized over [mesh+function, reconstruction_method].
6 *
7 * All iterative VR uses Jacobi iteration (SORInstead=false) so that
8 * golden values are deterministic across MPI partitionings (np=1,2,4).
9 *
10 * Error metric: L1 pointwise error at 6th-degree quadrature points,
11 * divided by domain volume (so the dimension matches the field function).
12 *
13 * Meshes:
14 * - Uniform_3x3_wall (9 quads, wall BC) -- polynomial exactness tests
15 * - IV10_10 (100 quads, periodic) + bisections -> 400, 1600 cells
16 * - IV10U_10 (322 tris, periodic) + bisections -> ~1288, ~5152 cells
17 *
18 * Golden values captured from commit c774b89 on dev/harry_refac1.
19 */
20
21#define DOCTEST_CONFIG_IMPLEMENT
22#include "doctest.h"
23
25
26#include <array>
27#include <cmath>
28#include <iostream>
29#include <iomanip>
30#include <functional>
31#include <nlohmann/json.hpp>
32#include <map>
33#include <string>
34#include <vector>
35
36using namespace DNDS;
37using namespace DNDS::Geom;
38
39static constexpr int g_dim = 2;
40static constexpr int g_nv = 1;
42
43static MPIInfo g_mpi;
44
45// ===================================================================
46// Mesh builder
47// ===================================================================
48static std::string meshPath(const std::string &name)
49{
50 std::string f(__FILE__);
51 for (int i = 0; i < 4; i++)
52 {
53 auto pos = f.rfind('/');
54 if (pos == std::string::npos)
55 pos = f.rfind('\\');
56 if (pos != std::string::npos)
57 f = f.substr(0, pos);
58 }
59 return f + "/data/mesh/" + name;
60}
61
62static ssp<UnstructuredMesh> buildMesh(
63 const std::string &file, bool periodic,
64 DNDS::real Lx, DNDS::real Ly, int nBisect = 0)
65{
66 auto mesh = std::make_shared<UnstructuredMesh>(g_mpi, g_dim);
67 UnstructuredMeshSerialRW reader(mesh, 0);
68
69 if (periodic)
70 {
71 tPoint zero{0, 0, 0};
72 mesh->SetPeriodicGeometry(
73 tPoint{Lx, 0, 0}, zero, zero,
74 tPoint{0, Ly, 0}, zero, zero,
75 tPoint{0, 0, 0}, zero, zero);
76 }
77
78 reader.ReadFromCGNSSerial(meshPath(file));
79 reader.Deduplicate1to1Periodic();
80 reader.BuildCell2Cell();
81
83 pOpt.metisSeed = 42;
84 reader.MeshPartitionCell2Cell(pOpt);
85 reader.PartitionReorderToMeshCell2Cell();
86
87 mesh->RecoverNode2CellAndNode2Bnd();
88 mesh->RecoverCell2CellAndBnd2Cell();
89 mesh->BuildGhostPrimary();
90 mesh->AdjGlobal2LocalPrimary();
91
92 // Bisect nBisect times via elevation+bisection
93 for (int ib = 0; ib < nBisect; ib++)
94 {
95 auto meshO2 = std::make_shared<UnstructuredMesh>(g_mpi, g_dim);
96 meshO2->BuildO2FromO1Elevation(*mesh);
97 meshO2->RecoverNode2CellAndNode2Bnd();
98 meshO2->RecoverCell2CellAndBnd2Cell();
99 meshO2->BuildGhostPrimary();
100 meshO2->AdjGlobal2LocalPrimary();
101
102 meshO2->AdjLocal2GlobalPrimary();
103 auto meshBis = std::make_shared<UnstructuredMesh>(g_mpi, g_dim);
104 meshBis->BuildBisectO1FormO2(*meshO2);
105 meshBis->RecoverNode2CellAndNode2Bnd();
106 meshBis->RecoverCell2CellAndBnd2Cell();
107 meshBis->BuildGhostPrimary();
108 meshBis->AdjGlobal2LocalPrimary();
109 mesh = meshBis;
110 }
111
112 mesh->InterpolateFace();
113 mesh->AssertOnFaces();
114 return mesh;
115}
116
117// ===================================================================
118// VR builder
119// ===================================================================
120enum class RecMethod
121{
122 GaussGreen, // explicit 2nd-order Gauss-Green gradient
123 VFV_P1_HQM, // iterative VR, maxOrder=1, HQM weights
124 VFV_P2_HQM, // iterative VR, maxOrder=2, HQM weights
125 VFV_P3_HQM, // iterative VR, maxOrder=3, HQM weights
129};
130
131static const char *recMethodName(RecMethod m)
132{
133 switch (m)
134 {
135 case RecMethod::GaussGreen: return "GG";
136 case RecMethod::VFV_P1_HQM: return "P1-HQM";
137 case RecMethod::VFV_P2_HQM: return "P2-HQM";
138 case RecMethod::VFV_P3_HQM: return "P3-HQM";
139 case RecMethod::VFV_P1_Default: return "P1-Def";
140 case RecMethod::VFV_P2_Default: return "P2-Def";
141 case RecMethod::VFV_P3_Default: return "P3-Def";
142 }
143 return "?";
144}
145
146static int recMethodOrder(RecMethod m)
147{
148 switch (m)
149 {
150 case RecMethod::GaussGreen: return 1;
152 case RecMethod::VFV_P1_Default: return 1;
154 case RecMethod::VFV_P2_Default: return 2;
156 case RecMethod::VFV_P3_Default: return 3;
157 }
158 return 1;
159}
160
161static ssp<tVR> buildVR(ssp<UnstructuredMesh> mesh, RecMethod method)
162{
163 auto vr = std::make_shared<tVR>(g_mpi, mesh);
164
165 CFV::VRSettings defaultSettings(g_dim);
166 nlohmann::ordered_json j;
167 defaultSettings.WriteIntoJson(j);
168
169 int order = recMethodOrder(method);
170 j["maxOrder"] = order;
171 j["intOrder"] = std::max(order + 2, 5);
172 j["cacheDiffBase"] = true;
173 // Force Jacobi iteration for deterministic results across MPI partitions
174 j["SORInstead"] = false;
175 j["jacobiRelax"] = 1.0;
176
177 bool isHQM = (method == RecMethod::VFV_P1_HQM ||
178 method == RecMethod::VFV_P2_HQM ||
179 method == RecMethod::VFV_P3_HQM);
180
181 if (method == RecMethod::GaussGreen)
182 {
183 j["subs2ndOrder"] = 1; // Gauss-Green
184 // weights don't matter for explicit GG
185 }
186 else if (isHQM)
187 {
188 j["subs2ndOrder"] = 0; // full VFV
189 j["functionalSettings"]["dirWeightScheme"] = "HQM_OPT";
190 j["functionalSettings"]["geomWeightScheme"] = "HQM_SD";
191 j["functionalSettings"]["geomWeightPower"] = 0.5;
192 j["functionalSettings"]["geomWeightBias"] = 1.0;
193 }
194 else
195 {
196 j["subs2ndOrder"] = 0; // full VFV
197 j["functionalSettings"]["dirWeightScheme"] = "Factorial";
198 j["functionalSettings"]["geomWeightScheme"] = "GWNone";
199 }
200
201 vr->parseSettings(j);
202 if (mesh->isPeriodic)
203 vr->SetPeriodicTransformations(); // identity for scalar
204 vr->ConstructMetrics();
205 vr->ConstructBaseAndWeight();
206 vr->ConstructRecCoeff();
207 return vr;
208}
209
210// ===================================================================
211// Boundary callback: Dirichlet (returns exact value for wall tests)
212// ===================================================================
213using ScalarFunc = std::function<DNDS::real(const tPoint &)>;
214
215static tVR::TFBoundary<g_nv> makeDirichletBC(const ScalarFunc &f)
216{
217 return [f](const auto &, const auto &, DNDS::index, DNDS::index, int,
218 const tPoint &, const tPoint &pPhy, t_index)
219 {
220 return Eigen::Vector<DNDS::real, g_nv>{f(pPhy)};
221 };
222}
223
224static tVR::TFBoundary<g_nv> g_zeroBC =
225 [](const auto &, const auto &, DNDS::index, DNDS::index, int,
226 const tPoint &, const tPoint &, t_index)
227{ return Eigen::Vector<DNDS::real, g_nv>::Zero(); };
228
229// ===================================================================
230// Core: run reconstruction and measure L1 error at 6th-degree quadrature
231// points, normalized by domain volume.
232// Returns the error and optionally prints iteration progress.
233// ===================================================================
234static DNDS::real runTest(
235 ssp<tVR> vr,
236 RecMethod method,
237 const ScalarFunc &exactFunc,
238 const tVR::TFBoundary<g_nv> &bc,
239 int maxIters,
240 DNDS::real convTol, // convergence threshold on increment; 0 = no early stop
241 bool printProgress)
242{
243 auto mesh = vr->mesh;
244
245 // --- Allocate arrays ---
247 vr->BuildUDof(u, 1);
248
249 // --- Set cell-averaged DOFs via quadrature ---
250 for (DNDS::index iCell = 0; iCell < mesh->NumCell(); iCell++)
251 {
252 auto qCell = vr->GetCellQuad(iCell);
253 Eigen::Vector<DNDS::real, g_nv> uc;
254 uc.setZero();
255 qCell.IntegrationSimple(
256 uc,
257 [&](auto &vInc, int iG)
258 {
259 vInc(0) = exactFunc(vr->GetCellQuadraturePPhys(iCell, iG)) *
260 vr->GetCellJacobiDet(iCell, iG);
261 });
262 u[iCell] = uc / vr->GetCellVol(iCell);
263 }
264 u.trans.startPersistentPull();
265 u.trans.waitPersistentPull();
266
267 // --- Reconstruct ---
268 if (method == RecMethod::GaussGreen)
269 {
270 // Explicit Gauss-Green: produces gradient, not polynomial coefficients
272 vr->BuildUGrad(uGrad, 1);
273 vr->DoReconstruction2ndGrad<g_nv>(uGrad, u, bc, 1 /*GG method*/);
274 uGrad.trans.startPersistentPull();
275 uGrad.trans.waitPersistentPull();
276
277 // Measure L1 error: u_rec(x) = u_mean + grad^T * (x - x_bary)
278 static const auto Seq012 = Eigen::seq(Eigen::fix<0>, Eigen::fix<g_dim - 1>);
279 DNDS::real errLocal = 0.0;
280 for (DNDS::index iCell = 0; iCell < mesh->NumCell(); iCell++)
281 {
282 auto qCell = vr->GetCellQuad(iCell);
283 DNDS::real errCell = 0.0;
284 qCell.IntegrationSimple(
285 errCell,
286 [&](DNDS::real &vInc, int iG)
287 {
288 tPoint pPhy = vr->GetCellQuadraturePPhys(iCell, iG);
289 DNDS::real uRecVal = u[iCell](0) +
290 (uGrad[iCell].transpose() *
291 (pPhy - vr->GetCellBary(iCell))(Seq012))(0);
292 DNDS::real uExact = exactFunc(pPhy);
293 vInc = std::abs(uRecVal - uExact) * vr->GetCellJacobiDet(iCell, iG);
294 });
295 errLocal += errCell;
296 }
297 DNDS::real errGlobal = 0.0;
298 MPI::Allreduce(&errLocal, &errGlobal, 1, DNDS_MPI_REAL, MPI_SUM, g_mpi.comm);
299 return errGlobal / vr->GetGlobalVol();
300 }
301 else
302 {
303 // Iterative VFV reconstruction
304 CFV::tURec<g_nv> uRec, uRecNew;
305 vr->BuildURec(uRec, 1);
306 vr->BuildURec(uRecNew, 1);
307
308 DNDS::real lastInc = veryLargeReal;
309 for (int iter = 0; iter < maxIters; iter++)
310 {
311 vr->DoReconstructionIter<g_nv>(
312 uRec, uRecNew, u, bc, /*putIntoNew=*/true);
313
314 // Compute increment norm
315 DNDS::real incLocal = 0.0;
316 for (DNDS::index iCell = 0; iCell < mesh->NumCell(); iCell++)
317 incLocal += (uRecNew[iCell] - uRec[iCell]).array().square().sum();
318 DNDS::real incGlobal = 0.0;
319 MPI::Allreduce(&incLocal, &incGlobal, 1, DNDS_MPI_REAL, MPI_SUM, g_mpi.comm);
320 incGlobal = std::sqrt(incGlobal / mesh->NumCellGlobal());
321
322 std::swap(uRec, uRecNew);
323 uRec.trans.startPersistentPull();
324 uRec.trans.waitPersistentPull();
325
326 lastInc = incGlobal;
327 if (printProgress && g_mpi.rank == 0 && (iter < 5 || (iter + 1) % 20 == 0))
328 std::cout << " iter " << iter + 1 << " inc = "
329 << std::scientific << incGlobal << std::endl;
330
331 if (convTol > 0 && incGlobal < convTol)
332 {
333 if (printProgress && g_mpi.rank == 0)
334 std::cout << " converged at iter " << iter + 1 << std::endl;
335 break;
336 }
337 }
338
339 // Measure L1 error at cell quadrature points (using VR's intOrder, which >= 6)
340 DNDS::real errLocal = 0.0;
341 for (DNDS::index iCell = 0; iCell < mesh->NumCell(); iCell++)
342 {
343 auto qCell = vr->GetCellQuad(iCell);
344 DNDS::real errCell = 0.0;
345 qCell.IntegrationSimple(
346 errCell,
347 [&](DNDS::real &vInc, int iG)
348 {
349 Eigen::VectorXd baseVal =
350 vr->GetIntPointDiffBaseValue(
351 iCell, -1, -1, iG, std::array<int, 1>{0}, 1) *
352 uRec[iCell];
353 DNDS::real uRecVal = baseVal(0) + u[iCell](0);
354 DNDS::real uExact = exactFunc(vr->GetCellQuadraturePPhys(iCell, iG));
355 vInc = std::abs(uRecVal - uExact) * vr->GetCellJacobiDet(iCell, iG);
356 });
357 errLocal += errCell;
358 }
359 DNDS::real errGlobal = 0.0;
360 MPI::Allreduce(&errLocal, &errGlobal, 1, DNDS_MPI_REAL, MPI_SUM, g_mpi.comm);
361 return errGlobal / vr->GetGlobalVol();
362 }
363}
364
365// ===================================================================
366// Prebuilt test data
367// ===================================================================
368
369// --- Wall mesh (polynomial exactness) ---
370static ssp<UnstructuredMesh> g_wall_mesh;
371
372// --- Periodic meshes for convergence (IV10 quad, IV10U tri) ---
373static ssp<UnstructuredMesh> g_iv10[3]; // bisect 0,1,2
374static ssp<UnstructuredMesh> g_iv10u[3]; // bisect 0,1,2
375
376// --- VR objects keyed by (mesh_ptr, method) ---
377// We build them on demand in the test runner to avoid combinatorial explosion
378// at startup. But for the wall mesh we prebuild a few.
379
380int main(int argc, char **argv)
381{
382 MPI_Init(&argc, &argv);
383 g_mpi.setWorld();
384
385 if (g_mpi.rank == 0)
386 std::cout << "=== Building meshes ===" << std::endl;
387
388 g_wall_mesh = buildMesh("Uniform_3x3_wall.cgns", false, 0, 0, 0);
389
390 for (int ib = 0; ib < 3; ib++)
391 {
392 if (g_mpi.rank == 0)
393 std::cout << " IV10_10 bisect=" << ib << std::endl;
394 g_iv10[ib] = buildMesh("IV10_10.cgns", true, 10, 10, ib);
395 }
396 for (int ib = 0; ib < 3; ib++)
397 {
398 if (g_mpi.rank == 0)
399 std::cout << " IV10U_10 bisect=" << ib << std::endl;
400 g_iv10u[ib] = buildMesh("IV10U_10.cgns", true, 10, 10, ib);
401 }
402
403 if (g_mpi.rank == 0)
404 std::cout << "=== Meshes built ===" << std::endl;
405
406 doctest::Context ctx;
407 ctx.applyCommandLine(argc, argv);
408 int res = ctx.run();
409
410 g_wall_mesh.reset();
411 for (auto &m : g_iv10)
412 m.reset();
413 for (auto &m : g_iv10u)
414 m.reset();
415 MPI_Finalize();
416 return res;
417}
418
419// ===================================================================
420// Test functions (periodic, L=10)
421// ===================================================================
422static const DNDS::real g_L = 10.0;
423static const DNDS::real g_k = 2.0 * pi / g_L;
424
425static ScalarFunc sinCos = [](const tPoint &p)
426{ return std::sin(g_k * p[0]) * std::cos(g_k * p[1]); };
427
428static ScalarFunc cosPlusCos = [](const tPoint &p)
429{ return std::cos(g_k * p[0]) + std::cos(g_k * p[1]); };
430
431// ===================================================================
432// POLYNOMIAL EXACTNESS on wall mesh (Uniform_3x3_wall, [-1,2]^2)
433// ===================================================================
434
435#define POLY_TEST(testName, method, polyFunc, polyDeg) \
436 TEST_CASE("Wall/" testName "/" #method) \
437 { \
438 ScalarFunc f = polyFunc; \
439 auto vr = buildVR(g_wall_mesh, RecMethod::method); \
440 auto bc = makeDirichletBC(f); \
441 DNDS::real err = runTest(vr, RecMethod::method, f, bc, 100, 1e-15, false); \
442 if (g_mpi.rank == 0) \
443 std::cout << "[Wall/" testName "/" #method "] err = " \
444 << std::scientific << err << std::endl; \
445 if (polyDeg == 0) \
446 CHECK(err < 1e-12); \
447 else \
448 CHECK(err < 10.0); /* golden-value regression checked separately */ \
449 }
450
451// Constant (degree 0): exact for all methods
452POLY_TEST("const", GaussGreen, [](const tPoint &) { return 1.0; }, 0)
453POLY_TEST("const", VFV_P1_HQM, [](const tPoint &) { return 1.0; }, 0)
454POLY_TEST("const", VFV_P2_HQM, [](const tPoint &) { return 1.0; }, 0)
455POLY_TEST("const", VFV_P3_HQM, [](const tPoint &) { return 1.0; }, 0)
456POLY_TEST("const", VFV_P1_Default, [](const tPoint &) { return 1.0; }, 0)
457
458// Linear (degree 1): exact for GG and p>=1
459POLY_TEST("linear", GaussGreen, [](const tPoint &p) { return p[0] + 2*p[1]; }, 1)
460POLY_TEST("linear", VFV_P1_HQM, [](const tPoint &p) { return p[0] + 2*p[1]; }, 1)
461POLY_TEST("linear", VFV_P2_HQM, [](const tPoint &p) { return p[0] + 2*p[1]; }, 1)
462POLY_TEST("linear", VFV_P3_HQM, [](const tPoint &p) { return p[0] + 2*p[1]; }, 1)
463POLY_TEST("linear", VFV_P1_Default, [](const tPoint &p) { return p[0] + 2*p[1]; }, 1)
464
465// Quadratic (degree 2): exact for p>=2
466POLY_TEST("quad", VFV_P2_HQM, [](const tPoint &p) { return p[0]*p[0] + p[1]*p[1]; }, 2)
467POLY_TEST("quad", VFV_P3_HQM, [](const tPoint &p) { return p[0]*p[0] + p[1]*p[1]; }, 2)
468POLY_TEST("quad", VFV_P2_Default, [](const tPoint &p) { return p[0]*p[0] + p[1]*p[1]; }, 2)
469
470// Cubic (degree 3): exact for p>=3
471POLY_TEST("cubic", VFV_P3_HQM, [](const tPoint &p) { return p[0]*p[0]*p[0] + p[0]*p[1]*p[1]; }, 3)
472POLY_TEST("cubic", VFV_P3_Default, [](const tPoint &p) { return p[0]*p[0]*p[0] + p[0]*p[1]*p[1]; }, 3)
473
474#undef POLY_TEST
475
476// ===================================================================
477// PERIODIC SMOOTH FUNCTIONS on IV10 (quad) and IV10U (tri) meshes
478// with convergence series (bisect 0,1,2).
479// Golden values -- to be filled after first acquisition run.
480// ===================================================================
481
483{
484 const char *meshName; // "IV10" or "IV10U"
485 ssp<UnstructuredMesh> *meshArray; // pointer to g_iv10 or g_iv10u
488 const char *funcName;
491 DNDS::real golden[3]; // golden L1/vol for bisect 0,1,2 (0 = not yet acquired)
492 bool checkConvergence; // whether to also check the iteration converges
493};
494
495// Golden values captured from commit c774b89 on dev/harry_refac1 (np=1).
496// Jacobi iteration ensures determinism across all np values.
497static PeriodicTestCase g_periodicTests[] = {
498 // IV10 (quad) + sin*cos
499 {"IV10", g_iv10, RecMethod::GaussGreen, sinCos, "sincos", 1, 0,
500 {1.5599270188e-02, 3.4233789068e-03, 7.8943623278e-04}, false},
501 {"IV10", g_iv10, RecMethod::VFV_P1_HQM, sinCos, "sincos", 200, 1e-14,
502 {4.6604402914e-02, 9.2961629825e-03, 1.5347312301e-03}, true},
503 {"IV10", g_iv10, RecMethod::VFV_P2_HQM, sinCos, "sincos", 200, 1e-14,
504 {3.0528143687e-03, 2.3057099673e-04, 2.4367733525e-05}, true},
505 {"IV10", g_iv10, RecMethod::VFV_P3_HQM, sinCos, "sincos", 200, 1e-14,
506 {1.9105219870e-03, 4.6701352192e-05, 1.4890868814e-06}, false},
507 {"IV10", g_iv10, RecMethod::VFV_P1_Default, sinCos, "sincos", 200, 1e-14,
508 {4.6604402914e-02, 9.2961629825e-03, 1.5347312301e-03}, false},
509 {"IV10", g_iv10, RecMethod::VFV_P3_Default, sinCos, "sincos", 200, 1e-14,
510 {1.8840503911e-03, 2.4995731251e-05, 8.3670061853e-07}, false},
511
512 // IV10U (tri) + sin*cos
513 {"IV10U", g_iv10u, RecMethod::GaussGreen, sinCos, "sincos", 1, 0,
514 {1.3027440876e-02, 3.8074751503e-03, 1.0853979656e-03}, false},
515 {"IV10U", g_iv10u, RecMethod::VFV_P1_HQM, sinCos, "sincos", 200, 1e-14,
516 {1.2040354507e-02, 2.2707678072e-03, 4.6833420900e-04}, false},
517 {"IV10U", g_iv10u, RecMethod::VFV_P2_HQM, sinCos, "sincos", 200, 1e-14,
518 {5.5617445849e-04, 5.9964034731e-05, 6.4337896956e-06}, false},
519 {"IV10U", g_iv10u, RecMethod::VFV_P3_HQM, sinCos, "sincos", 200, 1e-14,
520 {1.5878119293e-04, 9.8836538778e-06, 4.3407055649e-07}, false},
521
522 // IV10 (quad) + cos+cos
523 {"IV10", g_iv10, RecMethod::VFV_P3_HQM, cosPlusCos, "cos+cos", 200, 1e-14,
524 {5.5612138851e-04, 1.6082305163e-05, 6.6188225747e-07}, false},
525};
526
527static const int g_nPeriodicTests = sizeof(g_periodicTests) / sizeof(g_periodicTests[0]);
528
529TEST_CASE("Periodic reconstruction convergence series")
530{
531 for (int ti = 0; ti < g_nPeriodicTests; ti++)
532 {
533 auto &tc = g_periodicTests[ti];
534 std::string label = std::string(tc.meshName) + "/" + tc.funcName +
535 "/" + recMethodName(tc.method);
536
537 SUBCASE(label.c_str())
538 {
539 for (int ib = 0; ib < 3; ib++)
540 {
541 CAPTURE(ib);
542 auto mesh = tc.meshArray[ib];
543 auto vr = buildVR(mesh, tc.method);
544 bool print = (ib == 0 && tc.checkConvergence);
545 DNDS::real err = runTest(vr, tc.method, tc.func,
546 g_zeroBC, tc.maxIters, tc.convTol, print);
547
548 if (g_mpi.rank == 0)
549 std::cout << "[" << label << " bis=" << ib
550 << "] err = " << std::scientific
551 << std::setprecision(10) << err << std::endl;
552
553 if (tc.golden[ib] != 0.0)
554 CHECK(err == doctest::Approx(tc.golden[ib]).epsilon(1e-6));
555 else
556 CHECK(err >= 0.0); // acquisition: just check non-negative
557 }
558 }
559 }
560}
561
562// ===================================================================
563// Convergence check: selected VFV cases should converge
564// ===================================================================
565
566TEST_CASE("VFV P2 HQM converges on IV10 base mesh")
567{
568 auto vr = buildVR(g_iv10[0], RecMethod::VFV_P2_HQM);
569
571 vr->BuildUDof(u, 1);
572 CFV::tURec<g_nv> uRec, uRecNew;
573 vr->BuildURec(uRec, 1);
574 vr->BuildURec(uRecNew, 1);
575
576 for (DNDS::index iCell = 0; iCell < vr->mesh->NumCell(); iCell++)
577 {
578 auto qCell = vr->GetCellQuad(iCell);
579 Eigen::Vector<DNDS::real, g_nv> uc;
580 uc.setZero();
581 qCell.IntegrationSimple(
582 uc,
583 [&](auto &vInc, int iG)
584 {
585 vInc(0) = sinCos(vr->GetCellQuadraturePPhys(iCell, iG)) *
586 vr->GetCellJacobiDet(iCell, iG);
587 });
588 u[iCell] = uc / vr->GetCellVol(iCell);
589 }
590 u.trans.startPersistentPull();
591 u.trans.waitPersistentPull();
592
593 DNDS::real lastInc = veryLargeReal;
594 int convergedAt = 0;
595 for (int iter = 0; iter < 200; iter++)
596 {
597 vr->DoReconstructionIter<g_nv>(uRec, uRecNew, u, g_zeroBC, true);
598
599 DNDS::real incLocal = 0.0;
600 for (DNDS::index iCell = 0; iCell < vr->mesh->NumCell(); iCell++)
601 incLocal += (uRecNew[iCell] - uRec[iCell]).array().square().sum();
602 DNDS::real incGlobal = 0.0;
603 MPI::Allreduce(&incLocal, &incGlobal, 1, DNDS_MPI_REAL, MPI_SUM, g_mpi.comm);
604 incGlobal = std::sqrt(incGlobal / vr->mesh->NumCellGlobal());
605
606 std::swap(uRec, uRecNew);
607 uRec.trans.startPersistentPull();
608 uRec.trans.waitPersistentPull();
609
610 lastInc = incGlobal;
611 if (incGlobal < 1e-14)
612 {
613 convergedAt = iter + 1;
614 break;
615 }
616 }
617
618 if (g_mpi.rank == 0)
619 std::cout << "[Convergence P2-HQM IV10] " << convergedAt
620 << " iters, final inc = " << std::scientific << lastInc << std::endl;
621
622 CHECK(convergedAt > 0);
623 CHECK(convergedAt < 200);
624}
625
626// ===================================================================
627// LIMITER PROCEDURE TESTS
628//
629// After reconstruction, apply DoCalculateSmoothIndicator + DoLimiterWBAP_C
630// and measure the post-limiter L1 error. For scalar fields the
631// eigenvalue transform (FM/FMI) is identity.
632//
633// Golden value sentinel: 0.0 means "not yet acquired -- just print and
634// check non-negative".
635// ===================================================================
636
638
639/// Identity eigenvalue transform for nVarsFixed==1 (scalar).
640static tVR::tLimitBatch<g_nv> identityFM(
641 const Eigen::Vector<DNDS::real, g_nv> &,
642 const Eigen::Vector<DNDS::real, g_nv> &,
643 const tPoint &,
644 const Eigen::Ref<tVR::tLimitBatch<g_nv>> &data)
645{
646 return data;
647}
648
649/// Run reconstruction, then limiter, then measure L1 error.
650/// @param limiterKind 0 = WBAP_C, 1 = WBAP_3
651static DNDS::real runLimitedTest(
652 ssp<tVR> vr,
653 RecMethod method,
654 const ScalarFunc &exactFunc,
655 const tVR::TFBoundary<g_nv> &bc,
656 int maxIters,
657 DNDS::real convTol,
658 int limiterKind,
659 bool ifAll,
660 bool printProgress)
661{
662 auto mesh = vr->mesh;
663
664 // --- Allocate arrays ---
666 vr->BuildUDof(u, 1);
667
668 // --- Set cell-averaged DOFs via quadrature ---
669 for (DNDS::index iCell = 0; iCell < mesh->NumCell(); iCell++)
670 {
671 auto qCell = vr->GetCellQuad(iCell);
672 Eigen::Vector<DNDS::real, g_nv> uc;
673 uc.setZero();
674 qCell.IntegrationSimple(
675 uc,
676 [&](auto &vInc, int iG)
677 {
678 vInc(0) = exactFunc(vr->GetCellQuadraturePPhys(iCell, iG)) *
679 vr->GetCellJacobiDet(iCell, iG);
680 });
681 u[iCell] = uc / vr->GetCellVol(iCell);
682 }
683 u.trans.startPersistentPull();
684 u.trans.waitPersistentPull();
685
686 // --- Reconstruct (iterative VFV) ---
687 CFV::tURec<g_nv> uRec, uRecNew, uRecBuf;
688 vr->BuildURec(uRec, 1);
689 vr->BuildURec(uRecNew, 1);
690 vr->BuildURec(uRecBuf, 1);
691
692 for (int iter = 0; iter < maxIters; iter++)
693 {
694 vr->DoReconstructionIter<g_nv>(uRec, uRecNew, u, bc, true);
695 std::swap(uRec, uRecNew);
696 uRec.trans.startPersistentPull();
697 uRec.trans.waitPersistentPull();
698
699 DNDS::real incLocal = 0.0;
700 for (DNDS::index iCell = 0; iCell < mesh->NumCell(); iCell++)
701 incLocal += (uRecNew[iCell] - uRec[iCell]).array().square().sum();
702 DNDS::real incGlobal = 0.0;
703 MPI::Allreduce(&incLocal, &incGlobal, 1, DNDS_MPI_REAL, MPI_SUM, g_mpi.comm);
704 incGlobal = std::sqrt(incGlobal / mesh->NumCellGlobal());
705 if (convTol > 0 && incGlobal < convTol)
706 break;
707 }
708
709 // --- Smooth indicator ---
711 vr->BuildScalar(si);
712 vr->DoCalculateSmoothIndicator<g_nv, 1>(si, uRec, u, std::array<int, 1>{0});
713 si.trans.startPersistentPull();
714 si.trans.waitPersistentPull();
715
716 // --- Limiter ---
717 tVR::tFMEig<g_nv> fm = identityFM;
718 tVR::tFMEig<g_nv> fmi = identityFM;
719
720 if (limiterKind == 0)
721 vr->DoLimiterWBAP_C<g_nv>(u, uRec, uRecNew, uRecBuf, si, ifAll, fm, fmi, /*putIntoNew=*/true);
722 else
723 vr->DoLimiterWBAP_3<g_nv>(u, uRec, uRecNew, uRecBuf, si, ifAll, fm, fmi, /*putIntoNew=*/true);
724
725 // uRecNew now holds the limited reconstruction
726 auto &uRecLim = uRecNew;
727
728 // --- Measure L1 error ---
729 DNDS::real errLocal = 0.0;
730 for (DNDS::index iCell = 0; iCell < mesh->NumCell(); iCell++)
731 {
732 auto qCell = vr->GetCellQuad(iCell);
733 DNDS::real errCell = 0.0;
734 qCell.IntegrationSimple(
735 errCell,
736 [&](DNDS::real &vInc, int iG)
737 {
738 Eigen::VectorXd baseVal =
739 vr->GetIntPointDiffBaseValue(
740 iCell, -1, -1, iG, std::array<int, 1>{0}, 1) *
741 uRecLim[iCell];
742 DNDS::real uRecVal = baseVal(0) + u[iCell](0);
743 DNDS::real uExact = exactFunc(vr->GetCellQuadraturePPhys(iCell, iG));
744 vInc = std::abs(uRecVal - uExact) * vr->GetCellJacobiDet(iCell, iG);
745 });
746 errLocal += errCell;
747 }
748 DNDS::real errGlobal = 0.0;
749 MPI::Allreduce(&errLocal, &errGlobal, 1, DNDS_MPI_REAL, MPI_SUM, g_mpi.comm);
750 return errGlobal / vr->GetGlobalVol();
751}
752
754{
755 const char *meshName;
759 const char *funcName;
760 int limiterKind; // 0 = WBAP_C, 1 = WBAP_3
761 const char *limName;
762 bool ifAll; // limiter applied to all cells (no smooth-indicator skip)
763 DNDS::real golden[3]; // golden L1/vol for bisect 0,1,2
764 // 0.0 = not yet acquired (just print, CHECK >= 0)
765};
766
767static LimiterTestCase g_limiterTests[] = {
768 // IV10 (quad), P2-HQM, sincos, WBAP_C (ifAll=true so every cell is limited)
769 {"IV10", g_iv10, RecMethod::VFV_P2_HQM, sinCos, "sincos",
770 0, "CWBAP", true,
771 {6.6975633577e-02, 2.2647765241e-02, 9.2028443979e-03}},
772
773 // IV10 (quad), P3-HQM, sincos, WBAP_C (ifAll=true)
774 {"IV10", g_iv10, RecMethod::VFV_P3_HQM, sinCos, "sincos",
775 0, "CWBAP", true,
776 {7.1333971937e-02, 2.6226392939e-02, 1.1194383457e-02}},
777
778 // IV10U (tri), P2-HQM, sincos, WBAP_C (ifAll=true)
779 {"IV10U", g_iv10u, RecMethod::VFV_P2_HQM, sinCos, "sincos",
780 0, "CWBAP", true,
781 {3.5176301510e-02, 1.4925026476e-02, 6.9002847378e-03}},
782
783 // IV10 (quad), P2-HQM, sincos, WBAP_3 (ifAll=true)
784 {"IV10", g_iv10, RecMethod::VFV_P2_HQM, sinCos, "sincos",
785 1, "3WBAP", true,
786 {6.6488044323e-02, 2.2593150005e-02, 9.1963848358e-03}},
787
788 // IV10 (quad), P3-HQM, sincos, WBAP_3 (ifAll=true)
789 {"IV10", g_iv10, RecMethod::VFV_P3_HQM, sinCos, "sincos",
790 1, "3WBAP", true,
791 {7.1372867657e-02, 2.6697894435e-02, 1.1593632890e-02}},
792};
793
794static const int g_nLimiterTests = sizeof(g_limiterTests) / sizeof(g_limiterTests[0]);
795
796TEST_CASE("Limiter procedure: reconstruction + smooth indicator + WBAP limiter")
797{
798 for (int ti = 0; ti < g_nLimiterTests; ti++)
799 {
800 auto &tc = g_limiterTests[ti];
801 std::string label = std::string(tc.meshName) + "/" + tc.funcName +
802 "/" + recMethodName(tc.method) + "/" + tc.limName;
803
804 SUBCASE(label.c_str())
805 {
806 for (int ib = 0; ib < 3; ib++)
807 {
808 CAPTURE(ib);
809 auto mesh = tc.meshArray[ib];
810 auto vr = buildVR(mesh, tc.method);
811 DNDS::real err = runLimitedTest(
812 vr, tc.method, tc.func, g_zeroBC,
813 200, 1e-14, tc.limiterKind, tc.ifAll, false);
814
815 if (g_mpi.rank == 0)
816 std::cout << "[" << label << " bis=" << ib
817 << "] limited err = " << std::scientific
818 << std::setprecision(10) << err << std::endl;
819
820 if (tc.golden[ib] != 0.0)
821 CHECK(err == doctest::Approx(tc.golden[ib]).epsilon(1e-6));
822 else
823 CHECK(err >= 0.0);
824 }
825 }
826 }
827}
Primary solver state container: an ArrayEigenMatrix pair with MPI-collective vector-space operations.
Definition ArrayDOF.hpp:172
The VR class that provides any information needed in high-order CFV.
int main()
Definition dummy.cpp:3
int32_t t_index
Definition Geometric.hpp:6
Eigen::Vector3d tPoint
Definition Geometric.hpp:9
the host side operators are provided as implemented
DNDS_CONSTANT const real pi
π in double precision (matches DNDS_E_PI macro).
Definition Defines.hpp:199
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
double real
Canonical floating-point scalar used throughout DNDSR (double precision).
Definition Defines.hpp:105
DNDS_CONSTANT const real veryLargeReal
Catch-all upper bound ("practically infinity") for physical scalars.
Definition Defines.hpp:190
const MPI_Datatype DNDS_MPI_REAL
MPI datatype matching real (= MPI_REAL8).
Definition MPI.hpp:92
TTrans trans
Ghost-communication engine bound to father and son.
A means to translate nlohmann json into c++ primitive data types and back; and stores then during com...
Lightweight bundle of an MPI communicator and the calling rank's coordinates.
Definition MPI.hpp:215
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
ssp< UnstructuredMesh > * meshArray
ssp< UnstructuredMesh > * meshArray
real err
CHECK(result.size()==3)
double order
Definition test_ODE.cpp:257
auto res
Definition test_ODE.cpp:196
TEST_CASE("3D: VFV P2 HQM error < P1 on sinCos3D")
std::function< DNDS::real(const tPoint &)> ScalarFunc
#define POLY_TEST(testName, method, polyFunc, polyDeg)
CFV::VariationalReconstruction< g_dim > tVR