DNDSR 0.1.0.dev1+gcd065ad
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
38#include "DNDS/JsonUtil.hpp"
39#include "DNDS/ConfigParam.hpp"
41#include "DNDS/CsvLog.hpp"
42#include "DNDS/ObjectPool.hpp"
43#include "Solver/Linear.hpp"
44#include "Geom/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
53#include "DNDS/JsonUtil.hpp"
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 DNDS_FIELD(dtImplicit, "Max implicit time step; 1e100 for steady",
199 DNDS_FIELD(dtImplicitMin, "Minimum implicit time step",
201 DNDS_FIELD(nTimeStep, "Max number of time steps",
203 DNDS_FIELD(steadyQuit, "Quit on steady convergence");
204 DNDS_FIELD(useRestart, "Initialize from restart file");
205 DNDS_FIELD(useImplicitPP, "Use implicit positivity-preserving");
206 DNDS_FIELD(rhsFPPMode, "RHS flux-positivity-preserving mode");
207 DNDS_FIELD(rhsFPPScale, "RHS FPP scaling factor");
208 DNDS_FIELD(rhsFPPRelax, "RHS FPP relaxation factor",
209 DNDS::Config::range(0.0, 1.0));
210 DNDS_FIELD(incrementPPRelax, "Increment PP relaxation factor",
211 DNDS::Config::range(0.0, 1.0));
212 DNDS_FIELD(odeCode, "ODE integrator code");
213 DNDS_FIELD(tEnd, "End time for unsteady simulation");
214 DNDS_FIELD(odeSetting1, "ODE parameter 1");
215 DNDS_FIELD(odeSetting2, "ODE parameter 2");
216 DNDS_FIELD(odeSetting3, "ODE parameter 3");
217 DNDS_FIELD(odeSetting4, "ODE parameter 4");
218 config.field_json(&T::odeSettingsExtra, "odeSettingsExtra", "Extra ODE integrator settings (opaque JSON)");
219 DNDS_FIELD(partitionMeshOnly, "Only partition mesh, then exit");
220 DNDS_FIELD(dtIncreaseLimit, "Max factor for dt increase per step",
222 DNDS_FIELD(dtIncreaseAfterCount,"Steps before dt increase allowed",
224 DNDS_FIELD(dtCFLLimitScale, "CFL-based dt limit scale");
225 DNDS_FIELD(useDtPPLimit, "Use positivity-preserving dt limiter");
226 DNDS_FIELD(dtPPLimitRelax, "PP dt limiter relaxation",
227 DNDS::Config::range(0.0, 1.0));
228 DNDS_FIELD(dtPPLimitScale, "PP dt limiter scale",
230 }
232 {
233 return odeCode == 401 || (odeCode >= 411 && odeCode <= 413);
234 }
236
237 /**
238 * @brief Implicit reconstruction control parameters.
239 *
240 * Controls the iterative reconstruction solve within each time step:
241 * explicit vs implicit reconstruction, linear solver (SOR or GMRES/PCG),
242 * convergence thresholds, and gradient zeroing strategies.
243 */
245 {
246 bool useExplicit = false;
248 bool zeroGrads = false;
249 int recLinearScheme = 0; // 0 for original SOR, 1 for GMRES
250 int nGmresSpace = 5;
251 int nGmresIter = 10;
260 bool storeRecInc = false;
261 bool dampRecIncDTau = false;
265 {
266 DNDS_FIELD(useExplicit, "Use explicit reconstruction (no implicit)");
267 DNDS_FIELD(nInternalRecStep, "Number of internal reconstruction sub-steps",
269 DNDS_FIELD(zeroGrads, "Zero gradients before reconstruction");
270 DNDS_FIELD(recLinearScheme, "Reconstruction linear solver: 0=SOR, 1=GMRES");
271 DNDS_FIELD(nGmresSpace, "GMRES Krylov subspace size for reconstruction",
273 DNDS_FIELD(nGmresIter, "GMRES iterations for reconstruction",
275 DNDS_FIELD(gmresRecScale, "GMRES reconstruction scaling mode");
276 DNDS_FIELD(fpcgResetScheme, "FPCG reset scheme");
277 DNDS_FIELD(fpcgResetThres, "FPCG reset threshold",
278 DNDS::Config::range(0.0, 1.0));
279 DNDS_FIELD(fpcgResetReport, "FPCG reset report frequency");
280 DNDS_FIELD(fpcgMaxPHistory, "FPCG max preconditioner history",
282 DNDS_FIELD(recThreshold, "Reconstruction convergence threshold",
284 DNDS_FIELD(nRecConsolCheck, "Console output interval for reconstruction",
286 DNDS_FIELD(nRecMultiplyForZeroedGrad, "Rec iteration multiplier when grad is zeroed",
288 DNDS_FIELD(storeRecInc, "Store reconstruction increment");
289 DNDS_FIELD(dampRecIncDTau, "Damp reconstruction increment by dTau");
290 DNDS_FIELD(zeroRecForSteps, "Zero reconstruction for N outer steps",
292 DNDS_FIELD(zeroRecForStepsInternal, "Zero reconstruction for N internal steps",
294 }
296
297 /**
298 * @brief Output control parameters.
299 *
300 * Controls console output frequency and formatting, data file output
301 * intervals (VTK/HDF5), restart checkpoint intervals, time-average output,
302 * and logging precision.
303 */
305 {
308 int consoleOutputMode = 0; // 0 for basic, 1 for wall force out
311 std::vector<std::string> consoleMainOutputFormat{
312 "=== Step {termBold}[{step:4d}] ",
313 "res {termBold}{termRed}{resRel:.3e}{termReset} ",
314 "t,dT,dTaumin,CFL,nFix {termGreen}[{tSimu:.3e},{curDtImplicit:.3e},{curDtMin:.3e},{CFLNow:.3e},[alphaInc({nLimInc},{alphaMinInc}), betaRec({nLimBeta},{minBeta}), alphaRes({nLimAlpha},{minAlpha})]]{termReset} ",
315 "Time[{telapsed:.3f}] recTime[{trec:.3f}] rhsTime[{trhs:.3f}] commTime[{tcomm:.3f}] limTime[{tLim:.3f}] limTimeA[{tLimiterA:.3f}] limTimeB[{tLimiterB:.3f}]"};
316 std::vector<std::string> consoleMainOutputFormatInternal{
317 "\t Internal === Step [{step:4d},{iStep:2d},{iter:4d}] ",
318 "res {termRed}{resRel:.3e}{termReset} ",
319 "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} ",
320 "Time[{telapsedM:.3f}] recTime[{trecM:.3f}] rhsTime[{trhsM:.3f}] commTime[{tcommM:.3f}] limTime[{tLimM:.3f}] limTimeA[{tLimiterA:.3f}] limTimeB[{tLimiterB:.3f}]"};
321 std::vector<std::string> logfileOutputTitles{
322 "step", "iStep", "iterAll", "iter", "tSimu",
323 "res", "curDtImplicit", "curDtMin", "CFLNow",
324 "nLimInc", "alphaMinInc",
325 "nLimBeta", "minBeta",
326 "nLimAlpha", "minAlpha",
327 "tWall", "telapsed", "trec", "trhs", "tcomm", "tLim", "tLimiterA", "tLimiterB",
328 "fluxWall", "CL", "CD", "AoA"};
330 bool dataOutAtInit = false;
331 bool restartOutAtInit = false;
332 int nDataOut = 10000;
333 int nDataOutC = 50;
334 int nDataOutInternal = 10000;
336 int nRestartOut = INT_MAX;
337 int nRestartOutC = INT_MAX;
338 int nRestartOutInternal = INT_MAX;
339 int nRestartOutCInternal = INT_MAX;
340 int nTimeAverageOut = INT_MAX;
341 int nTimeAverageOutC = INT_MAX;
344 bool useCollectiveTimer = false;
345
347 {
348 DNDS_FIELD(nConsoleCheck, "Console output interval (outer steps)",
350 DNDS_FIELD(nConsoleCheckInternal, "Console output interval (internal steps)",
352 DNDS_FIELD(consoleOutputMode, "Console output mode: 0=basic, 1=wall force");
353 DNDS_FIELD(consoleOutputEveryFix, "Force console output every N fix iterations");
354 // nPrecisionConsole intentionally not registered: was not serialized
355 // in the original macro and existing config files lack this key.
356 DNDS_FIELD(consoleMainOutputFormat, "Format strings for console output lines");
357 DNDS_FIELD(consoleMainOutputFormatInternal, "Format strings for internal-step console output");
358 DNDS_FIELD(logfileOutputTitles, "Column titles for CSV log file");
359 DNDS_FIELD(nPrecisionLog, "Log file floating-point precision",
361 DNDS_FIELD(dataOutAtInit, "Output data at initialization");
362 DNDS_FIELD(restartOutAtInit, "Output restart at initialization");
363 DNDS_FIELD(nDataOut, "Data output interval (outer steps)",
365 DNDS_FIELD(nDataOutC, "Data output interval (convergence steps)",
367 DNDS_FIELD(nDataOutInternal, "Data output interval (internal steps)",
369 DNDS_FIELD(nDataOutCInternal, "Data output interval (internal convergence)");
370 DNDS_FIELD(nRestartOut, "Restart output interval (outer steps)");
371 DNDS_FIELD(nRestartOutC, "Restart output interval (convergence steps)");
372 DNDS_FIELD(nRestartOutInternal, "Restart output interval (internal steps)");
373 DNDS_FIELD(nRestartOutCInternal, "Restart output interval (internal convergence)");
374 DNDS_FIELD(nTimeAverageOut, "Time-average output interval (outer steps)");
375 DNDS_FIELD(nTimeAverageOutC, "Time-average output interval (convergence steps)");
376 DNDS_FIELD(tDataOut, "Output data at simulation time interval");
377 DNDS_FIELD(lazyCoverDataOutput, "Overwrite previous data output files");
378 DNDS_FIELD(useCollectiveTimer, "Use collective MPI timer for profiling");
379 }
381
382 /**
383 * @brief Implicit CFL number control parameters.
384 *
385 * Controls the CFL number, local vs global time stepping, CFL ramping
386 * schedule, local dTau smoothing, and RANS under-relaxation.
387 */
389 {
390 real CFL = 10;
391 int nForceLocalStartStep = INT_MAX;
392 int nCFLRampStart = INT_MAX;
393 int nCFLRampLength = INT_MAX;
395 bool useLocalDt = true;
396 int nSmoothDTau = 0;
399 {
400 DNDS_FIELD(CFL, "CFL number for implicit time stepping",
402 DNDS_FIELD(nForceLocalStartStep, "Step to force local time stepping",
404 DNDS_FIELD(nCFLRampStart, "Step to begin CFL ramping",
406 DNDS_FIELD(nCFLRampLength, "Number of steps for CFL ramp",
408 DNDS_FIELD(CFLRampEnd, "CFL value at end of ramp");
409 DNDS_FIELD(useLocalDt, "Use local (vs uniform) time step");
410 DNDS_FIELD(nSmoothDTau, "Smoothing passes for local dTau",
412 DNDS_FIELD(RANSRelax, "RANS equation under-relaxation factor",
413 DNDS::Config::range(0.0, 1.0));
414 }
416
417 /**
418 * @brief Convergence monitoring parameters.
419 *
420 * Controls inner-loop iteration counts, convergence thresholds,
421 * residual norm type, and CL-driver (lift-coefficient) adaptation.
422 */
424 {
431 int resBaseType = 0; // 1 to use rhs as base in unsteady
432 int mergeMultiResidual = 0; // for HM3, merge stage residuals
433 int normOrd = 1; // use LN norm
434 bool useVolWiseResidual = false;
435 bool useCLDriver = false;
437 {
438 DNDS_FIELD(nTimeStepInternal, "Max internal iterations per time step (0 = unlimited)",
440 DNDS_FIELD(nTimeStepInternalMin, "Min internal iterations per time step",
442 DNDS_FIELD(nAnchorUpdate, "Anchor update frequency",
444 DNDS_FIELD(nAnchorUpdateStart, "Step to start anchor updates",
446 DNDS_FIELD(rhsThresholdInternal, "RHS convergence threshold for internal loop",
448 DNDS_FIELD(res_base, "Residual base value (0=auto)");
449 DNDS_FIELD(resBaseType, "Residual base type: 0=default, 1=rhs-based for unsteady");
450 DNDS_FIELD(mergeMultiResidual, "Merge stage residuals for multi-stage: 0=off");
451 DNDS_FIELD(normOrd, "Residual norm order (1=L1, 2=L2, etc.)",
453 DNDS_FIELD(useVolWiseResidual, "Volume-weighted residual");
454 DNDS_FIELD(useCLDriver, "Enable CL-driven AoA adaptation");
455 }
457
458 /**
459 * @brief Data I/O control parameters.
460 *
461 * Controls mesh file paths and formats, mesh preprocessing (elevation,
462 * bisection, wall distance), output formats (Tecplot, VTK, HDF5),
463 * restart serialization, and coordinate transformations.
464 */
466 {
467 bool uniqueStamps = true;
470
471 int meshElevation = 0; // 0 = noOp, 1 = O1->O2
472 int meshElevationInternalSmoother = 0; // 0 = local interpolation, 1 = coupled
473 int meshElevationIter = 1000; // -1 to use handle all the nodes
480 int meshElevationBoundaryMode = 0; // 0: only wall bc; 1: wall + invis vall
481
483 int meshReorderCells = 0; // 0: natural; 1: reorder
486
487 int meshFormat = 0;
489 std::string meshFile = "data/mesh/NACA0012_WIDE_H3.cgns";
490 std::string outPltName = "data/out/debugData_";
491 std::string outLogName = "";
492 std::string outRestartName = "";
493
494 int outPltMode = 0; // 0 = serial, 1 = dist plt
495 int readMeshMode = 0; // 0 = serial cgns, 1 = dist json
497 bool outPltVTKFormat = false;
498 bool outPltVTKHDFFormat = false;
499 bool outAtPointData = true;
500 bool outAtCellData = true;
502 std::string vtuFloatEncodeMode = "binary";
503 int hdfChunkSize = 32768;
505 bool hdfCollOnData = false;
506 bool hdfCollOnMeta = true;
507 bool outVolumeData = true;
508 bool outBndData = false;
509
510 std::vector<std::string> outCellScalarNames{};
511
512 bool serializerSaveURec = false;
513
515
516 int rectifyNearPlane = 0; // 1: x 2: y 4: z
518
521 std::string meshPartitionedReaderType = "JSON";
522
523 const std::string &getOutLogName()
524 {
525 return outLogName.empty() ? outPltName : outLogName;
526 }
527
528 const std::string &getOutRestartName()
529 {
530 return outRestartName.empty() ? outPltName : outRestartName;
531 }
532
534 {
535 DNDS_FIELD(uniqueStamps, "Use unique output stamps per run");
536 DNDS_FIELD(meshRotZ, "Mesh rotation around Z axis (degrees)");
537 DNDS_FIELD(meshScale, "Mesh coordinate scaling factor",
539 DNDS_FIELD(meshElevation, "Mesh order elevation: 0=none, 1=O1->O2");
540 DNDS_FIELD(meshElevationInternalSmoother,"Elevation smoother: 0=local, 1=coupled");
541 DNDS_FIELD(meshElevationIter, "Elevation smoothing iterations");
542 DNDS_FIELD(meshElevationNSearch, "Elevation RBF neighbor search count",
544 DNDS_FIELD(meshElevationRBFRadius, "Elevation RBF support radius",
546 DNDS_FIELD(meshElevationRBFPower, "Elevation RBF power parameter");
547 DNDS_FIELD(meshElevationRBFKernel, "Elevation RBF kernel type");
548 DNDS_FIELD(meshElevationMaxIncludedAngle,"Max included angle for elevation (degrees)");
549 DNDS_FIELD(meshElevationRefDWall, "Reference wall distance for elevation",
551 DNDS_FIELD(meshElevationBoundaryMode, "Elevation boundary: 0=wall only, 1=wall+inviscid");
552 DNDS_FIELD(meshDirectBisect, "Mesh bisection passes",
554 DNDS_FIELD(meshReorderCells, "Reorder cells: 0=natural, 1=reorder");
555 DNDS_FIELD(meshBuildWallDist, "Build wall distance field: 0=no, 1=yes");
556 DNDS_FIELD(meshWallDistOptions, "Wall distance computation options");
557 DNDS_FIELD(meshFormat, "Mesh format code");
558 DNDS_FIELD(meshPartitionOptions, "Mesh partitioning options");
559 DNDS_FIELD(meshFile, "Input mesh file path");
560 DNDS_FIELD(outPltName, "Output plot file base name");
561 DNDS_FIELD(outLogName, "Output log file base name (empty=use outPltName)");
562 DNDS_FIELD(outRestartName, "Output restart file base name (empty=use outPltName)");
563 DNDS_FIELD(outPltMode, "Output mode: 0=serial, 1=distributed");
564 DNDS_FIELD(readMeshMode, "Read mesh mode: 0=serial CGNS, 1=distributed JSON");
565 DNDS_FIELD(outPltTecplotFormat, "Output in Tecplot format");
566 DNDS_FIELD(outPltVTKFormat, "Output in VTK XML format");
567 DNDS_FIELD(outPltVTKHDFFormat, "Output in VTK HDF format");
568 DNDS_FIELD(outAtPointData, "Output point-interpolated data");
569 DNDS_FIELD(outAtCellData, "Output cell-centered data");
570 DNDS_FIELD(nASCIIPrecision, "ASCII output floating-point precision",
572 DNDS_FIELD(vtuFloatEncodeMode, "VTU float encoding: binary or ascii",
573 DNDS::Config::enum_values({"binary", "ascii"}));
574 DNDS_FIELD(hdfChunkSize, "HDF5 chunk size for output",
576 DNDS_FIELD(hdfDeflateLevel, "HDF5 deflate compression level",
577 DNDS::Config::range(0, 9));
578 DNDS_FIELD(hdfCollOnMeta, "HDF5 collective I/O on metadata");
579 DNDS_FIELD(hdfCollOnData, "HDF5 collective I/O on data arrays");
580 DNDS_FIELD(outVolumeData, "Output volume data");
581 DNDS_FIELD(outBndData, "Output boundary data");
582 DNDS_FIELD(outCellScalarNames, "Additional cell scalar names to output");
583 DNDS_FIELD(serializerSaveURec, "Save reconstruction in restart");
584 DNDS_FIELD(allowAsyncPrintData, "Allow asynchronous data output");
585 DNDS_FIELD(rectifyNearPlane, "Rectify nodes near planes: bitmask 1=x,2=y,4=z");
586 DNDS_FIELD(rectifyNearPlaneThres, "Threshold for near-plane rectification",
588 config.field_section(&T::restartWriter, "restartWriter", "Restart file serializer settings");
589 config.field_section(&T::meshPartitionedWriter, "meshPartitionedWriter", "Partitioned mesh serializer settings");
590 DNDS_FIELD(meshPartitionedReaderType, "Partitioned mesh reader type",
591 DNDS::Config::enum_values({"JSON", "H5"}));
592 }
594
595 /**
596 * @brief Periodic boundary geometry definitions.
597 *
598 * Defines up to 3 periodic translation vectors and rotation parameters
599 * for periodic boundary matching.
600 */
602 {
603 Eigen::Vector<real, -1> PeriodicTranslation1;
604 Eigen::Vector<real, -1> PeriodicTranslation2;
605 Eigen::Vector<real, -1> PeriodicTranslation3;
606 Eigen::Vector<real, -1> PeriodicRotationCent1;
607 Eigen::Vector<real, -1> PeriodicRotationCent2;
608 Eigen::Vector<real, -1> PeriodicRotationCent3;
614 {
615 PeriodicTranslation1.resize(3);
616 PeriodicTranslation2.resize(3);
617 PeriodicTranslation3.resize(3);
618 PeriodicTranslation1 << 1, 0, 0;
619 PeriodicTranslation2 << 0, 1, 0;
620 PeriodicTranslation3 << 0, 0, 1;
621 PeriodicRotationCent1.setZero(3);
622 PeriodicRotationCent2.setZero(3);
623 PeriodicRotationCent3.setZero(3);
627 }
628
630 {
631 DNDS_FIELD(PeriodicTranslation1, "Periodic translation vector for pair 1");
632 DNDS_FIELD(PeriodicTranslation2, "Periodic translation vector for pair 2");
633 DNDS_FIELD(PeriodicTranslation3, "Periodic translation vector for pair 3");
634 DNDS_FIELD(PeriodicRotationCent1, "Rotation center for periodic pair 1");
635 DNDS_FIELD(PeriodicRotationCent2, "Rotation center for periodic pair 2");
636 DNDS_FIELD(PeriodicRotationCent3, "Rotation center for periodic pair 3");
637 DNDS_FIELD(PeriodicRotationEulerAngles1, "Rotation Euler angles (deg) for periodic pair 1");
638 DNDS_FIELD(PeriodicRotationEulerAngles2, "Rotation Euler angles (deg) for periodic pair 2");
639 DNDS_FIELD(PeriodicRotationEulerAngles3, "Rotation Euler angles (deg) for periodic pair 3");
640 DNDS_FIELD(periodicTolerance, "Tolerance for periodic node matching",
642 }
644
645 /**
646 * @brief Slope limiter control parameters.
647 *
648 * Controls whether limiters are active, which limiter variant to use
649 * (WBAP, CWBAP), positivity-preserving reconstruction limiting, and
650 * partial/preserved limiting strategies.
651 */
653 {
654 bool useLimiter = true;
655 bool usePPRecLimiter = true;
656 bool useViscousLimited = true;
659 int nPartialLimiterStart = INT_MAX;
661 bool preserveLimited = false;
663
665 {
666 DNDS_FIELD(useLimiter, "Enable slope limiter");
667 DNDS_FIELD(usePPRecLimiter, "Enable positivity-preserving reconstruction limiter");
668 DNDS_FIELD(useViscousLimited, "Apply limiter to viscous reconstruction");
669 DNDS_FIELD(smoothIndicatorProcedure, "Smooth indicator procedure index");
670 DNDS_FIELD(limiterProcedure, "Limiter variant: 0=WBAP (V2), 1=CWBAP (V3)");
671 DNDS_FIELD(nPartialLimiterStart, "Time step to begin partial limiting");
672 DNDS_FIELD(nPartialLimiterStartLocal, "Time step to begin local partial limiting");
673 DNDS_FIELD(preserveLimited, "Preserve limited reconstruction across steps");
674 DNDS_FIELD(ppRecLimiterCompressToMean, "PP limiter compresses toward cell mean");
675 }
677
678 /**
679 * @brief Linear solver control parameters.
680 *
681 * Controls the implicit linear solver: preconditioner type (Jacobi/GS/ILU),
682 * SGS iterations, GMRES settings, multi-grid options, and per-level
683 * coarse-grid solver configurations.
684 */
686 {
687 int jacobiCode = 1; // 0 for jacobi, 1 for gs, 2 for ilu
688 int sgsIter = 0;
689 int sgsWithRec = 0;
690 int gmresCode = 0; // 0 for lusgs, 1 for gmres, 2 for lusgs started gmres
691 int gmresScale = 0; // 0 for no scaling, 1 use refU, 2 use mean value
692 int nGmresSpace = 10;
693 int nGmresIter = 2;
697 int multiGridLP = 0;
702 {
703 int jacobiCode = 0; // 0 for jacobi, 1 for gs, 2 for ilu
704 int sgsIter = 0;
705 int gmresCode = 0; // 0 for lusgs, 1 for gmres, 2 for lusgs started gmres
706 int gmresScale = 0; // 0 for no scaling, 1 use refU, 2 use mean value
707 int nGmresIter = 2;
713
715 {
716 DNDS_FIELD(jacobiCode, "Preconditioner: 0=jacobi, 1=GS, 2=ILU");
717 DNDS_FIELD(sgsIter, "SGS iteration count",
719 DNDS_FIELD(gmresCode, "Linear solver: 0=LUSGS, 1=GMRES, 2=LUSGS+GMRES");
720 DNDS_FIELD(gmresScale, "GMRES scaling: 0=none, 1=refU, 2=mean");
721 DNDS_FIELD(nGmresIter, "GMRES iterations",
723 DNDS_FIELD(nSgsConsoleCheck, "Console check interval for SGS",
725 DNDS_FIELD(nGmresConsoleCheck, "Console check interval for GMRES",
727 DNDS_FIELD(multiGridNIter, "Multi-grid pre-smooth iterations (-1=auto)");
728 DNDS_FIELD(multiGridNIterPost, "Multi-grid post-smooth iterations",
730 DNDS_FIELD(centralSmoothInputResidual, "Smooth input residual on coarse grid");
731 }
732 };
733 std::map<std::string, CoarseGridLinearSolverControl> coarseGridLinearSolverControlList{
736 };
738
740 {
741 DNDS_FIELD(jacobiCode, "Preconditioner: 0=jacobi, 1=GS, 2=ILU");
742 DNDS_FIELD(sgsIter, "SGS iteration count",
744 DNDS_FIELD(sgsWithRec, "Couple SGS with reconstruction");
745 DNDS_FIELD(gmresCode, "Linear solver: 0=LUSGS, 1=GMRES, 2=LUSGS+GMRES");
746 DNDS_FIELD(gmresScale, "GMRES scaling: 0=none, 1=refU, 2=mean");
747 DNDS_FIELD(nGmresSpace, "GMRES Krylov subspace size",
749 DNDS_FIELD(nGmresIter, "GMRES outer iterations",
751 DNDS_FIELD(nSgsConsoleCheck, "Console check interval for SGS",
753 DNDS_FIELD(nGmresConsoleCheck, "Console check interval for GMRES",
755 DNDS_FIELD(initWithLastURecInc, "Initialize GMRES with last uRec increment");
756 DNDS_FIELD(multiGridLP, "Multi-grid levels (0=off)");
757 DNDS_FIELD(multiGridLPInnerNIter, "Multi-grid inner iterations",
759 DNDS_FIELD(multiGridLPStartIter, "Iteration to enable multi-grid",
761 DNDS_FIELD(multiGridLPInnerNSee, "Multi-grid inner console see interval",
763 config.template field_map_of<CoarseGridLinearSolverControl>(
764 &T::coarseGridLinearSolverControlList,
765 "coarseGridLinearSolverControlList",
766 "Per-level coarse grid linear solver settings");
767 config.field_section(&T::directPrecControl, "directPrecControl",
768 "Direct preconditioner settings");
769 }
771
772 /**
773 * @brief Restart checkpoint state.
774 *
775 * Records the step indices and file paths for restart continuation,
776 * including support for reading restart files from a different solver
777 * (with dimension remapping).
778 */
780 {
781 int iStep = -1;
783 int odeCodePrev = -1;
784 std::string lastRestartFile = "";
785 std::string otherRestartFile = "";
786 std::vector<int> otherRestartStoreDim;
788 {
789 DNDS_FIELD(iStep, "Restart step index");
790 DNDS_FIELD(iStepInternal, "Restart internal step index");
791 DNDS_FIELD(odeCodePrev, "Previous ODE code for restart");
792 DNDS_FIELD(lastRestartFile, "Path to last restart file");
793 DNDS_FIELD(otherRestartFile, "Path to alternate restart file");
794 DNDS_FIELD(otherRestartStoreDim, "Dimension mapping for alternate restart");
795 }
797 {
798 otherRestartStoreDim.resize(1);
799 for (int i = 0; i < otherRestartStoreDim.size(); i++)
801 }
803
804 /// @brief Time-averaging control for unsteady simulations.
806 {
807 bool enabled = false;
808
810 {
811 DNDS_FIELD(enabled, "Enable time-averaging of solution fields");
812 }
814
815 /// @brief Miscellaneous solver options (axisymmetric mode, passive scalar freezing, rec matrix output).
816 struct Others
817 {
820 bool printRecMatrix = false;
822
824 {
825 DNDS_FIELD(nFreezePassiveInner, "Freeze passive scalars for N inner steps",
827 DNDS_FIELD(axisSymmetric, "Axisymmetric mode: 0=off");
828 DNDS_FIELD(printRecMatrix, "Print reconstruction matrix to file");
829 config.field_section(&T::recMatrixWriter, "recMatrixWriter",
830 "Serializer for reconstruction matrix output");
831 }
833
834 EulerEvaluatorSettings<model> eulerSettings; ///< Physics settings passed to the EulerEvaluator.
835 CFV::VRSettings vfvSettings; ///< Variational reconstruction settings.
836 nlohmann::ordered_json bcSettings = nlohmann::ordered_json::array(); ///< Boundary condition definitions (JSON array).
837 std::map<std::string, std::string> bcNameMapping; ///< Mapping from mesh BC names to solver BC type names.
838
840 {
841 config.field_section(&T::timeMarchControl, "timeMarchControl", "Time marching settings");
842 config.field_section(&T::implicitReconstructionControl, "implicitReconstructionControl", "Implicit reconstruction settings");
843 config.field_section(&T::outputControl, "outputControl", "Output settings");
844 config.field_section(&T::implicitCFLControl, "implicitCFLControl", "Implicit CFL settings");
845 config.field_section(&T::convergenceControl, "convergenceControl", "Convergence monitoring settings");
846 config.field_section(&T::dataIOControl, "dataIOControl", "Data I/O and restart settings");
847 config.field_section(&T::boundaryDefinition, "boundaryDefinition", "Periodic boundary geometry");
848 config.field_section(&T::limiterControl, "limiterControl", "Slope limiter settings");
849 config.field_section(&T::linearSolverControl, "linearSolverControl", "Linear solver settings");
850 config.field_section(&T::timeAverageControl, "timeAverageControl", "Time averaging settings");
851 config.field_section(&T::others, "others", "Miscellaneous settings");
852 config.field_section(&T::restartState, "restartState", "Restart checkpoint state");
853 config.field_section(&T::eulerSettings, "eulerSettings", "Euler evaluator physics settings");
854 config.field_section(&T::vfvSettings, "vfvSettings", "Variational reconstruction settings");
855 config.field_json_schema(&T::bcSettings, "bcSettings",
856 "Boundary condition settings (per-BC array)",
857 []() { return bcSettingsSchema(); });
858 DNDS_FIELD(bcNameMapping, "Boundary name to type mapping");
859
860 config.check("bcSettings must be a JSON array", [](const T &s)
861 {
862 return s.bcSettings.is_array();
863 });
864 }
865
866 /// @brief Backward-compatible bidirectional JSON read/write.
867 ///
868 /// Delegates to the auto-generated from_json / to_json from
869 /// DNDS_DECLARE_CONFIG. Kept for call-site compatibility.
870 void
871 ReadWriteJson(nlohmann::ordered_json &jsonObj, int nVars, bool read)
872 {
873 (void)nVars; // nVars is already stored in eulerSettings via Configuration(nVars)
874 if (read)
875 from_json(jsonObj, *this);
876 else
877 to_json(jsonObj, *this);
878 }
879
881 {
882 }
883
884 Configuration(int nVars)
886 {
888 }
889
890 } config = Configuration{};
891
892 public:
893 /**
894 * @brief Construct an EulerSolver with MPI communicator information.
895 *
896 * Initializes the runtime variable count, output field sizes, and default
897 * configuration. Does not read mesh or allocate solver arrays.
898 *
899 * @param nmpi MPI communicator information.
900 * @param n_nVars Runtime number of conserved variables (default from model).
901 */
902 EulerSolver(const MPIInfo &nmpi, int n_nVars = getNVars(model)) : nVars(n_nVars), mpi(nmpi)
903 {
904 if (getNVars(model) == DynamicSize)
905 DNDS_assert_info(nVars >= getDim_Fixed(model) + 2, "nVars too small");
906 else
907 DNDS_assert_info(nVars == getNVars(model), "do not change the nVars for this model");
908 nOUTS = nVars + 4;
909 nOUTSPoint = nVars + 2;
910 nOUTSBnd = nVars + 2 + nVars + 3 + 1 + 3; // Uprim + (T,M) + F + Ft + faceZone + Norm
911
912 config = Configuration(nVars); //* important to initialize using nVars
913 }
914
915 /// @brief Destructor. Waits for all async output futures to complete.
917 {
918 int nBad{0};
919 do
920 {
921 nBad = 0;
922 for (auto &f : outFuture)
923 if (f.valid() && f.wait_for(std::chrono::microseconds(10)) != std::future_status::ready)
924 nBad++;
925 for (auto &f : outBndFuture)
926 if (f.valid() && f.wait_for(std::chrono::microseconds(10)) != std::future_status::ready)
927 nBad++;
928 if (outSeqFuture.valid() && outSeqFuture.wait_for(std::chrono::microseconds(10)) != std::future_status::ready)
929 nBad++;
930 } while (nBad);
931 }
932
933 /**
934 * @brief Load or write solver configuration from/to a JSON file.
935 *
936 * When read=true, parses the JSON file, optionally merges a patch file,
937 * applies key/value overrides, populates the Configuration struct and
938 * creates the boundary condition handler. When read=false, serializes the
939 * current configuration to a JSON file.
940 *
941 * @param jsonName Path to the JSON configuration file.
942 * @param read true=read from file, false=write to file.
943 * @param jsonMergeName Optional path to a JSON patch file to merge.
944 * @param overwriteKeys JSON pointer paths to override (e.g., "/timeMarchControl/CFL").
945 * @param overwriteValues Values corresponding to overwriteKeys.
946 */
947 void ConfigureFromJson(const std::string &jsonName, bool read = false, const std::string &jsonMergeName = "",
948 const std::vector<std::string> &overwriteKeys = {}, const std::vector<std::string> &overwriteValues = {})
949 {
950 if (read)
951 {
952 auto fIn = std::ifstream(jsonName);
953 DNDS_assert_info(fIn, "config file not existent");
954 gSetting = nlohmann::ordered_json::parse(fIn, nullptr, true, true);
955
956 if (read && !jsonMergeName.empty())
957 {
958 fIn = std::ifstream(jsonMergeName);
959 DNDS_assert_info(fIn, "config file patch not existent");
960 auto gSettingAdd = nlohmann::ordered_json::parse(fIn, nullptr, true, true);
961 gSetting.merge_patch(gSettingAdd);
962 }
963 DNDS_assert(overwriteKeys.size() == overwriteValues.size());
964 for (size_t i = 0; i < overwriteKeys.size(); i++)
965 {
966
967 auto key = nlohmann::ordered_json::json_pointer(overwriteKeys[i].c_str());
968 try
969 {
970 std::string valString =
971 fmt::format(R"({{
972 "__val_entry": {}
973 }})",
974 overwriteValues[i]);
975 auto valDoc = nlohmann::ordered_json::parse(valString, nullptr, true, true);
976 if (mpi.rank == 0)
977 log() << "JSON: overwrite key: " << key << std::endl
978 << "JSON: overwrite val: " << valDoc["__val_entry"] << std::endl;
979 gSetting[key] = valDoc["__val_entry"];
980 }
981 catch (const nlohmann::ordered_json::parse_error &e)
982 {
983 try
984 {
985 std::string valString =
986 fmt::format(R"({{
987"__val_entry": "{}"
988}})",
989 overwriteValues[i]);
990 auto valDoc = nlohmann::ordered_json::parse(valString, nullptr, true, true);
991 if (mpi.rank == 0)
992 log() << "JSON: overwrite key: " << key << std::endl
993 << "JSON: overwrite val: " << valDoc["__val_entry"] << std::endl;
994 gSetting[key] = valDoc["__val_entry"];
995 }
996 catch (const std::exception &e)
997 {
998 std::cerr << e.what() << "\n";
999 std::cerr << overwriteValues[i] << "\n";
1000 DNDS_assert(false);
1001 }
1002 }
1003 catch (const std::exception &e)
1004 {
1005 std::cerr << e.what() << "\n";
1006 std::cerr << overwriteValues[i] << "\n";
1007 DNDS_assert(false);
1008 }
1009 }
1010 config.ReadWriteJson(gSetting, nVars, read);
1011 // create from json the pBCHandler
1012 pBCHandler = std::make_shared<BoundaryHandler<model>>(nVars);
1013 from_json(config.bcSettings, *pBCHandler);
1014 gSetting["bcSettings"] = *pBCHandler;
1015 PrintConfig(true);
1016 if (mpi.rank == 0)
1017 log() << "JSON: read value:" << std::endl
1018 << std::setw(4) << gSetting << std::endl;
1019 }
1020 else
1021 {
1022 gSetting = nlohmann::ordered_json::object();
1023 config.ReadWriteJson(gSetting, nVars, read);
1024 if (pBCHandler) // todo: add example pBCHandler
1025 gSetting["bcSettings"] = *pBCHandler;
1026 if (mpi.rank == 0) // single call for output
1027 {
1028 std::filesystem::path outFile{jsonName};
1029 std::filesystem::create_directories(outFile.parent_path() / ".");
1030 auto fIn = std::ofstream(jsonName);
1031 DNDS_assert(fIn);
1032 fIn << std::setw(4) << gSetting;
1033 }
1034 MPI::Barrier(mpi.comm); // no go until output done
1035 }
1036
1037 if (mpi.rank == 0)
1038 log() << "JSON: Parse " << (read ? "read" : "write")
1039 << " Done ===" << std::endl;
1040 }
1041
1042 /**
1043 * @brief Read the mesh and initialize the full solver pipeline.
1044 *
1045 * Performs the complete initialization sequence:
1046 * 1. Read mesh (serial CGNS or distributed)
1047 * 2. Apply coordinate transformations (rotation, scaling, rectification)
1048 * 3. Partition the mesh across MPI ranks (Metis/ParMetis)
1049 * 4. Apply order elevation (O1->O2) and mesh bisection
1050 * 5. Build wall distance field (if requested)
1051 * 6. Construct the variational reconstruction (VFV)
1052 * 7. Create the EulerEvaluator and initialize FV infrastructure
1053 * 8. Allocate DOF, reconstruction, and Jacobian arrays
1054 * 9. Load restart data (if configured)
1055 * 10. Set initial conditions (if not restarting)
1056 */
1057 void ReadMeshAndInitialize();
1058
1059 /**
1060 * @brief Write the current configuration to a JSON log file.
1061 *
1062 * Extracts live settings from the VFV, evaluator, and BC handler,
1063 * serializes the full configuration to JSON, and writes it alongside
1064 * compile-time defines and commit information.
1065 *
1066 * @param updateCommit If true, include the git commit hash in the output.
1067 */
1068 void PrintConfig(bool updateCommit = false)
1069 {
1070 /***********************************************************/
1071 // if these objects are existent, extract settings from them
1072 if (vfv)
1073 config.vfvSettings = static_cast<const CFV::VRSettings &>(vfv->getSettings());
1074 if (pEval)
1075 config.eulerSettings = pEval->settings;
1076 if (pBCHandler)
1077 config.bcSettings = *pBCHandler;
1078 /***********************************************************/
1079 config.ReadWriteJson(gSetting, nVars, false);
1080 if (mpi.rank == 0)
1081 {
1082 std::string logConfigFileName = config.dataIOControl.getOutLogName() + "_" + output_stamp + ".config.json";
1083 std::filesystem::path outFile{logConfigFileName};
1084 std::filesystem::create_directories(outFile.parent_path() / ".");
1085 std::ofstream logConfig(logConfigFileName);
1086 DNDS_assert(logConfig);
1087 gSetting["___Compile_Time_Defines"] = DNDS_Defines_state;
1088 gSetting["___Runtime_PartitionNumber"] = mpi.size;
1089 gSetting["___Commit_ID"] = GetSetVersionName();
1090 // if (updateCommit)
1091 // {
1092 // std::ifstream commitIDFile("commitID.txt");
1093 // if (commitIDFile)
1094 // {
1095 // std::string commitHash;
1096 // commitIDFile >> commitHash;
1097 // gSetting["___Commit_ID"] = commitHash;
1098 // }
1099 // }
1100 logConfig << std::setw(4) << gSetting;
1101 logConfig.close();
1102 }
1103 }
1104 /// @brief Read a restart file and populate u (and optionally uRec) from it.
1105 void ReadRestart(std::string fname);
1106
1107 /// @brief Read a restart file from a different solver/model, remapping variable dimensions.
1108 void ReadRestartOtherSolver(std::string fname, const std::vector<int> &dimStore);
1109
1110 /// @brief Write the current solution state to a restart file.
1111 void PrintRestart(std::string fname);
1112
1113 using tAdditionalCellScalarList = tCellScalarList; ///< Type alias for additional output scalar list.
1114
1115 /// @brief Output mode selector for PrintData.
1117 {
1118 PrintDataLatest = 0, ///< Output the current (latest) solution.
1119 PrintDataTimeAverage = 1, ///< Output the time-averaged solution.
1120 };
1121
1122 /**
1123 * @brief Write solution data to VTK/HDF5/Tecplot output files.
1124 *
1125 * Gathers distributed data to serial (or writes distributed), converts
1126 * conservative to primitive variables, appends additional scalars (residual,
1127 * limiter flags), and writes the output in the configured format(s).
1128 *
1129 * @param fname Output file base name.
1130 * @param fnameSeries Series file name (for VTK time series).
1131 * @param odeResidualF Callback returning the ODE residual for each cell.
1132 * @param additionalCellScalars Additional per-cell scalar fields to output.
1133 * @param eval Reference to the evaluator (for output field computation).
1134 * @param TSimu Current simulation time (-1 for steady).
1135 * @param mode PrintDataLatest or PrintDataTimeAverage.
1136 */
1137 void PrintData(const std::string &fname, const std::string &fnameSeries,
1138 const tCellScalarFGet &odeResidualF,
1139 tAdditionalCellScalarList &additionalCellScalars,
1140 TEval &eval, real TSimu = -1.0, PrintDataMode mode = PrintDataLatest);
1141
1142 /// @brief Serialize the solution to a SerializerBase (currently unused standalone path).
1143 void WriteSerializer(Serializer::SerializerBaseSSP serializerP, const std::string &name) // currently not using
1144 {
1145 auto cwd = serializerP->GetCurrentPath();
1146 serializerP->CreatePath(name);
1147 serializerP->GoToPath(name);
1148
1149 u.WriteSerialize(serializerP, "meanValue");
1150
1151 serializerP->GoToPath(cwd);
1152
1153 nlohmann::ordered_json configJson;
1154 config.ReadWriteJson(configJson, nVars, false);
1155 serializerP->WriteString("lastConfig", configJson.dump());
1156
1158 {
1159 serializerP->WriteInt("hasReconstructionValue", 1);
1160 uRec.WriteSerialize(serializerP, "recValue");
1161 }
1162 else
1163 serializerP->WriteInt("hasReconstructionValue", 0);
1164 }
1165
1166 /// @brief Fill a log value map entry for arithmetic types (scalars).
1167 template <class TVal>
1168 std::enable_if_t<std::is_arithmetic_v<std::remove_reference_t<TVal>>>
1169 FillLogValue(tLogSimpleDIValueMap &v_map, const std::string &name, TVal &&val)
1170 {
1171 v_map[name] = val;
1172 }
1173
1174 /// @brief Fill a log value map entry for non-arithmetic types (vectors → per-component entries).
1175 template <class TVal>
1176 std::enable_if_t<!std::is_arithmetic_v<std::remove_reference_t<TVal>>>
1177 FillLogValue(tLogSimpleDIValueMap &v_map, const std::string &name, TVal &&val)
1178 {
1179 // std::vector<std::string> logfileOutputTitles{
1180 // "step", "iStep", "iter", "tSimu",
1181 // "res", "curDtImplicit", "curDtMin", "CFLNow",
1182 // "nLimInc", "alphaMinInc",
1183 // "nLimBeta", "minBeta",
1184 // "nLimAlpha", "minAlpha",
1185 // "tWall", "telapsed", "trec", "trhs", "tcomm", "tLim", "tLimiterA", "tLimiterB",
1186 // "fluxWall", "CL", "CD", "AoA"};
1187 if (name == "res" || name == "fluxWall")
1188 for (int i = 0; i < nVars; i++)
1189 v_map[name + std::to_string(i)] = val[i];
1190 else
1191 v_map[name] = 0;
1192 }
1193
1194 /// @brief Initialize the CSV error logger and value map for convergence monitoring.
1195 /// @return Tuple of (CsvLog writer, value map with all column names initialized).
1196 std::tuple<std::unique_ptr<CsvLog>, tLogSimpleDIValueMap> LogErrInitialize()
1197 {
1199 TU initVec;
1200 initVec.setZero(nVars);
1201 std::vector<std::string> realNames;
1202 for (auto name : config.outputControl.logfileOutputTitles)
1203 if (name == "res" || name == "fluxWall")
1204 {
1205 FillLogValue(v_map, name, initVec);
1206 for (int i = 0; i < nVars; i++)
1207 realNames.push_back(name + std::to_string(i));
1208 }
1209 else
1210 FillLogValue(v_map, name, 0.), realNames.push_back(name);
1211
1212 // std::string logErrFileName = config.dataIOControl.getOutLogName() + "_" + output_stamp + ".log";
1213 // std::filesystem::path outFile{logErrFileName};
1214 // std::filesystem::create_directories(outFile.parent_path() / ".");
1215 // auto pOs = std::make_unique<std::ofstream>(logErrFileName);
1216 // DNDS_assert_info(*pOs, "logErr file [" + logErrFileName + "] did not open");
1217
1218 return std::make_tuple(std::make_unique<DNDS::CsvLog>(
1219 realNames,
1221 ".log",
1222 100000),
1223 v_map);
1224 }
1225
1226 /**
1227 * @brief Run the main implicit time-marching loop.
1228 *
1229 * Executes the outer time-step loop with configurable ODE integrators
1230 * (backward Euler, BDF2, SDIRK, etc.). Each step involves:
1231 * 1. CFL ramping and time-step computation
1232 * 2. Variational reconstruction (implicit or explicit)
1233 * 3. Slope limiting and positivity-preserving enforcement
1234 * 4. RHS evaluation
1235 * 5. Implicit linear solve (LU-SGS, GMRES, or direct)
1236 * 6. Solution update with increment compression
1237 * 7. Convergence monitoring and data output
1238 *
1239 * Supports steady-state convergence checks, CL-driver AoA adaptation,
1240 * time-averaging, and restart checkpointing.
1241 */
1243
1244 /**
1245 * @brief Mutable state bundle for the time-marching loop.
1246 *
1247 * Holds all transient state that persists across time steps within
1248 * RunImplicitEuler: evaluator, ODE integrator, linear solvers, timing
1249 * statistics, residual bases, convergence trackers, output counters,
1250 * and CFL/dt history. The DNDS_EULERSOLVER_RUNNINGENV_GET_REF_LIST
1251 * macro provides convenient local aliases for all members.
1252 */
1254 {
1256 std::tuple<std::unique_ptr<CsvLog>, tLogSimpleDIValueMap> logErr;
1258 std::unique_ptr<tGMRES_u> gmres;
1259 std::unique_ptr<tGMRES_uRec> gmresRec;
1260 std::unique_ptr<tPCG_uRec> pcgRec;
1261 std::unique_ptr<tPCG_uRec> pcgRec1;
1262 double tstart = 0;
1263 double tstartInternal = 0;
1264 std::map<std::string, ScalarStatistics> tInternalStats;
1265 int stepCount = 0;
1271 int nextStepOut = -1;
1277
1279 bool ifOutT = false;
1282 std::vector<real> curDtImplicitHistory;
1283 int step = 0;
1284 int iterAll = 0;
1285 bool gradIsZero = true;
1286
1293
1295
1297
1298#define DNDS_EULERSOLVER_RUNNINGENV_GET_REF(name) auto &name = runningEnvironment.name
1299
1300#define DNDS_EULERSOLVER_RUNNINGENV_GET_REF_LIST \
1301 auto &eval = *runningEnvironment.pEval; \
1302 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(logErr); \
1303 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(ode); \
1304 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(gmres); \
1305 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(gmresRec); \
1306 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(pcgRec); \
1307 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(pcgRec1); \
1308 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(tstart); \
1309 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(tstartInternal); \
1310 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(tInternalStats); \
1311 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(stepCount); \
1312 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(resBaseC); \
1313 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(resBaseCInternal); \
1314 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(tSimu); \
1315 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(tAverage); \
1316 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(nextTout); \
1317 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(nextStepOut); \
1318 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(nextStepOutC); \
1319 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(nextStepRestart); \
1320 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(nextStepRestartC); \
1321 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(nextStepOutAverage); \
1322 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(nextStepOutAverageC); \
1323 \
1324 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(CFLNow); \
1325 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(ifOutT); \
1326 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(curDtMin); \
1327 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(curDtImplicit); \
1328 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(curDtImplicitHistory); \
1329 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(step); \
1330 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(iterAll); \
1331 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(gradIsZero); \
1332 \
1333 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(nLimBeta); \
1334 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(nLimAlpha); \
1335 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(minAlpha); \
1336 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(minBeta); \
1337 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(nLimInc); \
1338 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(alphaMinInc); \
1339 \
1340 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(dtIncreaseCounter); \
1341 \
1342 DNDS_EULERSOLVER_RUNNINGENV_GET_REF(addOutList);
1343
1345 };
1346 /// @brief Populate a RunningEnvironment with allocated solvers, loggers, and initial state.
1348
1349 /**
1350 * @brief Solve the implicit linear system at a given grid level.
1351 *
1352 * Dispatches to LU-SGS, GMRES, or LU-SGS-preconditioned GMRES depending on
1353 * configuration. Used within the ODE integrator's inner loop.
1354 *
1355 * @warning Explicit template instantiation does not exist; inlined only.
1356 */
1358 real alphaDiag, real t,
1359 TDof &cres, TDof &cx, TDof &cxInc, TRec &uRecC, TRec uRecIncC,
1360 JacobianDiagBlock<nVarsFixed> &JDC, tGMRES_u &gmres, int gridLevel);
1361 /**
1362 * @brief Apply the preconditioner (SGS sweeps or ILU) to a right-hand side vector.
1363 *
1364 * @warning Explicit template instantiation does not exist; inlined only.
1365 */
1366 void doPrecondition(real alphaDiag, real t,
1367 TDof &crhs, TDof &cx, TDof &cxInc, TDof &uTemp,
1368 JacobianDiagBlock<nVarsFixed> &JDC, TU &sgsRes, bool &inputIsZero, bool &hasLUDone, int gridLevel);
1369
1370 /// @brief Convergence/termination check functor called after each inner iteration.
1371 /// @return true if the inner loop should stop (converged or max iterations reached).
1372 bool functor_fstop(int iter, ArrayDOFV<nVarsFixed> &cres, int iStep, RunningEnvironment &env);
1373 /// @brief Main outer-loop functor: performs one full time step (reconstruction, RHS, linear solve, update).
1374 /// @return true if the outer loop should continue.
1376
1377 /// @name Accessors
1378 /// @{
1379 auto getMPI() const { return mpi; } ///< Get MPI communicator info.
1380 auto getMesh() const { return mesh; } ///< Get shared pointer to the mesh.
1381 auto getVFV() const { return vfv; } ///< Get shared pointer to the VFV reconstruction.
1382 /// @}
1383
1384 /// @name Test accessors (for unit testing the evaluator pipeline)
1385 /// @{
1386 auto &getU() { return u; } ///< Get mutable reference to the DOF array.
1387 auto &getURec() { return uRec; } ///< Get mutable reference to the reconstruction array.
1388 auto &getURecNew() { return uRecNew; } ///< Get mutable reference to the new reconstruction array.
1389 auto &getURecLimited() { return uRecLimited; } ///< Get mutable reference to the limited reconstruction array.
1390 auto &getEval() { return *pEval; } ///< Get mutable reference to the evaluator.
1391 auto &getConfiguration() { return config; } ///< Get mutable reference to the configuration.
1392 auto &getJSource() { return JSource; } ///< Get mutable reference to the source Jacobian.
1393 auto &getBetaPP() { return betaPP; } ///< Get mutable reference to the PP beta array.
1394 auto &getAlphaPP() { return alphaPP; } ///< Get mutable reference to the PP alpha array.
1395 auto &getDTauTmp() { return dTauTmp; } ///< Get mutable reference to the dTau temporary.
1396 auto &getIfUseLimiter() { return ifUseLimiter; }///< Get mutable reference to the limiter flag array.
1397 auto &getBCHandler() { return pBCHandler; } ///< Get mutable reference to the BC handler.
1398 /// @}
1399 };
1400}
1401
1402#define DNDS_EULERSOLVER_INS_EXTERN(model, ext) \
1403 namespace DNDS::Euler \
1404 { \
1405 ext template void EulerSolver<model>::RunImplicitEuler(); \
1406 ext template void EulerSolver<model>::InitializeRunningEnvironment( \
1407 EulerSolver<model>::RunningEnvironment &env); \
1408 }
1409
1417
1418#define DNDS_EULERSOLVER_PRINTDATA_INS_EXTERN(model, ext) \
1419 namespace DNDS::Euler \
1420 { \
1421 ext template void EulerSolver<model>::PrintData( \
1422 const std::string &fname, const std::string &fnameSeries, \
1423 const tCellScalarFGet &odeResidualF, \
1424 tAdditionalCellScalarList &additionalCellScalars, \
1425 TEval &eval, real tSimu, \
1426 PrintDataMode mode); \
1427 ext template void EulerSolver<model>::PrintRestart( \
1428 std::string fname); \
1429 ext template void EulerSolver<model>::ReadRestart( \
1430 std::string fname); \
1431 ext template void EulerSolver<model>::ReadRestartOtherSolver( \
1432 std::string fname, const std::vector<int> &dimStore); \
1433 }
1434
1442
1443#define DNDS_EULERSOLVER_INIT_INS_EXTERN(model, ext) \
1444 namespace DNDS::Euler \
1445 { \
1446 ext template void EulerSolver<model>::ReadMeshAndInitialize(); \
1447 ext template bool EulerSolver<model>::functor_fstop( \
1448 int iter, ArrayDOFV<nVarsFixed> &cres, int iStep, RunningEnvironment &env); \
1449 ext template bool EulerSolver<model>::functor_fmainloop( \
1450 RunningEnvironment &env); \
1451 }
1452
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:113
#define DNDS_assert(expr)
Debug-only assertion (compiled out when DNDS_NDEBUG is defined). Prints the expression + file/line + ...
Definition Errors.hpp:108
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:248
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:277
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: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
std::ostream & log()
Return the current DNDSR log stream (either std::cout or the installed file).
Definition Defines.cpp:49
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:190
std::string GetSetVersionName(const std::string &ver)
Read/set the build version string accessible from code.
Definition Defines.cpp:172
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:215
Config-backed factory selecting between JSON and HDF5 serializers.
Eigen::Matrix wrapper that hides begin/end from fmt.
Definition EigenUtil.hpp:62