21#define DOCTEST_CONFIG_IMPLEMENT
31#include <nlohmann/json.hpp>
39static constexpr int g_dim = 2;
40static constexpr int g_nv = 1;
48static std::string meshPath(
const std::string &name)
50 std::string f(__FILE__);
51 for (
int i = 0; i < 4; i++)
53 auto pos = f.rfind(
'/');
54 if (pos == std::string::npos)
56 if (pos != std::string::npos)
59 return f +
"/data/mesh/" + name;
63 const std::string &file,
bool periodic,
66 auto mesh = std::make_shared<UnstructuredMesh>(g_mpi, g_dim);
72 mesh->SetPeriodicGeometry(
73 tPoint{Lx, 0, 0}, zero, zero,
74 tPoint{0, Ly, 0}, zero, zero,
75 tPoint{0, 0, 0}, zero, zero);
78 reader.ReadFromCGNSSerial(meshPath(file));
79 reader.Deduplicate1to1Periodic();
80 reader.BuildCell2Cell();
84 reader.MeshPartitionCell2Cell(pOpt);
85 reader.PartitionReorderToMeshCell2Cell();
87 mesh->RecoverNode2CellAndNode2Bnd();
88 mesh->RecoverCell2CellAndBnd2Cell();
89 mesh->BuildGhostPrimary();
90 mesh->AdjGlobal2LocalPrimary();
93 for (
int ib = 0; ib < nBisect; ib++)
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();
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();
112 mesh->InterpolateFace();
113 mesh->AssertOnFaces();
131static const char *recMethodName(
RecMethod m)
163 auto vr = std::make_shared<tVR>(g_mpi, mesh);
166 nlohmann::ordered_json j;
167 defaultSettings.WriteIntoJson(j);
169 int order = recMethodOrder(method);
170 j[
"maxOrder"] =
order;
171 j[
"intOrder"] = std::max(
order + 2, 5);
172 j[
"cacheDiffBase"] =
true;
174 j[
"SORInstead"] =
false;
175 j[
"jacobiRelax"] = 1.0;
183 j[
"subs2ndOrder"] = 1;
188 j[
"subs2ndOrder"] = 0;
189 j[
"functionalSettings"][
"dirWeightScheme"] =
"HQM_OPT";
190 j[
"functionalSettings"][
"geomWeightScheme"] =
"HQM_SD";
191 j[
"functionalSettings"][
"geomWeightPower"] = 0.5;
192 j[
"functionalSettings"][
"geomWeightBias"] = 1.0;
196 j[
"subs2ndOrder"] = 0;
197 j[
"functionalSettings"][
"dirWeightScheme"] =
"Factorial";
198 j[
"functionalSettings"][
"geomWeightScheme"] =
"GWNone";
201 vr->parseSettings(j);
202 if (mesh->isPeriodic)
203 vr->SetPeriodicTransformations();
204 vr->ConstructMetrics();
205 vr->ConstructBaseAndWeight();
206 vr->ConstructRecCoeff();
215static tVR::TFBoundary<g_nv> makeDirichletBC(
const ScalarFunc &f)
220 return Eigen::Vector<DNDS::real, g_nv>{f(pPhy)};
224static tVR::TFBoundary<g_nv> g_zeroBC =
227{
return Eigen::Vector<DNDS::real, g_nv>::Zero(); };
238 const tVR::TFBoundary<g_nv> &bc,
243 auto mesh = vr->mesh;
250 for (
DNDS::index iCell = 0; iCell < mesh->NumCell(); iCell++)
252 auto qCell = vr->GetCellQuad(iCell);
253 Eigen::Vector<DNDS::real, g_nv> uc;
255 qCell.IntegrationSimple(
257 [&](
auto &vInc,
int iG)
259 vInc(0) = exactFunc(vr->GetCellQuadraturePPhys(iCell, iG)) *
260 vr->GetCellJacobiDet(iCell, iG);
262 u[iCell] = uc / vr->GetCellVol(iCell);
264 u.
trans.startPersistentPull();
265 u.
trans.waitPersistentPull();
272 vr->BuildUGrad(uGrad, 1);
273 vr->DoReconstruction2ndGrad<g_nv>(uGrad, u, bc, 1 );
274 uGrad.trans.startPersistentPull();
275 uGrad.trans.waitPersistentPull();
278 static const auto Seq012 = Eigen::seq(Eigen::fix<0>, Eigen::fix<g_dim - 1>);
280 for (
DNDS::index iCell = 0; iCell < mesh->NumCell(); iCell++)
282 auto qCell = vr->GetCellQuad(iCell);
284 qCell.IntegrationSimple(
288 tPoint pPhy = vr->GetCellQuadraturePPhys(iCell, iG);
290 (uGrad[iCell].transpose() *
291 (pPhy - vr->GetCellBary(iCell))(Seq012))(0);
293 vInc = std::abs(uRecVal - uExact) * vr->GetCellJacobiDet(iCell, iG);
299 return errGlobal / vr->GetGlobalVol();
305 vr->BuildURec(uRec, 1);
306 vr->BuildURec(uRecNew, 1);
309 for (
int iter = 0; iter < maxIters; iter++)
311 vr->DoReconstructionIter<g_nv>(
312 uRec, uRecNew, u, bc,
true);
316 for (
DNDS::index iCell = 0; iCell < mesh->NumCell(); iCell++)
317 incLocal += (uRecNew[iCell] - uRec[iCell]).array().square().sum();
320 incGlobal = std::sqrt(incGlobal / mesh->NumCellGlobal());
322 std::swap(uRec, uRecNew);
323 uRec.trans.startPersistentPull();
324 uRec.trans.waitPersistentPull();
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;
331 if (convTol > 0 && incGlobal < convTol)
333 if (printProgress && g_mpi.
rank == 0)
334 std::cout <<
" converged at iter " << iter + 1 << std::endl;
341 for (
DNDS::index iCell = 0; iCell < mesh->NumCell(); iCell++)
343 auto qCell = vr->GetCellQuad(iCell);
345 qCell.IntegrationSimple(
349 Eigen::VectorXd baseVal =
350 vr->GetIntPointDiffBaseValue(
351 iCell, -1, -1, iG, std::array<int, 1>{0}, 1) *
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);
361 return errGlobal / vr->GetGlobalVol();
382 MPI_Init(&argc, &argv);
386 std::cout <<
"=== Building meshes ===" << std::endl;
388 g_wall_mesh = buildMesh(
"Uniform_3x3_wall.cgns",
false, 0, 0, 0);
390 for (
int ib = 0; ib < 3; ib++)
393 std::cout <<
" IV10_10 bisect=" << ib << std::endl;
394 g_iv10[ib] = buildMesh(
"IV10_10.cgns",
true, 10, 10, ib);
396 for (
int ib = 0; ib < 3; ib++)
399 std::cout <<
" IV10U_10 bisect=" << ib << std::endl;
400 g_iv10u[ib] = buildMesh(
"IV10U_10.cgns",
true, 10, 10, ib);
404 std::cout <<
"=== Meshes built ===" << std::endl;
406 doctest::Context ctx;
407 ctx.applyCommandLine(argc, argv);
411 for (
auto &m : g_iv10)
413 for (
auto &m : g_iv10u)
426{
return std::sin(g_k * p[0]) * std::cos(g_k * p[1]); };
429{
return std::cos(g_k * p[0]) + std::cos(g_k * p[1]); };
435#define POLY_TEST(testName, method, polyFunc, polyDeg) \
436 TEST_CASE("Wall/" testName "/" #method) \
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; \
446 CHECK(err < 1e-12); \
500 {1.5599270188e-02, 3.4233789068e-03, 7.8943623278e-04},
false},
502 {4.6604402914e-02, 9.2961629825e-03, 1.5347312301e-03},
true},
504 {3.0528143687e-03, 2.3057099673e-04, 2.4367733525e-05},
true},
506 {1.9105219870e-03, 4.6701352192e-05, 1.4890868814e-06},
false},
508 {4.6604402914e-02, 9.2961629825e-03, 1.5347312301e-03},
false},
510 {1.8840503911e-03, 2.4995731251e-05, 8.3670061853e-07},
false},
514 {1.3027440876e-02, 3.8074751503e-03, 1.0853979656e-03},
false},
516 {1.2040354507e-02, 2.2707678072e-03, 4.6833420900e-04},
false},
518 {5.5617445849e-04, 5.9964034731e-05, 6.4337896956e-06},
false},
520 {1.5878119293e-04, 9.8836538778e-06, 4.3407055649e-07},
false},
524 {5.5612138851e-04, 1.6082305163e-05, 6.6188225747e-07},
false},
527static const int g_nPeriodicTests =
sizeof(g_periodicTests) /
sizeof(g_periodicTests[0]);
531 for (
int ti = 0; ti < g_nPeriodicTests; ti++)
533 auto &tc = g_periodicTests[ti];
534 std::string label = std::string(tc.meshName) +
"/" + tc.
funcName +
535 "/" + recMethodName(tc.method);
537 SUBCASE(label.c_str())
539 for (
int ib = 0; ib < 3; ib++)
542 auto mesh = tc.meshArray[ib];
543 auto vr = buildVR(mesh, tc.method);
544 bool print = (ib == 0 && tc.checkConvergence);
546 g_zeroBC, tc.maxIters, tc.convTol, print);
549 std::cout <<
"[" << label <<
" bis=" << ib
550 <<
"] err = " << std::scientific
551 << std::setprecision(10) <<
err << std::endl;
553 if (tc.golden[ib] != 0.0)
554 CHECK(
err == doctest::Approx(tc.golden[ib]).epsilon(1e-6));
573 vr->BuildURec(uRec, 1);
574 vr->BuildURec(uRecNew, 1);
576 for (
DNDS::index iCell = 0; iCell < vr->mesh->NumCell(); iCell++)
578 auto qCell = vr->GetCellQuad(iCell);
579 Eigen::Vector<DNDS::real, g_nv> uc;
581 qCell.IntegrationSimple(
583 [&](
auto &vInc,
int iG)
585 vInc(0) = sinCos(vr->GetCellQuadraturePPhys(iCell, iG)) *
586 vr->GetCellJacobiDet(iCell, iG);
588 u[iCell] = uc / vr->GetCellVol(iCell);
590 u.
trans.startPersistentPull();
591 u.
trans.waitPersistentPull();
595 for (
int iter = 0; iter < 200; iter++)
597 vr->DoReconstructionIter<g_nv>(uRec, uRecNew, u, g_zeroBC,
true);
600 for (
DNDS::index iCell = 0; iCell < vr->mesh->NumCell(); iCell++)
601 incLocal += (uRecNew[iCell] - uRec[iCell]).array().square().sum();
604 incGlobal = std::sqrt(incGlobal / vr->mesh->NumCellGlobal());
606 std::swap(uRec, uRecNew);
607 uRec.trans.startPersistentPull();
608 uRec.trans.waitPersistentPull();
611 if (incGlobal < 1e-14)
613 convergedAt = iter + 1;
619 std::cout <<
"[Convergence P2-HQM IV10] " << convergedAt
620 <<
" iters, final inc = " << std::scientific << lastInc << std::endl;
622 CHECK(convergedAt > 0);
623 CHECK(convergedAt < 200);
640static tVR::tLimitBatch<g_nv> identityFM(
641 const Eigen::Vector<DNDS::real, g_nv> &,
642 const Eigen::Vector<DNDS::real, g_nv> &,
644 const Eigen::Ref<tVR::tLimitBatch<g_nv>> &data)
655 const tVR::TFBoundary<g_nv> &bc,
662 auto mesh = vr->mesh;
669 for (
DNDS::index iCell = 0; iCell < mesh->NumCell(); iCell++)
671 auto qCell = vr->GetCellQuad(iCell);
672 Eigen::Vector<DNDS::real, g_nv> uc;
674 qCell.IntegrationSimple(
676 [&](
auto &vInc,
int iG)
678 vInc(0) = exactFunc(vr->GetCellQuadraturePPhys(iCell, iG)) *
679 vr->GetCellJacobiDet(iCell, iG);
681 u[iCell] = uc / vr->GetCellVol(iCell);
683 u.
trans.startPersistentPull();
684 u.
trans.waitPersistentPull();
688 vr->BuildURec(uRec, 1);
689 vr->BuildURec(uRecNew, 1);
690 vr->BuildURec(uRecBuf, 1);
692 for (
int iter = 0; iter < maxIters; iter++)
694 vr->DoReconstructionIter<g_nv>(uRec, uRecNew, u, bc,
true);
695 std::swap(uRec, uRecNew);
696 uRec.trans.startPersistentPull();
697 uRec.trans.waitPersistentPull();
700 for (
DNDS::index iCell = 0; iCell < mesh->NumCell(); iCell++)
701 incLocal += (uRecNew[iCell] - uRec[iCell]).array().square().sum();
704 incGlobal = std::sqrt(incGlobal / mesh->NumCellGlobal());
705 if (convTol > 0 && incGlobal < convTol)
712 vr->DoCalculateSmoothIndicator<g_nv, 1>(si, uRec, u, std::array<int, 1>{0});
713 si.trans.startPersistentPull();
714 si.trans.waitPersistentPull();
717 tVR::tFMEig<g_nv> fm = identityFM;
718 tVR::tFMEig<g_nv> fmi = identityFM;
720 if (limiterKind == 0)
721 vr->DoLimiterWBAP_C<g_nv>(u, uRec, uRecNew, uRecBuf, si, ifAll, fm, fmi,
true);
723 vr->DoLimiterWBAP_3<g_nv>(u, uRec, uRecNew, uRecBuf, si, ifAll, fm, fmi,
true);
726 auto &uRecLim = uRecNew;
730 for (
DNDS::index iCell = 0; iCell < mesh->NumCell(); iCell++)
732 auto qCell = vr->GetCellQuad(iCell);
734 qCell.IntegrationSimple(
738 Eigen::VectorXd baseVal =
739 vr->GetIntPointDiffBaseValue(
740 iCell, -1, -1, iG, std::array<int, 1>{0}, 1) *
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);
750 return errGlobal / vr->GetGlobalVol();
771 {6.6975633577e-02, 2.2647765241e-02, 9.2028443979e-03}},
776 {7.1333971937e-02, 2.6226392939e-02, 1.1194383457e-02}},
781 {3.5176301510e-02, 1.4925026476e-02, 6.9002847378e-03}},
786 {6.6488044323e-02, 2.2593150005e-02, 9.1963848358e-03}},
791 {7.1372867657e-02, 2.6697894435e-02, 1.1593632890e-02}},
794static const int g_nLimiterTests =
sizeof(g_limiterTests) /
sizeof(g_limiterTests[0]);
796TEST_CASE(
"Limiter procedure: reconstruction + smooth indicator + WBAP limiter")
798 for (
int ti = 0; ti < g_nLimiterTests; ti++)
800 auto &tc = g_limiterTests[ti];
801 std::string label = std::string(tc.meshName) +
"/" + tc.
funcName +
802 "/" + recMethodName(tc.method) +
"/" + tc.limName;
804 SUBCASE(label.c_str())
806 for (
int ib = 0; ib < 3; ib++)
809 auto mesh = tc.meshArray[ib];
810 auto vr = buildVR(mesh, tc.method);
812 vr, tc.method, tc.func, g_zeroBC,
813 200, 1e-14, tc.limiterKind, tc.ifAll,
false);
816 std::cout <<
"[" << label <<
" bis=" << ib
817 <<
"] limited err = " << std::scientific
818 << std::setprecision(10) <<
err << std::endl;
820 if (tc.golden[ib] != 0.0)
821 CHECK(
err == doctest::Approx(tc.golden[ib]).epsilon(1e-6));
Primary solver state container: an ArrayEigenMatrix pair with MPI-collective vector-space operations.
The VR class that provides any information needed in high-order CFV.
the host side operators are provided as implemented
DNDS_CONSTANT const real pi
π in double precision (matches DNDS_E_PI macro).
int64_t index
Global row / DOF index type (signed 64-bit; handles multi-billion-cell meshes).
std::shared_ptr< T > ssp
Shortened alias for std::shared_ptr used pervasively in DNDSR.
double real
Canonical floating-point scalar used throughout DNDSR (double precision).
DNDS_CONSTANT const real veryLargeReal
Catch-all upper bound ("practically infinity") for physical scalars.
const MPI_Datatype DNDS_MPI_REAL
MPI datatype matching real (= MPI_REAL8).
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.
int rank
This rank's 0-based index within comm (-1 until initialised).
MPI_Comm comm
The underlying MPI communicator handle.
void setWorld()
Initialise the object to MPI_COMM_WORLD. Requires MPI_Init to have run.
ssp< UnstructuredMesh > * meshArray
ssp< UnstructuredMesh > * meshArray
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