DNDSR 0.2.1
Distributed Numeric Data Structure for CFV
Loading...
Searching...
No Matches
EulerSolver.hpp
Go to the documentation of this file.
1/**
2 * @file EulerSolver.hpp
3 * @brief Top-level solver orchestrator for compressible Navier-Stokes / Euler simulations.
4 *
5 * Provides the EulerSolver class template which owns and coordinates all solver
6 * components: mesh infrastructure, variational reconstruction, the EulerEvaluator
7 * spatial discretization, DOF arrays, ODE integrators, linear solvers (LU-SGS / GMRES),
8 * and I/O subsystems.
9 *
10 * Responsibilities include:
11 * - JSON-based configuration loading, merging, and validation
12 * - Mesh reading (serial CGNS or distributed), partitioning, order elevation, and bisection
13 * - Solver initialization (evaluator, DOF allocation, restart loading)
14 * - Implicit time-marching loop (dual time stepping, CFL ramping, convergence monitoring)
15 * - VTK/HDF5/Tecplot data output and restart I/O
16 * - Time-averaging for unsteady statistics
17 *
18 * Supported model specializations (via EulerModel enum):
19 * NS, NS_2D, NS_SA, NS_SA_3D, NS_2EQ, NS_2EQ_3D, NS_3D
20 *
21 * @see EulerEvaluator.hpp Spatial discretization and flux evaluation
22 * @see Solver/ODE.hpp ODE integrator interface
23 * @see Solver/Linear.hpp GMRES and preconditioner interface
24 */
25#pragma once
26
27#include <iomanip>
28#include <functional>
29#include <tuple>
30#include <filesystem>
31#include <mutex>
32#include <future>
33
34// #ifndef __DNDS_REALLY_COMPILING__
35// #define __DNDS_REALLY_COMPILING__
36// #define __DNDS_REALLY_COMPILING__HEADER_ON__
37// #endif
41#include "DNDS/CsvLog.hpp"
42#include "DNDS/ObjectPool.hpp"
43#include "Solver/Linear.hpp"
44#include "Geom/Mesh/Mesh.hpp"
46#include "Gas.hpp"
47#include "EulerEvaluator.hpp"
48#include "EulerBC.hpp"
49// #ifdef __DNDS_REALLY_COMPILING__HEADER_ON__
50// #undef __DNDS_REALLY_COMPILING__
51// #endif
52
54
55#include "Solver/ODE.hpp"
56#include "Solver/Linear.hpp"
57
58namespace DNDS::Euler
59{
60
61 /**
62 * @brief Top-level solver orchestrator for compressible Navier-Stokes / Euler equations.
63 *
64 * Owns the complete solver pipeline: mesh, variational reconstruction, spatial
65 * evaluator, DOF arrays, Jacobian data, ODE integrator, linear solvers, and I/O.
66 * Provides the main time-marching loop (RunImplicitEuler) which drives steady-state
67 * or unsteady simulations with implicit dual-time stepping, CFL ramping, and
68 * convergence monitoring.
69 *
70 * @tparam model EulerModel enum selecting the equation set and spatial dimension
71 * (NS, NS_2D, NS_SA, NS_SA_3D, NS_2EQ, NS_2EQ_3D, NS_3D).
72 */
73 template <EulerModel model>
75 {
76 int nVars = getNVars(model); ///< Runtime number of conserved variables.
77
78 public:
79 typedef EulerEvaluator<model> TEval; ///< Evaluator type for this model.
80 static const int nVarsFixed = TEval::nVarsFixed; ///< Compile-time number of conserved variables.
81
82 static const int dim = TEval::dim; ///< Spatial dimension (2 or 3).
83 // static const int gdim = TEval::gdim;
84 static const int gDim = TEval::gDim; ///< Geometric dimension of the mesh.
85 static const int I4 = TEval::I4; ///< Energy equation index (= dim + 1).
86
87 typedef typename TEval::TU TU; ///< Conservative variable vector type.
88 typedef typename TEval::TDiffU TDiffU; ///< Gradient of conserved variables type.
89 typedef typename TEval::TJacobianU TJacobianU; ///< Flux Jacobian matrix type.
90 typedef typename TEval::TVec TVec; ///< Spatial vector type.
91 typedef typename TEval::TMat TMat; ///< Spatial matrix type.
92 typedef typename TEval::TDof TDof; ///< Cell-centered DOF array type.
93 typedef typename TEval::TRec TRec; ///< Reconstruction coefficient array type.
94 typedef typename TEval::TScalar TScalar; ///< Scalar reconstruction coefficient array type.
95 typedef typename TEval::TVFV TVFV; ///< Variational reconstruction type.
96 typedef typename TEval::TpVFV TpVFV; ///< Shared pointer to VFV type.
97
98 using tGMRES_u = Linear::GMRES_LeftPreconditioned<TDof>; ///< GMRES solver type for conservative DOFs.
99 using tGMRES_uRec = Linear::GMRES_LeftPreconditioned<TRec>; ///< GMRES solver type for reconstruction coefficients.
101
102 private:
103 MPIInfo mpi; ///< MPI communicator and rank information.
104 ssp<Geom::UnstructuredMesh> mesh, meshBnd; ///< Volume mesh and (optional) boundary surface mesh.
105 TpVFV vfv; // ! gDim -> 3 for intellisense ///< Variational reconstruction object.
106 ssp<Geom::UnstructuredMeshSerialRW> reader, readerBnd; ///< Mesh reader for volume and boundary meshes.
107 ssp<EulerEvaluator<model>> pEval; ///< Spatial evaluator instance.
108
109 ArrayDOFV<nVarsFixed> u, uIncBufODE, wAveraged, uAveraged; ///< DOF arrays: solution, ODE increment buffer, time-averaged fields.
110 ObjectPool<ArrayDOFV<nVarsFixed>> uPool; ///< Object pool for temporary DOF arrays (used by ODE integrators).
111 ArrayRECV<nVarsFixed> uRec, uRecLimited, uRecNew, uRecNew1, uRecOld, uRec1, uRecInc, uRecInc1, uRecB, uRecB1; ///< Reconstruction arrays (current, limited, new, old, increment, etc.).
112 JacobianDiagBlock<nVarsFixed> JD, JD1, JDTmp, JSource, JSource1, JSourceTmp; ///< Diagonal Jacobian blocks for implicit methods.
113 ssp<JacobianLocalLU<nVarsFixed>> JLocalLU; ///< Local LU factorization for direct preconditioner.
114 ArrayDOFV<1> alphaPP, alphaPP1, betaPP, betaPP1, alphaPP_tmp, dTauTmp; ///< Positivity-preserving limiter scalars and time-step buffer.
115
116 int nOUTS = {-1}; ///< Number of output scalars per cell in volume output.
117 int nOUTSPoint{-1}; ///< Number of output scalars per node in point output.
118 int nOUTSBnd{-1}; ///< Number of output scalars per face in boundary output.
119 // rho u v w p T M ifUseLimiter RHS
120 ssp<ArrayEigenVector<Eigen::Dynamic>> outDist; ///< Distributed cell output array.
121 ssp<ArrayEigenVector<Eigen::Dynamic>> outSerial; ///< Serial (gathered) cell output array.
122 ArrayTransformerType<ArrayEigenVector<Eigen::Dynamic>>::Type outDist2SerialTrans; ///< Transformer for distributed-to-serial cell output.
123
124 ssp<ArrayEigenVector<Eigen::Dynamic>> outDistPoint; ///< Distributed node output array.
125 ssp<ArrayEigenVector<Eigen::Dynamic>> outGhostPoint; ///< Ghost-node output array.
126 ssp<ArrayEigenVector<Eigen::Dynamic>> outSerialPoint; ///< Serial (gathered) node output array.
127 ArrayTransformerType<ArrayEigenVector<Eigen::Dynamic>>::Type outDist2SerialTransPoint; ///< Transformer for distributed-to-serial node output.
128 ArrayPair<ArrayEigenVector<Eigen::Dynamic>> outDistPointPair; ///< Array pair for async node output.
129 static const int maxOutFutures{3}; ///< Maximum number of concurrent async output futures.
130 std::mutex outArraysMutex; ///< Mutex protecting output arrays during async writes.
131 std::array<std::future<void>, maxOutFutures> outFuture; ///< Futures for async volume output. Mind the order, relies on the arrays and the mutex.
132
133 ssp<ArrayEigenVector<Eigen::Dynamic>> outDistBnd; ///< Distributed boundary output array.
134 // ssp<ArrayEigenVector<Eigen::Dynamic>> outGhostBnd;
135 ssp<ArrayEigenVector<Eigen::Dynamic>> outSerialBnd; ///< Serial (gathered) boundary output array.
136 ArrayTransformerType<ArrayEigenVector<Eigen::Dynamic>>::Type outDist2SerialTransBnd; ///< Transformer for distributed-to-serial boundary output.
137 // ArrayPair<ArrayEigenVector<Eigen::Dynamic>> outDistBndPair;
138 std::mutex outBndArraysMutex; ///< Mutex protecting boundary output arrays during async writes.
139 std::array<std::future<void>, maxOutFutures> outBndFuture; ///< Futures for async boundary output. Mind the order, relies on the arrays and the mutex.
140 std::future<void> outSeqFuture; ///< Future for sequential (non-parallel) output operations.
141
142 // std::vector<uint32_t> ifUseLimiter;
143 CFV::tScalarPair ifUseLimiter; ///< Per-cell flag indicating whether the limiter was active.
144
145 ssp<BoundaryHandler<model>> pBCHandler; ///< Boundary condition handler (shared with evaluator).
146
147 public:
148 nlohmann::ordered_json gSetting; ///< Full JSON configuration (read from file, may be modified at runtime).
149 std::string output_stamp = ""; ///< Unique stamp appended to output filenames for this run.
150
151 /**
152 * @brief Complete solver configuration, serializable to/from JSON.
153 *
154 * Aggregates all sub-configurations: time marching, reconstruction, output,
155 * CFL control, convergence, data I/O, boundary definitions, limiters, linear
156 * solver, restart state, time averaging, evaluator settings, VFV settings,
157 * and boundary condition definitions. Each sub-struct uses DNDS_DECLARE_CONFIG
158 * for automatic JSON serialization.
159 */
161 {
162
163 /**
164 * @brief Time marching control parameters.
165 *
166 * Controls physical/pseudo time step sizes, ODE integrator selection,
167 * positivity-preserving options, and dt ramping strategies.
168 */
170 {
173 int nTimeStep = 1000000;
174 bool steadyQuit = false;
175 bool useRestart = false;
176 bool useImplicitPP = false;
177 int rhsFPPMode = 0;
181 int odeCode = 0;
187 nlohmann::ordered_json odeSettingsExtra;
188 bool partitionMeshOnly = false;
192 bool useDtPPLimit = false;
196 {
197 // clang-format off
198 DNDS_FIELD(dtImplicit, "Max implicit time step; 1e100 for steady",
200 DNDS_FIELD(dtImplicitMin, "Minimum implicit time step",
202 DNDS_FIELD(nTimeStep, "Max number of time steps",
204 DNDS_FIELD(steadyQuit, "Quit on steady convergence");
205 DNDS_FIELD(useRestart, "Initialize from restart file");
206 DNDS_FIELD(useImplicitPP, "Use implicit positivity-preserving");
207 DNDS_FIELD(rhsFPPMode, "RHS flux-positivity-preserving mode");
208 DNDS_FIELD(rhsFPPScale, "RHS FPP scaling factor");
209 DNDS_FIELD(rhsFPPRelax, "RHS FPP relaxation factor",
210 DNDS::Config::range(0.0, 1.0));
211 DNDS_FIELD(incrementPPRelax, "Increment PP relaxation factor",
212 DNDS::Config::range(0.0, 1.0));
213 DNDS_FIELD(odeCode, "ODE integrator code");
214 DNDS_FIELD(tEnd, "End time for unsteady simulation");
215 DNDS_FIELD(odeSetting1, "ODE parameter 1");
216 DNDS_FIELD(odeSetting2, "ODE parameter 2");
217 DNDS_FIELD(odeSetting3, "ODE parameter 3");
218 DNDS_FIELD(odeSetting4, "ODE parameter 4");
219 config.field_json(&T::odeSettingsExtra, "odeSettingsExtra", "Extra ODE integrator settings (opaque JSON)");
220 DNDS_FIELD(partitionMeshOnly, "Only partition mesh, then exit");
221 DNDS_FIELD(dtIncreaseLimit, "Max factor for dt increase per step",
223 DNDS_FIELD(dtIncreaseAfterCount,"Steps before dt increase allowed",
225 DNDS_FIELD(dtCFLLimitScale, "CFL-based dt limit scale");
226 DNDS_FIELD(useDtPPLimit, "Use positivity-preserving dt limiter");
227 DNDS_FIELD(dtPPLimitRelax, "PP dt limiter relaxation",
228 DNDS::Config::range(0.0, 1.0));
229 DNDS_FIELD(dtPPLimitScale, "PP dt limiter scale",
231 // clang-format on
232 }
234 {
235 return odeCode == 401 || (odeCode >= 411 && odeCode <= 413);
236 }
238
239 /**
240 * @brief Implicit reconstruction control parameters.
241 *
242 * Controls the iterative reconstruction solve within each time step:
243 * explicit vs implicit reconstruction, linear solver (SOR or GMRES/PCG),
244 * convergence thresholds, and gradient zeroing strategies.
245 */
247 {
248 bool useExplicit = false;
250 bool zeroGrads = false;
251 int recLinearScheme = 0; // 0 for original SOR, 1 for GMRES
252 int nGmresSpace = 5;
253 int nGmresIter = 10;
262 bool storeRecInc = false;
263 bool dampRecIncDTau = false;
267 {
268 // clang-format off
269 DNDS_FIELD(useExplicit, "Use explicit reconstruction (no implicit)");
270 DNDS_FIELD(nInternalRecStep, "Number of internal reconstruction sub-steps",
272 DNDS_FIELD(zeroGrads, "Zero gradients before reconstruction");
273 DNDS_FIELD(recLinearScheme, "Reconstruction linear solver: 0=SOR, 1=GMRES");
274 DNDS_FIELD(nGmresSpace, "GMRES Krylov subspace size for reconstruction",
276 DNDS_FIELD(nGmresIter, "GMRES iterations for reconstruction",
278 DNDS_FIELD(gmresRecScale, "GMRES reconstruction scaling mode");
279 DNDS_FIELD(fpcgResetScheme, "FPCG reset scheme");
280 DNDS_FIELD(fpcgResetThres, "FPCG reset threshold",
281 DNDS::Config::range(0.0, 1.0));
282 DNDS_FIELD(fpcgResetReport, "FPCG reset report frequency");
283 DNDS_FIELD(fpcgMaxPHistory, "FPCG max preconditioner history",
285 DNDS_FIELD(recThreshold, "Reconstruction convergence threshold",
287 DNDS_FIELD(nRecConsolCheck, "Console output interval for reconstruction",
289 DNDS_FIELD(nRecMultiplyForZeroedGrad, "Rec iteration multiplier when grad is zeroed",
291 DNDS_FIELD(storeRecInc, "Store reconstruction increment");
292 DNDS_FIELD(dampRecIncDTau, "Damp reconstruction increment by dTau");
293 DNDS_FIELD(zeroRecForSteps, "Zero reconstruction for N outer steps",
295 DNDS_FIELD(zeroRecForStepsInternal, "Zero reconstruction for N internal steps",
297 // clang-format on
298 }
300
301 /**
302 * @brief Output control parameters.
303 *
304 * Controls console output frequency and formatting, data file output
305 * intervals (VTK/HDF5), restart checkpoint intervals, time-average output,
306 * and logging precision.
307 */
309 {
312 int consoleOutputMode = 0; // 0 for basic, 1 for wall force out
315 std::vector<std::string> consoleMainOutputFormat{
316 "=== Step {termBold}[{step:4d}] ",
317 "res {termBold}{termRed}{resRel:.3e}{termReset} ",
318 "t,dT,dTaumin,CFL,nFix {termGreen}[{tSimu:.3e},{curDtImplicit:.3e},{curDtMin:.3e},{CFLNow:.3e},[alphaInc({nLimInc},{alphaMinInc}), betaRec({nLimBeta},{minBeta}), alphaRes({nLimAlpha},{minAlpha})]]{termReset} ",
319 "Time[{telapsed:.3f}] recTime[{trec:.3f}] rhsTime[{trhs:.3f}] commTime[{tcomm:.3f}] limTime[{tLim:.3f}] limTimeA[{tLimiterA:.3f}] limTimeB[{tLimiterB:.3f}]"};
320 std::vector<std::string> consoleMainOutputFormatInternal{
321 "\t Internal === Step [{step:4d},{iStep:2d},{iter:4d}] ",
322 "res {termRed}{resRel:.3e}{termReset} ",
323 "t,dT,dTaumin,CFL,nFix {termGreen}[{tSimu:.3e},{curDtImplicit:.3e},{curDtMin:.3e},{CFLNow:.3e},[alphaInc({nLimInc},{alphaMinInc:.3g}), betaRec({nLimBeta},{minBeta:.3g}), alphaRes({nLimAlpha},{minAlpha:.3g})]]{termReset} ",
324 "Time[{telapsedM:.3f}] recTime[{trecM:.3f}] rhsTime[{trhsM:.3f}] commTime[{tcommM:.3f}] limTime[{tLimM:.3f}] limTimeA[{tLimiterA:.3f}] limTimeB[{tLimiterB:.3f}]"};
325 std::vector<std::string> logfileOutputTitles{
326 "step", "iStep", "iterAll", "iter", "tSimu",
327 "res", "curDtImplicit", "curDtMin", "CFLNow",
328 "nLimInc", "alphaMinInc",
329 "nLimBeta", "minBeta",
330 "nLimAlpha", "minAlpha",
331 "tWall", "telapsed", "trec", "trhs", "tcomm", "tLim", "tLimiterA", "tLimiterB",
332 "fluxWall", "CL", "CD", "AoA"};
334 bool dataOutAtInit = false;
335 bool restartOutAtInit = false;
336 int nDataOut = 10000;
337 int nDataOutC = 50;
338 int nDataOutInternal = 10000;
340 int nRestartOut = INT_MAX;
341 int nRestartOutC = INT_MAX;
342 int nRestartOutInternal = INT_MAX;
343 int nRestartOutCInternal = INT_MAX;
344 int nTimeAverageOut = INT_MAX;
345 int nTimeAverageOutC = INT_MAX;
348 bool useCollectiveTimer = false;
349
351 {
352 // clang-format off
353 DNDS_FIELD(nConsoleCheck, "Console output interval (outer steps)",
355 DNDS_FIELD(nConsoleCheckInternal, "Console output interval (internal steps)",
357 DNDS_FIELD(consoleOutputMode, "Console output mode: 0=basic, 1=wall force");
358 DNDS_FIELD(consoleOutputEveryFix, "Force console output every N fix iterations");
359 // nPrecisionConsole intentionally not registered: was not serialized
360 // in the original macro and existing config files lack this key.
361 DNDS_FIELD(consoleMainOutputFormat, "Format strings for console output lines");
362 DNDS_FIELD(consoleMainOutputFormatInternal, "Format strings for internal-step console output");
363 DNDS_FIELD(logfileOutputTitles, "Column titles for CSV log file");
364 DNDS_FIELD(nPrecisionLog, "Log file floating-point precision",
366 DNDS_FIELD(dataOutAtInit, "Output data at initialization");
367 DNDS_FIELD(restartOutAtInit, "Output restart at initialization");
368 DNDS_FIELD(nDataOut, "Data output interval (outer steps)",
370 DNDS_FIELD(nDataOutC, "Data output interval (convergence steps)",
372 DNDS_FIELD(nDataOutInternal, "Data output interval (internal steps)",
374 DNDS_FIELD(nDataOutCInternal, "Data output interval (internal convergence)");
375 DNDS_FIELD(nRestartOut, "Restart output interval (outer steps)");
376 DNDS_FIELD(nRestartOutC, "Restart output interval (convergence steps)");
377 DNDS_FIELD(nRestartOutInternal, "Restart output interval (internal steps)");
378 DNDS_FIELD(nRestartOutCInternal, "Restart output interval (internal convergence)");
379 DNDS_FIELD(nTimeAverageOut, "Time-average output interval (outer steps)");
380 DNDS_FIELD(nTimeAverageOutC, "Time-average output interval (convergence steps)");
381 DNDS_FIELD(tDataOut, "Output data at simulation time interval");
382 DNDS_FIELD(lazyCoverDataOutput, "Overwrite previous data output files");
383 DNDS_FIELD(useCollectiveTimer, "Use collective MPI timer for profiling");
384 // clang-format on
385 }
387
388 /**
389 * @brief Implicit CFL number control parameters.
390 *
391 * Controls the CFL number, local vs global time stepping, CFL ramping
392 * schedule, local dTau smoothing, and RANS under-relaxation.
393 */
395 {
396 real CFL = 10;
397 int nForceLocalStartStep = INT_MAX;
398 int nCFLRampStart = INT_MAX;
399 int nCFLRampLength = INT_MAX;
401 bool useLocalDt = true;
402 int nSmoothDTau = 0;
405 {
406 // clang-format off
407 DNDS_FIELD(CFL, "CFL number for implicit time stepping",
409 DNDS_FIELD(nForceLocalStartStep, "Step to force local time stepping",
411 DNDS_FIELD(nCFLRampStart, "Step to begin CFL ramping",
413 DNDS_FIELD(nCFLRampLength, "Number of steps for CFL ramp",
415 DNDS_FIELD(CFLRampEnd, "CFL value at end of ramp");
416 DNDS_FIELD(useLocalDt, "Use local (vs uniform) time step");
417 DNDS_FIELD(nSmoothDTau, "Smoothing passes for local dTau",
419 DNDS_FIELD(RANSRelax, "RANS equation under-relaxation factor",
420 DNDS::Config::range(0.0, 1.0));
421 // clang-format on
422 }
424
425 /**
426 * @brief Convergence monitoring parameters.
427 *
428 * Controls inner-loop iteration counts, convergence thresholds,
429 * residual norm type, and CL-driver (lift-coefficient) adaptation.
430 */
432 {
439 int resBaseType = 0; // 1 to use rhs as base in unsteady
440 int mergeMultiResidual = 0; // for HM3, merge stage residuals
441 int normOrd = 1; // use LN norm
442 bool useVolWiseResidual = false;
443 bool useCLDriver = false;
445 {
446 // clang-format off
447 DNDS_FIELD(nTimeStepInternal, "Max internal iterations per time step (0 = unlimited)",
449 DNDS_FIELD(nTimeStepInternalMin, "Min internal iterations per time step",
451 DNDS_FIELD(nAnchorUpdate, "Anchor update frequency",
453 DNDS_FIELD(nAnchorUpdateStart, "Step to start anchor updates",
455 DNDS_FIELD(rhsThresholdInternal, "RHS convergence threshold for internal loop",
457 DNDS_FIELD(res_base, "Residual base value (0=auto)");
458 DNDS_FIELD(resBaseType, "Residual base type: 0=default, 1=rhs-based for unsteady");
459 DNDS_FIELD(mergeMultiResidual, "Merge stage residuals for multi-stage: 0=off");
460 DNDS_FIELD(normOrd, "Residual norm order (1=L1, 2=L2, etc.)",
462 DNDS_FIELD(useVolWiseResidual, "Volume-weighted residual");
463 DNDS_FIELD(useCLDriver, "Enable CL-driven AoA adaptation");
464 // clang-format on
465 }
467
468 /**
469 * @brief Data I/O control parameters.
470 *
471 * Controls mesh file paths and formats, mesh preprocessing (elevation,
472 * bisection, wall distance), output formats (Tecplot, VTK, HDF5),
473 * restart serialization, and coordinate transformations.
474 */
476 {
477 bool uniqueStamps = true;
480
481 int meshElevation = 0; // 0 = noOp, 1 = O1->O2
482 int meshElevationInternalSmoother = 0; // 0 = local interpolation, 1 = coupled
483 int meshElevationIter = 1000; // -1 to use handle all the nodes
490 int meshElevationBoundaryMode = 0; // 0: only wall bc; 1: wall + invis vall
491
493 int meshReorderCells = 0; // 0: natural; 1: reorder
496
497 int meshFormat = 0;
499 std::string meshFile = "data/mesh/NACA0012_WIDE_H3.cgns";
500 std::string outPltName = "data/out/debugData_";
501 std::string outLogName = "";
502 std::string outRestartName = "";
503
504 int outPltMode = 0; // 0 = serial, 1 = dist plt
505 int readMeshMode = 0; // 0 = serial cgns, 1 = dist json
507 bool outPltVTKFormat = false;
508 bool outPltVTKHDFFormat = false;
509 bool outAtPointData = true;
510 bool outAtCellData = true;
512 std::string vtuFloatEncodeMode = "binary";
513 int hdfChunkSize = 32768;
515 bool hdfCollOnData = false;
516 bool hdfCollOnMeta = true;
517 bool outVolumeData = true;
518 bool outBndData = false;
519
520 std::vector<std::string> outCellScalarNames{};
521
522 bool serializerSaveURec = false;
523
525
526 int rectifyNearPlane = 0; // 1: x 2: y 4: z
528
531 std::string meshPartitionedReaderType = "JSON";
532
533 const std::string &getOutLogName()
534 {
535 return outLogName.empty() ? outPltName : outLogName;
536 }
537
538 const std::string &getOutRestartName()
539 {
540 return outRestartName.empty() ? outPltName : outRestartName;
541 }
542
544 {
545 // clang-format off
546 DNDS_FIELD(uniqueStamps, "Use unique output stamps per run");
547 DNDS_FIELD(meshRotZ, "Mesh rotation around Z axis (degrees)");
548 DNDS_FIELD(meshScale, "Mesh coordinate scaling factor",
550 DNDS_FIELD(meshElevation, "Mesh order elevation: 0=none, 1=O1->O2");
551 DNDS_FIELD(meshElevationInternalSmoother,"Elevation smoother: 0=local, 1=coupled");
552 DNDS_FIELD(meshElevationIter, "Elevation smoothing iterations");
553 DNDS_FIELD(meshElevationNSearch, "Elevation RBF neighbor search count",
555 DNDS_FIELD(meshElevationRBFRadius, "Elevation RBF support radius",
557 DNDS_FIELD(meshElevationRBFPower, "Elevation RBF power parameter");
558 DNDS_FIELD(meshElevationRBFKernel, "Elevation RBF kernel type");
559 DNDS_FIELD(meshElevationMaxIncludedAngle,"Max included angle for elevation (degrees)");
560 DNDS_FIELD(meshElevationRefDWall, "Reference wall distance for elevation",
562 DNDS_FIELD(meshElevationBoundaryMode, "Elevation boundary: 0=wall only, 1=wall+inviscid");
563 DNDS_FIELD(meshDirectBisect, "Mesh bisection passes",
565 DNDS_FIELD(meshReorderCells, "Reorder cells: 0=natural, 1=reorder");
566 DNDS_FIELD(meshBuildWallDist, "Build wall distance field: 0=no, 1=yes");
567 DNDS_FIELD(meshWallDistOptions, "Wall distance computation options");
568 DNDS_FIELD(meshFormat, "Mesh format code");
569 DNDS_FIELD(meshPartitionOptions, "Mesh partitioning options");
570 DNDS_FIELD(meshFile, "Input mesh file path");
571 DNDS_FIELD(outPltName, "Output plot file base name");
572 DNDS_FIELD(outLogName, "Output log file base name (empty=use outPltName)");
573 DNDS_FIELD(outRestartName, "Output restart file base name (empty=use outPltName)");
574 DNDS_FIELD(outPltMode, "Output mode: 0=serial, 1=distributed");
575 DNDS_FIELD(readMeshMode, "Read mesh mode: 0=serial CGNS, 1=distributed JSON");
576 DNDS_FIELD(outPltTecplotFormat, "Output in Tecplot format");
577 DNDS_FIELD(outPltVTKFormat, "Output in VTK XML format");
578 DNDS_FIELD(outPltVTKHDFFormat, "Output in VTK HDF format");
579 DNDS_FIELD(outAtPointData, "Output point-interpolated data");
580 DNDS_FIELD(outAtCellData, "Output cell-centered data");
581 DNDS_FIELD(nASCIIPrecision, "ASCII output floating-point precision",
583 DNDS_FIELD(vtuFloatEncodeMode, "VTU float encoding: binary or ascii",
584 DNDS::Config::enum_values({"binary", "ascii"}));
585 DNDS_FIELD(hdfChunkSize, "HDF5 chunk size for output",
587 DNDS_FIELD(hdfDeflateLevel, "HDF5 deflate compression level",
588 DNDS::Config::range(0, 9));
589 DNDS_FIELD(hdfCollOnMeta, "HDF5 collective I/O on metadata");
590 DNDS_FIELD(hdfCollOnData, "HDF5 collective I/O on data arrays");
591 DNDS_FIELD(outVolumeData, "Output volume data");
592 DNDS_FIELD(outBndData, "Output boundary data");
593 DNDS_FIELD(outCellScalarNames, "Additional cell scalar names to output");
594 DNDS_FIELD(serializerSaveURec, "Save reconstruction in restart");
595 DNDS_FIELD(allowAsyncPrintData, "Allow asynchronous data output");
596 DNDS_FIELD(rectifyNearPlane, "Rectify nodes near planes: bitmask 1=x,2=y,4=z");
597 DNDS_FIELD(rectifyNearPlaneThres, "Threshold for near-plane rectification",
599 config.field_section(&T::restartWriter, "restartWriter", "Restart file serializer settings");
600 config.field_section(&T::meshPartitionedWriter, "meshPartitionedWriter", "Partitioned mesh serializer settings");
601 DNDS_FIELD(meshPartitionedReaderType, "Partitioned mesh reader type",
602 DNDS::Config::enum_values({"JSON", "H5"}));
603 // clang-format on
604 }
606
607 /**
608 * @brief Periodic boundary geometry definitions.
609 *
610 * Defines up to 3 periodic translation vectors and rotation parameters
611 * for periodic boundary matching.
612 */
614 {
615 Eigen::Vector<real, -1> PeriodicTranslation1;
616 Eigen::Vector<real, -1> PeriodicTranslation2;
617 Eigen::Vector<real, -1> PeriodicTranslation3;
618 Eigen::Vector<real, -1> PeriodicRotationCent1;
619 Eigen::Vector<real, -1> PeriodicRotationCent2;
620 Eigen::Vector<real, -1> PeriodicRotationCent3;
626 {
627 PeriodicTranslation1.resize(3);
628 PeriodicTranslation2.resize(3);
629 PeriodicTranslation3.resize(3);
630 PeriodicTranslation1 << 1, 0, 0;
631 PeriodicTranslation2 << 0, 1, 0;
632 PeriodicTranslation3 << 0, 0, 1;
633 PeriodicRotationCent1.setZero(3);
634 PeriodicRotationCent2.setZero(3);
635 PeriodicRotationCent3.setZero(3);
639 }
640
642 {
643 // clang-format off
644 DNDS_FIELD(PeriodicTranslation1, "Periodic translation vector for pair 1");
645 DNDS_FIELD(PeriodicTranslation2, "Periodic translation vector for pair 2");
646 DNDS_FIELD(PeriodicTranslation3, "Periodic translation vector for pair 3");
647 DNDS_FIELD(PeriodicRotationCent1, "Rotation center for periodic pair 1");
648 DNDS_FIELD(PeriodicRotationCent2, "Rotation center for periodic pair 2");
649 DNDS_FIELD(PeriodicRotationCent3, "Rotation center for periodic pair 3");
650 DNDS_FIELD(PeriodicRotationEulerAngles1, "Rotation Euler angles (deg) for periodic pair 1");
651 DNDS_FIELD(PeriodicRotationEulerAngles2, "Rotation Euler angles (deg) for periodic pair 2");
652 DNDS_FIELD(PeriodicRotationEulerAngles3, "Rotation Euler angles (deg) for periodic pair 3");
653 DNDS_FIELD(periodicTolerance, "Tolerance for periodic node matching",
655 // clang-format on
656 }
658
659 /**
660 * @brief Slope limiter control parameters.
661 *
662 * Controls whether limiters are active, which limiter variant to use
663 * (WBAP, CWBAP), positivity-preserving reconstruction limiting, and
664 * partial/preserved limiting strategies.
665 */
667 {
668 bool useLimiter = true;
669 bool usePPRecLimiter = true;
670 bool useViscousLimited = true;
673 int nPartialLimiterStart = INT_MAX;
675 bool preserveLimited = false;
677
679 {
680 // clang-format off
681 DNDS_FIELD(useLimiter, "Enable slope limiter");
682 DNDS_FIELD(usePPRecLimiter, "Enable positivity-preserving reconstruction limiter");
683 DNDS_FIELD(useViscousLimited, "Apply limiter to viscous reconstruction");
684 DNDS_FIELD(smoothIndicatorProcedure, "Smooth indicator procedure index");
685 DNDS_FIELD(limiterProcedure, "Limiter variant: 0=WBAP (V2), 1=CWBAP (V3)");
686 DNDS_FIELD(nPartialLimiterStart, "Time step to begin partial limiting");
687 DNDS_FIELD(nPartialLimiterStartLocal, "Time step to begin local partial limiting");
688 DNDS_FIELD(preserveLimited, "Preserve limited reconstruction across steps");
689 DNDS_FIELD(ppRecLimiterCompressToMean, "PP limiter compresses toward cell mean");
690 // clang-format on
691 }
693
694 /**
695 * @brief Linear solver control parameters.
696 *
697 * Controls the implicit linear solver: preconditioner type (Jacobi/GS/ILU),
698 * SGS iterations, GMRES settings, multi-grid options, and per-level
699 * coarse-grid solver configurations.
700 */
702 {
703 int jacobiCode = 1; // 0 for jacobi, 1 for gs, 2 for ilu
704 int sgsIter = 0;
705 int sgsWithRec = 0;
706 int gmresCode = 0; // 0 for lusgs, 1 for gmres, 2 for lusgs started gmres
707 int gmresScale = 0; // 0 for no scaling, 1 use refU, 2 use mean value
708 int nGmresSpace = 10;
709 int nGmresIter = 2;
713 int multiGridLP = 0;
718 {
719 int jacobiCode = 0; // 0 for jacobi, 1 for gs, 2 for ilu
720 int sgsIter = 0;
721 int gmresCode = 0; // 0 for lusgs, 1 for gmres, 2 for lusgs started gmres
722 int gmresScale = 0; // 0 for no scaling, 1 use refU, 2 use mean value
723 int nGmresIter = 2;
729
731 {
732 // clang-format off
733 DNDS_FIELD(jacobiCode, "Preconditioner: 0=jacobi, 1=GS, 2=ILU");
734 DNDS_FIELD(sgsIter, "SGS iteration count",
736 DNDS_FIELD(gmresCode, "Linear solver: 0=LUSGS, 1=GMRES, 2=LUSGS+GMRES");
737 DNDS_FIELD(gmresScale, "GMRES scaling: 0=none, 1=refU, 2=mean");
738 DNDS_FIELD(nGmresIter, "GMRES iterations",
740 DNDS_FIELD(nSgsConsoleCheck, "Console check interval for SGS",
742 DNDS_FIELD(nGmresConsoleCheck, "Console check interval for GMRES",
744 DNDS_FIELD(multiGridNIter, "Multi-grid pre-smooth iterations (-1=auto)");
745 DNDS_FIELD(multiGridNIterPost, "Multi-grid post-smooth iterations",
747 DNDS_FIELD(centralSmoothInputResidual, "Smooth input residual on coarse grid");
748 // clang-format on
749 }
750 };
751 std::map<std::string, CoarseGridLinearSolverControl> coarseGridLinearSolverControlList{
754 };
756
758 {
759 // clang-format off
760 DNDS_FIELD(jacobiCode, "Preconditioner: 0=jacobi, 1=GS, 2=ILU");
761 DNDS_FIELD(sgsIter, "SGS iteration count",
763 DNDS_FIELD(sgsWithRec, "Couple SGS with reconstruction");
764 DNDS_FIELD(gmresCode, "Linear solver: 0=LUSGS, 1=GMRES, 2=LUSGS+GMRES");
765 DNDS_FIELD(gmresScale, "GMRES scaling: 0=none, 1=refU, 2=mean");
766 DNDS_FIELD(nGmresSpace, "GMRES Krylov subspace size",
768 DNDS_FIELD(nGmresIter, "GMRES outer iterations",
770 DNDS_FIELD(nSgsConsoleCheck, "Console check interval for SGS",
772 DNDS_FIELD(nGmresConsoleCheck, "Console check interval for GMRES",
774 DNDS_FIELD(initWithLastURecInc, "Initialize GMRES with last uRec increment");
775 DNDS_FIELD(multiGridLP, "Multi-grid levels (0=off)");
776 DNDS_FIELD(multiGridLPInnerNIter, "Multi-grid inner iterations",
778 DNDS_FIELD(multiGridLPStartIter, "Iteration to enable multi-grid",
780 DNDS_FIELD(multiGridLPInnerNSee, "Multi-grid inner console see interval",
782 config.template field_map_of<CoarseGridLinearSolverControl>(
783 &T::coarseGridLinearSolverControlList,
784 "coarseGridLinearSolverControlList",
785 "Per-level coarse grid linear solver settings");
786 config.field_section(&T::directPrecControl, "directPrecControl",
787 "Direct preconditioner settings");
788 // clang-format on
789 }
791
792 /**
793 * @brief Restart checkpoint state.
794 *
795 * Records the step indices and file paths for restart continuation,
796 * including support for reading restart files from a different solver
797 * (with dimension remapping).
798 */
800 {
801 int iStep = -1;
803 int odeCodePrev = -1;
804 std::string lastRestartFile = "";
805 std::string otherRestartFile = "";
806 std::vector<int> otherRestartStoreDim;
808 {
809 // clang-format off
810 DNDS_FIELD(iStep, "Restart step index");
811 DNDS_FIELD(iStepInternal, "Restart internal step index");
812 DNDS_FIELD(odeCodePrev, "Previous ODE code for restart");
813 DNDS_FIELD(lastRestartFile, "Path to last restart file");
814 DNDS_FIELD(otherRestartFile, "Path to alternate restart file");
815 DNDS_FIELD(otherRestartStoreDim, "Dimension mapping for alternate restart");
816 // clang-format on
817 }
819 {
820 otherRestartStoreDim.resize(1);
821 for (int i = 0; i < otherRestartStoreDim.size(); i++)
823 }
825
826 /// @brief Time-averaging control for unsteady simulations.
828 {
829 bool enabled = false;
830
832 {
833 // clang-format off
834 DNDS_FIELD(enabled, "Enable time-averaging of solution fields");
835 // clang-format on
836 }
838
839 /// @brief Miscellaneous solver options (axisymmetric mode, passive scalar freezing, rec matrix output).
840 struct Others
841 {
844 bool printRecMatrix = false;
846
848 {
849 // clang-format off
850 DNDS_FIELD(nFreezePassiveInner, "Freeze passive scalars for N inner steps",
852 DNDS_FIELD(axisSymmetric, "Axisymmetric mode: 0=off");
853 DNDS_FIELD(printRecMatrix, "Print reconstruction matrix to file");
854 config.field_section(&T::recMatrixWriter, "recMatrixWriter",
855 "Serializer for reconstruction matrix output");
856 // clang-format on
857 }
859
860 EulerEvaluatorSettings<model> eulerSettings; ///< Physics settings passed to the EulerEvaluator.
861 CFV::VRSettings vfvSettings; ///< Variational reconstruction settings.
862 nlohmann::ordered_json bcSettings = nlohmann::ordered_json::array(); ///< Boundary condition definitions (JSON array).
863 std::map<std::string, std::string> bcNameMapping; ///< Mapping from mesh BC names to solver BC type names.
864
866 {
867 // clang-format off
868 config.field_section(&T::timeMarchControl, "timeMarchControl", "Time marching settings");
869 config.field_section(&T::implicitReconstructionControl, "implicitReconstructionControl", "Implicit reconstruction settings");
870 config.field_section(&T::outputControl, "outputControl", "Output settings");
871 config.field_section(&T::implicitCFLControl, "implicitCFLControl", "Implicit CFL settings");
872 config.field_section(&T::convergenceControl, "convergenceControl", "Convergence monitoring settings");
873 config.field_section(&T::dataIOControl, "dataIOControl", "Data I/O and restart settings");
874 config.field_section(&T::boundaryDefinition, "boundaryDefinition", "Periodic boundary geometry");
875 config.field_section(&T::limiterControl, "limiterControl", "Slope limiter settings");
876 config.field_section(&T::linearSolverControl, "linearSolverControl", "Linear solver settings");
877 config.field_section(&T::timeAverageControl, "timeAverageControl", "Time averaging settings");
878 config.field_section(&T::others, "others", "Miscellaneous settings");
879 config.field_section(&T::restartState, "restartState", "Restart checkpoint state");
880 config.field_section(&T::eulerSettings, "eulerSettings", "Euler evaluator physics settings");
881 config.field_section(&T::vfvSettings, "vfvSettings", "Variational reconstruction settings");
882 config.field_json_schema(&T::bcSettings, "bcSettings",
883 "Boundary condition settings (per-BC array)",
884 []() { return bcSettingsSchema(); });
885 DNDS_FIELD(bcNameMapping, "Boundary name to type mapping");
886
887 config.check("bcSettings must be a JSON array", [](const T &s)
888 {
889 return s.bcSettings.is_array();
890 });
891 // clang-format on
892 }
893
894 /// @brief Backward-compatible bidirectional JSON read/write.
895 ///
896 /// Delegates to the auto-generated from_json / to_json from
897 /// DNDS_DECLARE_CONFIG. Kept for call-site compatibility.
898 void
899 ReadWriteJson(nlohmann::ordered_json &jsonObj, int nVars, bool read)
900 {
901 (void)nVars; // nVars is already stored in eulerSettings via Configuration(nVars)
902 if (read)
903 from_json(jsonObj, *this);
904 else
905 to_json(jsonObj, *this);
906 }
907
909 {
910 }
911
912 Configuration(int nVars)
914 {
916 }
917
918 } config = Configuration{};
919
920 public:
921 /**
922 * @brief Construct an EulerSolver with MPI communicator information.
923 *
924 * Initializes the runtime variable count, output field sizes, and default
925 * configuration. Does not read mesh or allocate solver arrays.
926 *
927 * @param nmpi MPI communicator information.
928 * @param n_nVars Runtime number of conserved variables (default from model).
929 */
930 EulerSolver(const MPIInfo &nmpi, int n_nVars = getNVars(model)) : nVars(n_nVars), mpi(nmpi)
931 {
932 if (getNVars(model) == DynamicSize)
933 DNDS_assert_info(nVars >= getDim_Fixed(model) + 2, "nVars too small");
934 else
935 DNDS_assert_info(nVars == getNVars(model), "do not change the nVars for this model");
936 nOUTS = nVars + 4;
937 nOUTSPoint = nVars + 2;
938 nOUTSBnd = nVars + 2 + nVars + 3 + 1 + 3; // Uprim + (T,M) + F + Ft + faceZone + Norm
939
940 config = Configuration(nVars); //* important to initialize using nVars
941 }
942
943 /// @brief Destructor. Waits for all async output futures to complete.
945 {
946 int nBad{0};
947 do
948 {
949 nBad = 0;
950 for (auto &f : outFuture)
951 if (f.valid() && f.wait_for(std::chrono::microseconds(10)) != std::future_status::ready)
952 nBad++;
953 for (auto &f : outBndFuture)
954 if (f.valid() && f.wait_for(std::chrono::microseconds(10)) != std::future_status::ready)
955 nBad++;
956 if (outSeqFuture.valid() && outSeqFuture.wait_for(std::chrono::microseconds(10)) != std::future_status::ready)
957 nBad++;
958 } while (nBad);
959 }
960
961 /**
962 * @brief Load or write solver configuration from/to a JSON file.
963 *
964 * When read=true, parses the JSON file, optionally merges a patch file,
965 * applies key/value overrides, populates the Configuration struct and
966 * creates the boundary condition handler. When read=false, serializes the
967 * current configuration to a JSON file.
968 *
969 * @param jsonName Path to the JSON configuration file.
970 * @param read true=read from file, false=write to file.
971 * @param jsonMergeName Optional path to a JSON patch file to merge.
972 * @param overwriteKeys JSON pointer paths to override (e.g., "/timeMarchControl/CFL").
973 * @param overwriteValues Values corresponding to overwriteKeys.
974 */
975 void ConfigureFromJson(const std::string &jsonName, bool read = false, const std::string &jsonMergeName = "",
976 const std::vector<std::string> &overwriteKeys = {}, const std::vector<std::string> &overwriteValues = {})
977 {
978 if (read)
979 {
980 auto fIn = std::ifstream(jsonName);
981 DNDS_assert_info(fIn, "config file not existent");
982 gSetting = nlohmann::ordered_json::parse(fIn, nullptr, true, true);
983
984 if (read && !jsonMergeName.empty())
985 {
986 fIn = std::ifstream(jsonMergeName);
987 DNDS_assert_info(fIn, "config file patch not existent");
988 auto gSettingAdd = nlohmann::ordered_json::parse(fIn, nullptr, true, true);
989 gSetting.merge_patch(gSettingAdd);
990 }
991 DNDS_assert(overwriteKeys.size() == overwriteValues.size());
992 for (size_t i = 0; i < overwriteKeys.size(); i++)
993 {
994
995 auto key = nlohmann::ordered_json::json_pointer(overwriteKeys[i].c_str());
996 try
997 {
998 std::string valString =
999 fmt::format(R"({{
1000 "__val_entry": {}
1001 }})",
1002 overwriteValues[i]);
1003 auto valDoc = nlohmann::ordered_json::parse(valString, nullptr, true, true);
1004 if (mpi.rank == 0)
1005 log() << "JSON: overwrite key: " << key << std::endl
1006 << "JSON: overwrite val: " << valDoc["__val_entry"] << std::endl;
1007 gSetting[key] = valDoc["__val_entry"];
1008 }
1009 catch (const nlohmann::ordered_json::parse_error &e)
1010 {
1011 try
1012 {
1013 std::string valString =
1014 fmt::format(R"({{
1015"__val_entry": "{}"
1016}})",
1017 overwriteValues[i]);
1018 auto valDoc = nlohmann::ordered_json::parse(valString, nullptr, true, true);
1019 if (mpi.rank == 0)
1020 log() << "JSON: overwrite key: " << key << std::endl
1021 << "JSON: overwrite val: " << valDoc["__val_entry"] << std::endl;
1022 gSetting[key] = valDoc["__val_entry"];
1023 }
1024 catch (const std::exception &e)
1025 {
1026 std::cerr << e.what() << "\n";
1027 std::cerr << overwriteValues[i] << "\n";
1028 DNDS_assert(false);
1029 }
1030 }
1031 catch (const std::exception &e)
1032 {
1033 std::cerr << e.what() << "\n";
1034 std::cerr << overwriteValues[i] << "\n";
1035 DNDS_assert(false);
1036 }
1037 }
1038 config.ReadWriteJson(gSetting, nVars, read);
1039 // create from json the pBCHandler
1040 pBCHandler = std::make_shared<BoundaryHandler<model>>(nVars);
1041 from_json(config.bcSettings, *pBCHandler);
1042 gSetting["bcSettings"] = *pBCHandler;
1043 PrintConfig(true);
1044 if (mpi.rank == 0)
1045 log() << "JSON: read value:" << std::endl
1046 << std::setw(4) << gSetting << std::endl;
1047 }
1048 else
1049 {
1050 gSetting = nlohmann::ordered_json::object();
1051 config.ReadWriteJson(gSetting, nVars, read);
1052 if (pBCHandler) // todo: add example pBCHandler
1053 gSetting["bcSettings"] = *pBCHandler;
1054 if (mpi.rank == 0) // single call for output
1055 {
1056 std::filesystem::path outFile{jsonName};
1057 std::filesystem::create_directories(outFile.parent_path() / ".");
1058 auto fIn = std::ofstream(jsonName);
1059 DNDS_assert(fIn);
1060 fIn << std::setw(4) << gSetting;
1061 }
1062 MPI::Barrier(mpi.comm); // no go until output done
1063 }
1064
1065 if (mpi.rank == 0)
1066 log() << "JSON: Parse " << (read ? "read" : "write")
1067 << " Done ===" << std::endl;
1068 }
1069
1070 /**
1071 * @brief Read the mesh and initialize the full solver pipeline.
1072 *
1073 * Performs the complete initialization sequence:
1074 * 1. Read mesh (serial CGNS or distributed)
1075 * 2. Apply coordinate transformations (rotation, scaling, rectification)
1076 * 3. Partition the mesh across MPI ranks (Metis/ParMetis)
1077 * 4. Apply order elevation (O1->O2) and mesh bisection
1078 * 5. Build wall distance field (if requested)
1079 * 6. Construct the variational reconstruction (VFV)
1080 * 7. Create the EulerEvaluator and initialize FV infrastructure
1081 * 8. Allocate DOF, reconstruction, and Jacobian arrays
1082 * 9. Load restart data (if configured)
1083 * 10. Set initial conditions (if not restarting)
1084 */
1085 void ReadMeshAndInitialize();
1086
1087 /**
1088 * @brief Write the current configuration to a JSON log file.
1089 *
1090 * Extracts live settings from the VFV, evaluator, and BC handler,
1091 * serializes the full configuration to JSON, and writes it alongside
1092 * compile-time defines and commit information.
1093 *
1094 * @param updateCommit If true, include the git commit hash in the output.
1095 */
1096 void PrintConfig(bool updateCommit = false)
1097 {
1098 /***********************************************************/
1099 // if these objects are existent, extract settings from them
1100 if (vfv)
1101 config.vfvSettings = static_cast<const CFV::VRSettings &>(vfv->getSettings());
1102 if (pEval)
1103 config.eulerSettings = pEval->settings;
1104 if (pBCHandler)
1105 config.bcSettings = *pBCHandler;
1106 /***********************************************************/
1107 config.ReadWriteJson(gSetting, nVars, false);
1108 if (mpi.rank == 0)
1109 {
1110 std::string logConfigFileName = config.dataIOControl.getOutLogName() + "_" + output_stamp + ".config.json";
1111 std::filesystem::path outFile{logConfigFileName};
1112 std::filesystem::create_directories(outFile.parent_path() / ".");
1113 std::ofstream logConfig(logConfigFileName);
1114 DNDS_assert(logConfig);
1115 gSetting["___Compile_Time_Defines"] = DNDS_Defines_state;
1116 gSetting["___Runtime_PartitionNumber"] = mpi.size;
1117 gSetting["___Commit_ID"] = GetSetVersionName();
1118 // if (updateCommit)
1119 // {
1120 // std::ifstream commitIDFile("commitID.txt");
1121 // if (commitIDFile)
1122 // {
1123 // std::string commitHash;
1124 // commitIDFile >> commitHash;
1125 // gSetting["___Commit_ID"] = commitHash;
1126 // }
1127 // }
1128 logConfig << std::setw(4) << gSetting;
1129 logConfig.close();
1130 }
1131 }
1132 /// @brief Read a restart file and populate u (and optionally uRec) from it.
1133 void ReadRestart(std::string fname);
1134
1135 /// @brief Read a restart file from a different solver/model, remapping variable dimensions.
1136 void ReadRestartOtherSolver(std::string fname, const std::vector<int> &dimStore);
1137
1138 /// @brief Write the current solution state to a restart file.
1139 void PrintRestart(std::string fname);
1140
1141 using tAdditionalCellScalarList = tCellScalarList; ///< Type alias for additional output scalar list.
1142
1143 /// @brief Output mode selector for PrintData.
1145 {
1146 PrintDataLatest = 0, ///< Output the current (latest) solution.
1147 PrintDataTimeAverage = 1, ///< Output the time-averaged solution.
1148 };
1149
1150 /**
1151 * @brief Write solution data to VTK/HDF5/Tecplot output files.
1152 *
1153 * Gathers distributed data to serial (or writes distributed), converts
1154 * conservative to primitive variables, appends additional scalars (residual,
1155 * limiter flags), and writes the output in the configured format(s).
1156 *
1157 * @param fname Output file base name.
1158 * @param fnameSeries Series file name (for VTK time series).
1159 * @param odeResidualF Callback returning the ODE residual for each cell.
1160 * @param additionalCellScalars Additional per-cell scalar fields to output.
1161 * @param eval Reference to the evaluator (for output field computation).
1162 * @param TSimu Current simulation time (-1 for steady).
1163 * @param mode PrintDataLatest or PrintDataTimeAverage.
1164 */
1165 void PrintData(const std::string &fname, const std::string &fnameSeries,
1166 const tCellScalarFGet &odeResidualF,
1167 tAdditionalCellScalarList &additionalCellScalars,
1168 TEval &eval, real TSimu = -1.0, PrintDataMode mode = PrintDataLatest);
1169
1170 /// @brief Serialize the solution to a SerializerBase (currently unused standalone path).
1171 void WriteSerializer(Serializer::SerializerBaseSSP serializerP, const std::string &name) // currently not using
1172 {
1173 auto cwd = serializerP->GetCurrentPath();
1174 serializerP->CreatePath(name);
1175 serializerP->GoToPath(name);
1176
1177 u.WriteSerialize(serializerP, "meanValue");
1178
1179 serializerP->GoToPath(cwd);
1180
1181 nlohmann::ordered_json configJson;
1182 config.ReadWriteJson(configJson, nVars, false);
1183 serializerP->WriteString("lastConfig", configJson.dump());
1184
1186 {
1187 serializerP->WriteInt("hasReconstructionValue", 1);
1188 uRec.WriteSerialize(serializerP, "recValue");
1189 }
1190 else
1191 serializerP->WriteInt("hasReconstructionValue", 0);
1192 }
1193
1194 /// @brief Fill a log value map entry for arithmetic types (scalars).
1195 template <class TVal>
1196 std::enable_if_t<std::is_arithmetic_v<std::remove_reference_t<TVal>>>
1197 FillLogValue(tLogSimpleDIValueMap &v_map, const std::string &name, TVal &&val)
1198 {
1199 v_map[name] = val;
1200 }
1201
1202 /// @brief Fill a log value map entry for non-arithmetic types (vectors → per-component entries).
1203 template <class TVal>
1204 std::enable_if_t<!std::is_arithmetic_v<std::remove_reference_t<TVal>>>
1205 FillLogValue(tLogSimpleDIValueMap &v_map, const std::string &name, TVal &&val)
1206 {
1207 // std::vector<std::string> logfileOutputTitles{
1208 // "step", "iStep", "iter", "tSimu",
1209 // "res", "curDtImplicit", "curDtMin", "CFLNow",
1210 // "nLimInc", "alphaMinInc",
1211 // "nLimBeta", "minBeta",
1212 // "nLimAlpha", "minAlpha",
1213 // "tWall", "telapsed", "trec", "trhs", "tcomm", "tLim", "tLimiterA", "tLimiterB",
1214 // "fluxWall", "CL", "CD", "AoA"};
1215 if (name == "res" || name == "fluxWall")
1216 for (int i = 0; i < nVars; i++)
1217 v_map[name + std::to_string(i)] = val[i];
1218 else
1219 v_map[name] = 0;
1220 }
1221
1222 /// @brief Initialize the CSV error logger and value map for convergence monitoring.
1223 /// @return Tuple of (CsvLog writer, value map with all column names initialized).
1224 std::tuple<std::unique_ptr<CsvLog>, tLogSimpleDIValueMap> LogErrInitialize()
1225 {
1227 TU initVec;
1228 initVec.setZero(nVars);
1229 std::vector<std::string> realNames;
1230 for (auto name : config.outputControl.logfileOutputTitles)
1231 if (name == "res" || name == "fluxWall")
1232 {
1233 FillLogValue(v_map, name, initVec);
1234 for (int i = 0; i < nVars; i++)
1235 realNames.push_back(name + std::to_string(i));
1236 }
1237 else
1238 FillLogValue(v_map, name, 0.), realNames.push_back(name);
1239
1240 // std::string logErrFileName = config.dataIOControl.getOutLogName() + "_" + output_stamp + ".log";
1241 // std::filesystem::path outFile{logErrFileName};
1242 // std::filesystem::create_directories(outFile.parent_path() / ".");
1243 // auto pOs = std::make_unique<std::ofstream>(logErrFileName);
1244 // DNDS_assert_info(*pOs, "logErr file [" + logErrFileName + "] did not open");
1245
1246 return std::make_tuple(std::make_unique<DNDS::CsvLog>(
1247 realNames,
1249 ".log",
1250 100000),
1251 v_map);
1252 }
1253
1254 /**
1255 * @brief Run the main implicit time-marching loop.
1256 *
1257 * Executes the outer time-step loop with configurable ODE integrators
1258 * (backward Euler, BDF2, SDIRK, etc.). Each step involves:
1259 * 1. CFL ramping and time-step computation
1260 * 2. Variational reconstruction (implicit or explicit)
1261 * 3. Slope limiting and positivity-preserving enforcement
1262 * 4. RHS evaluation
1263 * 5. Implicit linear solve (LU-SGS, GMRES, or direct)
1264 * 6. Solution update with increment compression
1265 * 7. Convergence monitoring and data output
1266 *
1267 * Supports steady-state convergence checks, CL-driver AoA adaptation,
1268 * time-averaging, and restart checkpointing.
1269 */
1271
1272 /**
1273 * @brief Mutable state bundle for the time-marching loop.
1274 *
1275 * Holds all transient state that persists across time steps within
1276 * RunImplicitEuler: evaluator, ODE integrator, linear solvers, timing
1277 * statistics, residual bases, convergence trackers, output counters,
1278 * and CFL/dt history. The DNDS_EULERSOLVER_RUNNINGENV_GET_REF_LIST
1279 * macro provides convenient local aliases for all members.
1280 */
1282 {
1284 std::tuple<std::unique_ptr<CsvLog>, tLogSimpleDIValueMap> logErr;
1286 std::unique_ptr<tGMRES_u> gmres;
1287 std::unique_ptr<tGMRES_uRec> gmresRec;
1288 std::unique_ptr<tPCG_uRec> pcgRec;
1289 std::unique_ptr<tPCG_uRec> pcgRec1;
1290 double tstart = 0;
1291 double tstartInternal = 0;
1292 std::map<std::string, ScalarStatistics> tInternalStats;
1293 int stepCount = 0;
1299 int nextStepOut = -1;
1305
1307 bool ifOutT = false;
1310 std::vector<real> curDtImplicitHistory;
1311 int step = 0;
1312 int iterAll = 0;
1313 bool gradIsZero = true;
1314
1321
1323
1325
1326#define DNDS_EULERSOLVER_RUNNINGENV_GET_REF(name) auto &name = runningEnvironment.name
1327
1328#define DNDS_EULERSOLVER_RUNNINGENV_GET_REF_LIST \
1329 auto &eval = *runningEnvironment.pEval; \
1330 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(logErr); \
1331 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(ode); \
1332 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(gmres); \
1333 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(gmresRec); \
1334 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(pcgRec); \
1335 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(pcgRec1); \
1336 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(tstart); \
1337 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(tstartInternal); \
1338 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(tInternalStats); \
1339 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(stepCount); \
1340 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(resBaseC); \
1341 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(resBaseCInternal); \
1342 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(tSimu); \
1343 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(tAverage); \
1344 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(nextTout); \
1345 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(nextStepOut); \
1346 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(nextStepOutC); \
1347 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(nextStepRestart); \
1348 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(nextStepRestartC); \
1349 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(nextStepOutAverage); \
1350 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(nextStepOutAverageC); \
1351 \
1352 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(CFLNow); \
1353 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(ifOutT); \
1354 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(curDtMin); \
1355 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(curDtImplicit); \
1356 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(curDtImplicitHistory); \
1357 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(step); \
1358 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(iterAll); \
1359 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(gradIsZero); \
1360 \
1361 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(nLimBeta); \
1362 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(nLimAlpha); \
1363 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(minAlpha); \
1364 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(minBeta); \
1365 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(nLimInc); \
1366 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(alphaMinInc); \
1367 \
1368 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(dtIncreaseCounter); \
1369 \
1370 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(addOutList);
1371
1373 };
1374 /// @brief Populate a RunningEnvironment with allocated solvers, loggers, and initial state.
1376
1377 /**
1378 * @brief Solve the implicit linear system at a given grid level.
1379 *
1380 * Dispatches to LU-SGS, GMRES, or LU-SGS-preconditioned GMRES depending on
1381 * configuration. Used within the ODE integrator's inner loop.
1382 *
1383 * @warning Explicit template instantiation does not exist; inlined only.
1384 */
1386 real alphaDiag, real t,
1387 TDof &cres, TDof &cx, TDof &cxInc, TRec &uRecC, TRec uRecIncC,
1388 JacobianDiagBlock<nVarsFixed> &JDC, tGMRES_u &gmres, int gridLevel);
1389 /**
1390 * @brief Apply the preconditioner (SGS sweeps or ILU) to a right-hand side vector.
1391 *
1392 * @warning Explicit template instantiation does not exist; inlined only.
1393 */
1394 void doPrecondition(real alphaDiag, real t,
1395 TDof &crhs, TDof &cx, TDof &cxInc, TDof &uTemp,
1396 JacobianDiagBlock<nVarsFixed> &JDC, TU &sgsRes, bool &inputIsZero, bool &hasLUDone, int gridLevel);
1397
1398 /// @brief Convergence/termination check functor called after each inner iteration.
1399 /// @return true if the inner loop should stop (converged or max iterations reached).
1400 bool functor_fstop(int iter, ArrayDOFV<nVarsFixed> &cres, int iStep, RunningEnvironment &env);
1401 /// @brief Main outer-loop functor: performs one full time step (reconstruction, RHS, linear solve, update).
1402 /// @return true if the outer loop should continue.
1404
1405 /// @name Accessors
1406 /// @{
1407 auto getMPI() const { return mpi; } ///< Get MPI communicator info.
1408 auto getMesh() const { return mesh; } ///< Get shared pointer to the mesh.
1409 auto getVFV() const { return vfv; } ///< Get shared pointer to the VFV reconstruction.
1410 /// @}
1411
1412 /// @name Test accessors (for unit testing the evaluator pipeline)
1413 /// @{
1414 auto &getU() { return u; } ///< Get mutable reference to the DOF array.
1415 auto &getURec() { return uRec; } ///< Get mutable reference to the reconstruction array.
1416 auto &getURecNew() { return uRecNew; } ///< Get mutable reference to the new reconstruction array.
1417 auto &getURecLimited() { return uRecLimited; } ///< Get mutable reference to the limited reconstruction array.
1418 auto &getEval() { return *pEval; } ///< Get mutable reference to the evaluator.
1419 auto &getConfiguration() { return config; } ///< Get mutable reference to the configuration.
1420 auto &getJSource() { return JSource; } ///< Get mutable reference to the source Jacobian.
1421 auto &getBetaPP() { return betaPP; } ///< Get mutable reference to the PP beta array.
1422 auto &getAlphaPP() { return alphaPP; } ///< Get mutable reference to the PP alpha array.
1423 auto &getDTauTmp() { return dTauTmp; } ///< Get mutable reference to the dTau temporary.
1424 auto &getIfUseLimiter() { return ifUseLimiter; } ///< Get mutable reference to the limiter flag array.
1425 auto &getBCHandler() { return pBCHandler; } ///< Get mutable reference to the BC handler.
1426 /// @}
1427 };
1428}
1429
1430#define DNDS_EULERSOLVER_INS_EXTERN(model, ext) \
1431 namespace DNDS::Euler \
1432 { \
1433 ext template void EulerSolver<model>::RunImplicitEuler(); \
1434 ext template void EulerSolver<model>::InitializeRunningEnvironment( \
1435 EulerSolver<model>::RunningEnvironment &env); \
1436 }
1437
1445
1446#define DNDS_EULERSOLVER_PRINTDATA_INS_EXTERN(model, ext) \
1447 namespace DNDS::Euler \
1448 { \
1449 ext template void EulerSolver<model>::PrintData( \
1450 const std::string &fname, const std::string &fnameSeries, \
1451 const tCellScalarFGet &odeResidualF, \
1452 tAdditionalCellScalarList &additionalCellScalars, \
1453 TEval &eval, real tSimu, \
1454 PrintDataMode mode); \
1455 ext template void EulerSolver<model>::PrintRestart( \
1456 std::string fname); \
1457 ext template void EulerSolver<model>::ReadRestart( \
1458 std::string fname); \
1459 ext template void EulerSolver<model>::ReadRestartOtherSolver( \
1460 std::string fname, const std::vector<int> &dimStore); \
1461 }
1462
1470
1471#define DNDS_EULERSOLVER_INIT_INS_EXTERN(model, ext) \
1472 namespace DNDS::Euler \
1473 { \
1474 ext template void EulerSolver<model>::ReadMeshAndInitialize(); \
1475 ext template bool EulerSolver<model>::functor_fstop( \
1476 int iter, ArrayDOFV<nVarsFixed> &cres, int iStep, RunningEnvironment &env); \
1477 ext template bool EulerSolver<model>::functor_fmainloop( \
1478 RunningEnvironment &env); \
1479 }
1480
pybind11-style configuration registration with macro-based field declaration and namespace-scoped tag...
#define DNDS_FIELD(name_, desc_,...)
Register a field inside a DNDS_DECLARE_CONFIG body.
Streaming CSV writer with auto-rotation and a dual int-or-double value type.
#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
Boundary condition types, handlers, integration recorders, and 1-D profile utilities for the compress...
Core finite-volume evaluator for compressible Navier-Stokes / Euler equations.
#define DNDS_EULERSOLVER_INIT_INS_EXTERN(model, ext)
#define DNDS_EULERSOLVER_PRINTDATA_INS_EXTERN(model, ext)
#define DNDS_EULERSOLVER_INS_EXTERN(model, ext)
Ideal-gas Riemann solvers, flux functions, and thermodynamic utilities for the compressible Euler / N...
JSON-to-Eigen conversion utilities and nlohmann_json helper macros.
Pre-allocated object pool with RAII checkout/return semantics.
Configurable factory that builds either a SerializerJSON or a SerializerH5 with all tunables exposed ...
The VR class that provides any information needed in high-order CFV.
MPI-parallel distributed array of per-cell Degrees-of-Freedom (conservative variable vectors).
Definition Euler.hpp:93
MPI-parallel distributed array of per-cell reconstruction coefficient matrices.
Definition Euler.hpp:401
Per-zone boundary condition handler for Euler/Navier-Stokes solvers.
Definition EulerBC.hpp:212
Core finite-volume evaluator for compressible Navier-Stokes / Euler equations.
static const int nVarsFixed
Number of conserved variables (compile-time fixed).
static const auto I4
Index of the energy equation (= dim+1).
ssp< CFV::VariationalReconstruction< gDim > > TpVFV
Shared pointer to the variational reconstruction object.
static const int gDim
Geometric dimension of the mesh.
static const int dim
Spatial dimension (2 or 3).
Top-level solver orchestrator for compressible Navier-Stokes / Euler equations.
void PrintData(const std::string &fname, const std::string &fnameSeries, const tCellScalarFGet &odeResidualF, tAdditionalCellScalarList &additionalCellScalars, TEval &eval, real TSimu=-1.0, PrintDataMode mode=PrintDataLatest)
Write solution data to VTK/HDF5/Tecplot output files.
auto & getU()
Get mutable reference to the DOF array.
TEval::TRec TRec
Reconstruction coefficient array type.
void ReadRestartOtherSolver(std::string fname, const std::vector< int > &dimStore)
Read a restart file from a different solver/model, remapping variable dimensions.
TEval::TVFV TVFV
Variational reconstruction type.
void WriteSerializer(Serializer::SerializerBaseSSP serializerP, const std::string &name)
Serialize the solution to a SerializerBase (currently unused standalone path).
void RunImplicitEuler()
Run the main implicit time-marching loop.
void ReadMeshAndInitialize()
Read the mesh and initialize the full solver pipeline.
TEval::TDiffU TDiffU
Gradient of conserved variables type.
EulerSolver(const MPIInfo &nmpi, int n_nVars=getNVars(model))
Construct an EulerSolver with MPI communicator information.
bool functor_fstop(int iter, ArrayDOFV< nVarsFixed > &cres, int iStep, RunningEnvironment &env)
Convergence/termination check functor called after each inner iteration.
void solveLinear(real alphaDiag, real t, TDof &cres, TDof &cx, TDof &cxInc, TRec &uRecC, TRec uRecIncC, JacobianDiagBlock< nVarsFixed > &JDC, tGMRES_u &gmres, int gridLevel)
Solve the implicit linear system at a given grid level.
static const int gDim
Geometric dimension of the mesh.
auto & getConfiguration()
Get mutable reference to the configuration.
auto & getAlphaPP()
Get mutable reference to the PP alpha array.
auto & getURecNew()
Get mutable reference to the new reconstruction array.
std::enable_if_t<!std::is_arithmetic_v< std::remove_reference_t< TVal > > > FillLogValue(tLogSimpleDIValueMap &v_map, const std::string &name, TVal &&val)
Fill a log value map entry for non-arithmetic types (vectors → per-component entries).
auto & getBCHandler()
Get mutable reference to the BC handler.
void InitializeRunningEnvironment(RunningEnvironment &env)
Populate a RunningEnvironment with allocated solvers, loggers, and initial state.
TEval::TVec TVec
Spatial vector type.
TEval::TpVFV TpVFV
Shared pointer to VFV type.
TEval::TU TU
Conservative variable vector type.
TEval::TScalar TScalar
Scalar reconstruction coefficient array type.
auto getMesh() const
Get shared pointer to the mesh.
auto & getEval()
Get mutable reference to the evaluator.
void doPrecondition(real alphaDiag, real t, TDof &crhs, TDof &cx, TDof &cxInc, TDof &uTemp, JacobianDiagBlock< nVarsFixed > &JDC, TU &sgsRes, bool &inputIsZero, bool &hasLUDone, int gridLevel)
Apply the preconditioner (SGS sweeps or ILU) to a right-hand side vector.
auto & getBetaPP()
Get mutable reference to the PP beta array.
auto & getDTauTmp()
Get mutable reference to the dTau temporary.
auto & getURec()
Get mutable reference to the reconstruction array.
void ReadRestart(std::string fname)
Read a restart file and populate u (and optionally uRec) from it.
EulerEvaluator< model > TEval
Evaluator type for this model.
~EulerSolver()
Destructor. Waits for all async output futures to complete.
struct DNDS::Euler::EulerSolver::Configuration config
TEval::TJacobianU TJacobianU
Flux Jacobian matrix type.
bool functor_fmainloop(RunningEnvironment &env)
Main outer-loop functor: performs one full time step (reconstruction, RHS, linear solve,...
nlohmann::ordered_json gSetting
Full JSON configuration (read from file, may be modified at runtime).
void PrintConfig(bool updateCommit=false)
Write the current configuration to a JSON log file.
std::enable_if_t< std::is_arithmetic_v< std::remove_reference_t< TVal > > > FillLogValue(tLogSimpleDIValueMap &v_map, const std::string &name, TVal &&val)
Fill a log value map entry for arithmetic types (scalars).
TEval::TMat TMat
Spatial matrix type.
auto & getURecLimited()
Get mutable reference to the limited reconstruction array.
static const int dim
Spatial dimension (2 or 3).
void ConfigureFromJson(const std::string &jsonName, bool read=false, const std::string &jsonMergeName="", const std::vector< std::string > &overwriteKeys={}, const std::vector< std::string > &overwriteValues={})
Load or write solver configuration from/to a JSON file.
auto getMPI() const
Get MPI communicator info.
std::string output_stamp
Unique stamp appended to output filenames for this run.
std::tuple< std::unique_ptr< CsvLog >, tLogSimpleDIValueMap > LogErrInitialize()
Initialize the CSV error logger and value map for convergence monitoring.
static const int nVarsFixed
Compile-time number of conserved variables.
static const int I4
Energy equation index (= dim + 1).
tCellScalarList tAdditionalCellScalarList
Type alias for additional output scalar list.
TEval::TDof TDof
Cell-centered DOF array type.
PrintDataMode
Output mode selector for PrintData.
@ PrintDataTimeAverage
Output the time-averaged solution.
@ PrintDataLatest
Output the current (latest) solution.
void PrintRestart(std::string fname)
Write the current solution state to a restart file.
auto & getJSource()
Get mutable reference to the source Jacobian.
auto & getIfUseLimiter()
Get mutable reference to the limiter flag array.
auto getVFV() const
Get shared pointer to the VFV reconstruction.
Per-cell diagonal or block-diagonal Jacobian storage for implicit time stepping.
Generic object pool: caches unique_ptr<T> instances and hands them out with RAII return-on-destructio...
DNDS::ArrayPair< DNDS::ArrayEigenVector< 1 > > tScalarPair
Definition VRDefines.hpp:87
EnumValuesTag enum_values(std::vector< std::string > vals)
Create an enum allowed-values tag.
RangeTag range(double min)
Create a minimum-only range constraint.
std::vector< std::tuple< std::string, const tCellScalarFGet > > tCellScalarList
List of (field_name, getter_function) pairs for cell-scalar output.
Definition EulerBC.hpp:706
std::function< real(index)> tCellScalarFGet
Function type returning a scalar value for a given cell index [0, NumCell).
Definition EulerBC.hpp:704
nlohmann::ordered_json bcSettingsSchema()
Build a JSON Schema (draft-07) for the bcSettings array.
Definition EulerBC.hpp:96
MPI_int Barrier(MPI_Comm comm)
Wrapper over MPI_Barrier.
Definition MPI.cpp:247
ssp< SerializerBase > SerializerBaseSSP
std::map< std::string, LogSimpleDIValue > tLogSimpleDIValueMap
Convenience alias: the title -> value map type passed to CsvLog::WriteLine.
Definition CsvLog.hpp:141
DNDS_CONSTANT const rowsize DynamicSize
Template parameter flag: "row width is set at runtime but uniform".
Definition Defines.hpp:282
void from_json(const nlohmann::ordered_json &j, host_device_vector< real > &v)
Definition JsonUtil.hpp:183
int64_t index
Global row / DOF index type (signed 64-bit; handles multi-billion-cell meshes).
Definition Defines.hpp:112
std::shared_ptr< T > ssp
Shortened alias for std::shared_ptr used pervasively in DNDSR.
Definition Defines.hpp:143
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
void to_json(nlohmann::ordered_json &j, const host_device_vector< real > &v)
Definition JsonUtil.hpp:177
DNDS_CONSTANT const real veryLargeReal
Catch-all upper bound ("practically infinity") for physical scalars.
Definition Defines.hpp:195
std::string GetSetVersionName(const std::string &ver)
Read/set the build version string accessible from code.
Definition Defines.cpp:175
void WriteSerialize(Serializer::SerializerBaseSSP serializerP, const std::string &name, bool includePIG=true, bool includeSon=true)
Writes the ArrayPair (father, optional son, optional ghost mapping).
A means to translate nlohmann json into c++ primitive data types and back; and stores then during com...
Master configuration struct for the compressible Euler/Navier-Stokes evaluator.
Geom::UnstructuredMesh::WallDistOptions meshWallDistOptions
Geom::UnstructuredMeshSerialRW::PartitionOptions meshPartitionOptions
std::map< std::string, CoarseGridLinearSolverControl > coarseGridLinearSolverControlList
Miscellaneous solver options (axisymmetric mode, passive scalar freezing, rec matrix output).
Serializer::SerializerFactory recMatrixWriter
Time-averaging control for unsteady simulations.
Complete solver configuration, serializable to/from JSON.
CFV::VRSettings vfvSettings
Variational reconstruction settings.
struct DNDS::Euler::EulerSolver::Configuration::DataIOControl dataIOControl
struct DNDS::Euler::EulerSolver::Configuration::ImplicitReconstructionControl implicitReconstructionControl
std::map< std::string, std::string > bcNameMapping
Mapping from mesh BC names to solver BC type names.
EulerEvaluatorSettings< model > eulerSettings
Physics settings passed to the EulerEvaluator.
struct DNDS::Euler::EulerSolver::Configuration::TimeAverageControl timeAverageControl
struct DNDS::Euler::EulerSolver::Configuration::Others others
void ReadWriteJson(nlohmann::ordered_json &jsonObj, int nVars, bool read)
Backward-compatible bidirectional JSON read/write.
struct DNDS::Euler::EulerSolver::Configuration::LinearSolverControl linearSolverControl
struct DNDS::Euler::EulerSolver::Configuration::TimeMarchControl timeMarchControl
struct DNDS::Euler::EulerSolver::Configuration::LimiterControl limiterControl
struct DNDS::Euler::EulerSolver::Configuration::BoundaryDefinition boundaryDefinition
struct DNDS::Euler::EulerSolver::Configuration::RestartState restartState
struct DNDS::Euler::EulerSolver::Configuration::ConvergenceControl convergenceControl
struct DNDS::Euler::EulerSolver::Configuration::OutputControl outputControl
nlohmann::ordered_json bcSettings
Boundary condition definitions (JSON array).
struct DNDS::Euler::EulerSolver::Configuration::ImplicitCFLControl implicitCFLControl
Mutable state bundle for the time-marching loop.
std::tuple< std::unique_ptr< CsvLog >, tLogSimpleDIValueMap > logErr
ssp< ODE::ImplicitDualTimeStep< ArrayDOFV< nVarsFixed >, ArrayDOFV< 1 > > > ode
Eigen::VectorFMTSafe< real, -1 > resBaseC
Eigen::VectorFMTSafe< real, -1 > resBaseCInternal
std::unique_ptr< tGMRES_uRec > gmresRec
std::map< std::string, ScalarStatistics > tInternalStats
Lightweight bundle of an MPI communicator and the calling rank's coordinates.
Definition MPI.hpp:231
Config-backed factory selecting between JSON and HDF5 serializers.
Eigen::Matrix wrapper that hides begin/end from fmt.
Definition EigenUtil.hpp:66
DistributedHex3D mesh