DNDSR 0.2.1
Distributed Numeric Data Structure for CFV
Loading...
Searching...
No Matches
EulerEvaluatorSettings.hpp
Go to the documentation of this file.
1/** @file EulerEvaluatorSettings.hpp
2 * @brief Complete solver configuration for the Euler/Navier-Stokes evaluator.
3 *
4 * Contains the EulerEvaluatorSettings template struct which aggregates all
5 * runtime parameters for the compressible flow solver:
6 * - Jacobian options (scalar vs. Roe, wall treatment).
7 * - Reconstruction and limiting parameters.
8 * - Riemann solver selection and tuning.
9 * - Wall-distance computation settings.
10 * - RANS turbulence model selection (SA, k-omega) and DES length scales.
11 * - Viscous flux and source term options.
12 * - Rotating reference frame (FrameConstRotation).
13 * - CL (lift-coefficient) driver configuration.
14 * - Region-based initial conditions (box, plane, exprtk).
15 * - Ideal gas thermodynamic properties.
16 *
17 * All settings use DNDS_DECLARE_CONFIG for automatic JSON schema generation
18 * and serialization/deserialization.
19 */
20#pragma once
23
24#include "Euler.hpp"
25#include "Gas.hpp"
26#include "CLDriver.hpp"
27#include <unordered_set>
28
29namespace DNDS::Euler
30{
31 /**
32 * @brief Master configuration struct for the compressible Euler/Navier-Stokes evaluator.
33 *
34 * Organizes all solver-tunable parameters into a single struct that supports
35 * JSON round-trip serialization via DNDS_DECLARE_CONFIG. After deserialization,
36 * the finalize() hook validates cross-field constraints (e.g. mutually exclusive
37 * Jacobian modes) and computes derived reference quantities.
38 *
39 * @tparam model The EulerModel tag (e.g. NS_SA_3D) that determines variable
40 * count, spatial dimension, and available turbulence model traits.
41 */
42 template <EulerModel model>
44 {
45 using Traits = EulerModelTraits<model>; ///< Compile-time model traits.
46 static const int nVarsFixed = getnVarsFixed(model); ///< Compile-time variable count.
47 static const int dim = getDim_Fixed(model); ///< Physical dimension (2 or 3).
48 static const int gDim = getGeomDim_Fixed(model); ///< Geometric dimension (may differ for axi-symmetric).
49 static const auto I4 = dim + 1; ///< Index of the energy equation in the state vector.
50
51 /// @name Jacobian Options
52 /// @{
53 bool useScalarJacobian = false; ///< Use scalar (diagonal) Jacobian approximation instead of block.
54 bool useRoeJacobian = false; ///< Use Roe-linearization-based Jacobian.
55 bool noRsOnWall = false; ///< Disable the Riemann solver on wall boundary faces.
56 bool noGRPOnWall = false; ///< Disable the Generalized Riemann Problem (GRP) on wall faces.
57 bool ignoreSourceTerm = false; ///< Completely ignore source terms (must be false when RANS or body forces are active).
58 /// @}
59
60 /// @name Reconstruction
61 /// @{
62 int direct2ndRecMethod = 1; ///< Direct 2nd-order reconstruction method selector.
63 int specialBuiltinInitializer = 0; ///< Index of a built-in special initializer (0 = none).
64 real uRecAlphaCompressPower = 1; ///< Alpha compression power for the reconstruction limiter.
65 real uRecBetaCompressPower = 1; ///< Beta compression power for the reconstruction limiter.
66 bool forceVolURecBeta = true; ///< Force volume-based beta in the reconstruction.
67 bool ppEpsIsRelaxed = false; ///< Use relaxed positivity-preserving epsilon.
68 /// @}
69
70 real RANSBottomLimit = 0.01; ///< Lower clamp for RANS turbulence variables.
71
72 /// @name Riemann Solver Configuration
73 /// @{
74 Gas::RiemannSolverType rsType = Gas::Roe; ///< Primary Riemann solver type.
75 Gas::RiemannSolverType rsTypeAux = Gas::UnknownRS; ///< Auxiliary Riemann solver type (UnknownRS = same as primary).
76 Gas::RiemannSolverType rsTypeWall = Gas::UnknownRS; ///< Wall-face Riemann solver type (UnknownRS = same as primary).
77 real rsFixScale = 1; ///< Entropy-fix scaling factor for the Riemann solver.
78 real rsIncFScale = 1; ///< Incremental flux scaling factor.
79 int rsMeanValueEig = 0; ///< Mean-value eigenvalue computation mode.
80 int rsRotateScheme = 0; ///< Riemann solver rotation scheme selector.
81 /// @}
82
83 /// @name Wall-Distance Computation
84 /// @{
85 real minWallDist = 1e-12; ///< Minimum wall distance clamp (avoids singularities).
86 int wallDistExection = 0; ///< Execution mode: 0 = parallel, 1 = serial.
87 real wallDistRefineMax = 1; ///< Maximum wall-distance refinement factor.
88 int wallDistScheme = 0; ///< Wall-distance computation scheme selector.
89 int wallDistCellLoadSize = 1024 * 32; ///< Cell batch size for wall-distance computation.
90 int wallDistIter = 1000; ///< Maximum iterations for the wall-distance solver.
91 int wallDistLinSolver = 0; ///< Linear solver: 0 = Jacobi, 1 = GMRES.
92 real wallDistResTol = 1e-4; ///< Residual tolerance for wall-distance convergence.
93 int wallDistIterStart = 100; ///< Starting iteration count for the wall-distance solver.
94 int wallDistPoissonP = 2; ///< Poisson equation power in the wall-distance PDE.
95 real wallDistDTauScale = 100.; ///< Pseudo-time step scaling for wall-distance solver.
96 int wallDistNJacobiSweep = 10; ///< Number of Jacobi sweeps per wall-distance iteration.
97 /// @}
98
99 /// @name RANS / DES Configuration
100 /// @{
101 real SADESScale = veryLargeReal; ///< SA-DES length scale (veryLargeReal effectively disables DES).
102 int SADESMode = 1; ///< SA-DES mode selector (1 = DDES, etc.).
103 /**
104 * @brief SA model variant selector.
105 *
106 * - **0 (default):** Current formulation — rotation correction uses cRot = 2.0
107 * with corrected strain-rate magnitude |S| = ||S_ij + S_ij^T|| / sqrt(2),
108 * SRotCor = cRot * min(0, |Omega| - |S|), and the implicit Jacobian source
109 * includes the negative part of production: min(P, 0) * 1.
110 * - **1 (legacy):** Pre-31578ce (dev/harry_ba3) formulation — rotation correction
111 * uses cRot = 1.0 with the Frobenius norm SS = ||S_ij + S_ij^T||,
112 * SRotCor = cRot * min(0, SS - S), and the implicit Jacobian source omits
113 * production entirely (P * 0).
114 */
115 int SAVersion = 0;
116 RANSModel ransModel = RANSModel::RANS_None; ///< RANS turbulence model (RANS_None, RANS_SA, RANS_KOWilcox, etc.).
117 int ransUseQCR = 0; ///< Enable QCR (Quadratic Constitutive Relation) correction.
118 int ransSARotCorrection = 1; ///< SA rotation/curvature correction mode.
119 int ransEigScheme = 0; ///< Eigenvalue computation scheme for RANS.
120 int ransForce2nd = 0; ///< Force 2nd-order accuracy for RANS variables.
121 int ransSource2nd = 0; ///< Enable 2nd-order RANS source term discretization.
122 /// @}
123
124 /// @name Viscous Flux and Source Options
125 /// @{
126 int source2nd = 0; ///< Enable 2nd-order source term discretization.
127 int usePrimGradInVisFlux = 0; ///< Use primitive-variable gradients in viscous flux.
128 int useSourceGradFixGG = 0; ///< Apply Green-Gauss gradient fix for source terms.
129 int nCentralSmoothStep = 0; ///< Number of central-difference smoothing steps.
130 real centralSmoothEps = 0.5; ///< Epsilon for central smoothing.
131 Eigen::Vector<real, 3> constMassForce = Eigen::Vector<real, 3>{0, 0, 0}; ///< Constant body force vector [fx, fy, fz].
132 /// @}
133 /**
134 * @brief Constant-rotation reference frame settings.
135 *
136 * When enabled, the solver transforms the governing equations into
137 * a non-inertial frame rotating at a constant angular velocity about
138 * the specified axis through the specified center.
139 */
141 {
142 bool enabled = false; ///< Enable the rotating frame.
143 Geom::tPoint axis = Geom::tPoint{0, 0, 1}; ///< Rotation axis (unit vector; normalized in finalize()).
144 Geom::tPoint center = Geom::tPoint{0, 0, 0}; ///< Center of rotation [x, y, z].
145 real rpm = 0; ///< Rotational speed in revolutions per minute.
146
147 /// @brief Compute angular velocity magnitude (rad/s) from RPM.
148 /// @return Omega = rpm * 2π / 60.
150 {
151 return rpm * (2 * pi / 60.);
152 }
153
154 /// @brief Return the angular velocity vector (axis * Omega).
155 /// @return 3-D omega vector aligned with the rotation axis.
157 {
158 return axis * Omega();
159 }
160
161 /**
162 * @brief Project a position vector onto the plane perpendicular to the axis.
163 * @param r Position vector in the absolute frame.
164 * @return Radial component of @p r (axis-normal projection).
165 */
167 {
168 return r - r.dot(axis) * axis;
169 }
170
171 /**
172 * @brief Build the local cylindrical (r, θ, z) coordinate frame at position @p r.
173 *
174 * Column 0 = radial unit vector, column 1 = tangential (axis × r̂),
175 * column 2 = axial (same as the rotation axis).
176 *
177 * @param r Position vector in the absolute frame.
178 * @return 3×3 matrix whose columns are the (r, θ, z) basis vectors.
179 */
181 {
182 Geom::tPoint rn = rVec(r).normalized();
183 Geom::tGPoint ret;
184 ret(EigenAll, 0) = rn;
185 ret(EigenAll, 2) = axis;
186 ret(EigenAll, 1) = axis.cross(rn);
187 return ret;
188 }
190 {
191 // clang-format off
192 DNDS_FIELD(enabled, "Enable constant-rotation reference frame");
193 DNDS_FIELD(axis, "Rotation axis (unit vector)");
194 DNDS_FIELD(center, "Rotation center coordinates");
195 DNDS_FIELD(rpm, "Rotational speed in RPM");
196 // clang-format on
197 }
198 } frameConstRotation; ///< Rotating reference frame configuration.
199 CLDriverSettings cLDriverSettings; ///< Lift-coefficient (CL) driver settings.
200 std::vector<std::string> cLDriverBCNames; ///< Boundary zone names for CL driver force integration.
201 Eigen::Vector<real, -1> farFieldStaticValue = Eigen::Vector<real, 5>{1, 0, 0, 0, 2.5}; ///< Far-field reference state vector (size = nVars).
202 /**
203 * @brief Axis-aligned box region for initial condition specification.
204 *
205 * Cells whose centroids lie within [x0,x1]×[y0,y1]×[z0,z1] are
206 * initialized to the state vector @c v.
207 */
209 {
210 real x0{0}, x1{0}, y0{0}, y1{0}, z0{0}, z1{0}; ///< Box bounds [min, max] per axis.
211 Eigen::Vector<real, -1> v; ///< Initial state vector (size = nVars).
212
214 {
215 // clang-format off
216 DNDS_FIELD(x0, "Box x-min");
217 DNDS_FIELD(x1, "Box x-max");
218 DNDS_FIELD(y0, "Box y-min");
219 DNDS_FIELD(y1, "Box y-max");
220 DNDS_FIELD(z0, "Box z-min");
221 DNDS_FIELD(z1, "Box z-max");
222 DNDS_FIELD(v, "Initial value vector (size = nVars)");
223 // clang-format on
224 }
225 };
226 std::vector<BoxInitializer> boxInitializers; ///< List of box-region initial condition specifiers.
227
228 /**
229 * @brief Half-space region for initial condition specification.
230 *
231 * Cells satisfying `a*x + b*y + c*z >= h` are initialized to the
232 * state vector @c v. The normal direction is (a, b, c).
233 */
235 {
236 real a{0}, b{0}, c{0}, h{0}; ///< Plane equation coefficients: a*x + b*y + c*z = h.
237 Eigen::Vector<real, -1> v; ///< Initial state vector (size = nVars).
238
240 {
241 // clang-format off
242 DNDS_FIELD(a, "Plane normal x-component");
243 DNDS_FIELD(b, "Plane normal y-component");
244 DNDS_FIELD(c, "Plane normal z-component");
245 DNDS_FIELD(h, "Plane offset");
246 DNDS_FIELD(v, "Initial value vector (size = nVars)");
247 // clang-format on
248 }
249 };
250 std::vector<PlaneInitializer> planeInitializers; ///< List of plane-region initial condition specifiers.
251
252 /**
253 * @brief Expression-based initial condition using the ExprTk library.
254 *
255 * Evaluates user-supplied mathematical expressions (one per line in @c exprs)
256 * to compute the initial state at each cell centroid. Lines are concatenated
257 * with newlines to form a single ExprTk program string.
258 */
260 {
261 std::vector<std::string> exprs; ///< ExprTk expression lines (concatenated with newlines).
262
264 {
265 // clang-format off
266 DNDS_FIELD(exprs, "Expression lines (concatenated with newlines)");
267 // clang-format on
268 }
269
270 /**
271 * @brief Concatenate all expression lines into a single ExprTk program string.
272 * @return The full expression string with newline separators.
273 */
274 std::string GetExpr() const
275 {
276 std::string ret;
277 for (auto &line : exprs)
278 ret += line + "\n";
279 return ret;
280 }
281 };
282 std::vector<ExprtkInitializer> exprtkInitializers; ///< List of ExprTk-based initial condition specifiers.
283
284 /**
285 * @brief Ideal gas thermodynamic property set.
286 *
287 * Stores gamma, gas constant, viscosity parameters, and Prandtl number.
288 * The heat capacity CpGas is a derived quantity recomputed automatically
289 * after deserialization via the post_read hook calling recomputeDerived().
290 */
292 {
293 real gamma = 1.4; ///< Ratio of specific heats (Cp/Cv).
294 real Rgas = 1; ///< Specific gas constant (J/(kg·K) in dimensional runs).
295 real muGas = 1; ///< Dynamic viscosity (or reference viscosity for Sutherland).
296 real prGas = 0.72; ///< Prandtl number.
297 real CpGas = Rgas * gamma / (gamma - 1); ///< Heat capacity at constant pressure (derived, not serialized).
298 real TRef = 273.15; ///< Reference temperature (K) for Sutherland's law.
299 real CSutherland = 110.4; ///< Sutherland constant (K).
300 int muModel = 1; ///< Viscosity model: 0 = constant, 1 = Sutherland, 2 = constant_nu.
301
303 {
304 // clang-format off
305 DNDS_FIELD(gamma, "Ratio of specific heats",
307 DNDS_FIELD(Rgas, "Specific gas constant",
309 DNDS_FIELD(muGas, "Dynamic viscosity",
311 DNDS_FIELD(prGas, "Prandtl number",
313 DNDS_FIELD(TRef, "Reference temperature (K)");
314 DNDS_FIELD(CSutherland, "Sutherland constant (K)");
315 DNDS_FIELD(muModel, "Viscosity model: 0=constant, 1=sutherland, 2=constant_nu");
316 // CpGas is derived: recomputed after deserialization
317 config.post_read([](T &s) { s.recomputeDerived(); });
318 // clang-format on
319 }
320
321 /// @brief Recompute derived quantities (CpGas) from gamma and Rgas after deserialization.
323 {
324 CpGas = Rgas * gamma / (gamma - 1);
325 }
326 } idealGasProperty; ///< Ideal gas thermodynamic property configuration.
327
328 /***************************************************************************************************/
329 // end of setting entries
330 /***************************************************************************************************/
331
332 int _nVars = 0; ///< Runtime nVars, not serialized. Set by ctor, preserved across from_json.
333 Eigen::Vector<real, -1> refU; ///< Reference conservative state (derived from farFieldStaticValue).
334 Eigen::Vector<real, -1> refUPrim; ///< Reference primitive state (derived from farFieldStaticValue).
335
337 {
338 // clang-format off
339 DNDS_FIELD(useScalarJacobian, "Use scalar Jacobian approximation");
340 DNDS_FIELD(useRoeJacobian, "Use Roe-based Jacobian");
341 DNDS_FIELD(noRsOnWall, "Disable Riemann solver on wall boundaries");
342 DNDS_FIELD(noGRPOnWall, "Disable GRP on wall boundaries");
343 DNDS_FIELD(ignoreSourceTerm, "Ignore source terms");
344 DNDS_FIELD(direct2ndRecMethod, "Direct 2nd-order reconstruction method");
345 DNDS_FIELD(specialBuiltinInitializer, "Special built-in initializer code");
346 DNDS_FIELD(uRecAlphaCompressPower, "uRec alpha compression power");
347 DNDS_FIELD(uRecBetaCompressPower, "uRec beta compression power");
348 DNDS_FIELD(forceVolURecBeta, "Force volume uRec beta");
349 DNDS_FIELD(ppEpsIsRelaxed, "Positivity-preserving epsilon is relaxed");
350 DNDS_FIELD(RANSBottomLimit, "RANS variable bottom limit",
352 config.field_alias(&T::rsType, "riemannSolverType",
353 "Riemann solver type");
354 config.field_alias(&T::rsTypeAux, "riemannSolverTypeAux",
355 "Auxiliary Riemann solver type");
356 config.field_alias(&T::rsTypeWall, "riemannSolverTypeWall",
357 "Wall Riemann solver type");
358 DNDS_FIELD(rsFixScale, "Riemann solver entropy fix scale");
359 DNDS_FIELD(rsIncFScale, "Riemann solver increment flux scale");
360 DNDS_FIELD(rsMeanValueEig, "Riemann solver mean-value eigenvalue mode");
361 DNDS_FIELD(rsRotateScheme, "Riemann solver rotation scheme");
362 DNDS_FIELD(minWallDist, "Minimum wall distance clamp",
364 DNDS_FIELD(wallDistExection, "Wall distance execution mode: 0=parallel, 1=serial");
365 DNDS_FIELD(wallDistRefineMax, "Wall distance max refinement");
366 DNDS_FIELD(wallDistScheme, "Wall distance computation scheme");
367 DNDS_FIELD(wallDistCellLoadSize, "Wall distance cell load batch size",
369 DNDS_FIELD(wallDistIter, "Wall distance solver iterations",
371 DNDS_FIELD(wallDistLinSolver, "Wall distance linear solver: 0=jacobi, 1=gmres");
372 DNDS_FIELD(wallDistResTol, "Wall distance residual tolerance",
374 DNDS_FIELD(wallDistIterStart, "Wall distance solver start iteration",
376 DNDS_FIELD(wallDistPoissonP, "Wall distance Poisson equation power");
377 DNDS_FIELD(wallDistDTauScale, "Wall distance pseudo-time scale",
379 DNDS_FIELD(wallDistNJacobiSweep, "Wall distance Jacobi sweep count",
381 DNDS_FIELD(SADESScale, "SA-DES length scale");
382 DNDS_FIELD(SADESMode, "SA-DES mode");
383 DNDS_FIELD(SAVersion, "SA variant: 0=current (cRot=2, corrected |S|, min(P,0) Jacobian), 1=legacy/pre-31578ce (cRot=1, Frobenius SS, P*0 Jacobian)");
384 DNDS_FIELD(ransModel, "RANS turbulence model");
385 DNDS_FIELD(ransUseQCR, "Use QCR correction for RANS");
386 DNDS_FIELD(ransSARotCorrection, "SA rotation correction");
387 DNDS_FIELD(ransEigScheme, "RANS eigenvalue scheme");
388 DNDS_FIELD(ransForce2nd, "Force 2nd-order RANS");
389 DNDS_FIELD(ransSource2nd, "RANS source 2nd-order");
390 DNDS_FIELD(source2nd, "Source term 2nd-order");
391 DNDS_FIELD(usePrimGradInVisFlux, "Use primitive gradient in viscous flux");
392 DNDS_FIELD(useSourceGradFixGG, "Use source gradient fix for Green-Gauss");
393 DNDS_FIELD(nCentralSmoothStep, "Central smoothing steps",
395 DNDS_FIELD(centralSmoothEps, "Central smoothing epsilon");
396 DNDS_FIELD(constMassForce, "Constant mass force vector (3D)");
397 config.field_section(&T::frameConstRotation, "frameConstRotation",
398 "Constant-rotation reference frame settings");
399 config.field_section(&T::cLDriverSettings, "cLDriverSettings",
400 "CL driver settings");
401 DNDS_FIELD(cLDriverBCNames, "Boundary names for CL driver force integration");
402 DNDS_FIELD(farFieldStaticValue, "Far-field static value vector (size = nVars)");
403 config.template field_array_of<BoxInitializer>(
404 &T::boxInitializers, "boxInitializers",
405 "Box region initializers");
406 config.template field_array_of<PlaneInitializer>(
407 &T::planeInitializers, "planeInitializers",
408 "Plane region initializers");
409 config.template field_array_of<ExprtkInitializer>(
410 &T::exprtkInitializers, "exprtkInitializers",
411 "Exprtk expression initializers");
412 config.field_section(&T::idealGasProperty, "idealGasProperty",
413 "Ideal gas thermodynamic properties");
414
415 // Cross-field checks
416 config.check("useScalarJacobian and useRoeJacobian are mutually exclusive",
417 [](const T &s) { return !(s.useScalarJacobian && s.useRoeJacobian); });
418 config.check("ransModel must not be RANS_Unknown",
419 [](const T &s) { return s.ransModel != RANS_Unknown; });
420
421 // Post-read hook: finalize derived quantities using stored _nVars
422 config.post_read([](T &s) { s.finalize(); });
423 // clang-format on
424 }
425
426 /// @brief Default constructor (used for schema emission; _nVars remains 0).
428
429 /**
430 * @brief Construct with a known variable count and set model-appropriate defaults.
431 *
432 * If the model includes SA or 2-equation RANS traits, the default ransModel
433 * is set accordingly. The farFieldStaticValue is sized to @p nVars and
434 * initialized to a default freestream state.
435 *
436 * @param nVars Number of conservative variables for this model.
437 */
438 EulerEvaluatorSettings(int nVars) : _nVars(nVars)
439 {
440 if constexpr (Traits::hasSA)
441 {
443 }
444 if constexpr (Traits::has2EQ)
445 {
447 }
448 farFieldStaticValue.setOnes(nVars);
449 DNDS_assert(nVars > I4);
450 farFieldStaticValue(0) = 1;
451 farFieldStaticValue(Eigen::seq(Eigen::fix<0>, Eigen::fix<I4>)).setZero();
452 farFieldStaticValue(I4) = 2.5;
453 }
454
455 /**
456 * @brief Post-deserialization finalization: cross-field validation and derived quantities.
457 *
458 * Checks dimensional consistency of farFieldStaticValue, boxInitializers, and
459 * planeInitializers against _nVars. Normalizes the rotation axis if the rotating
460 * frame is enabled. Computes refU and refUPrim from the far-field state and ideal
461 * gas properties for use as reference scales in the solver.
462 *
463 * Uses the stored _nVars set by the constructor. Called automatically by the
464 * post_read hook after from_json, or explicitly after copy-construction.
465 * If _nVars <= 0 (e.g. default-constructed for schema emission), this is a no-op.
466 */
467 void finalize()
468 {
469 int nVars = _nVars;
470 if (nVars <= 0)
471 return; // skip finalize if nVars not set (e.g. schema emission default-ctor)
472 DNDS_assert(constMassForce.size() == 3);
473 DNDS_assert(farFieldStaticValue.size() == nVars);
475 std::unordered_set<EulerModel>{NS_SA, NS_SA_3D, NS_2EQ, NS_2EQ_3D}.count(model))
477 "you have set source term, do not use ignoreSourceTerm! ");
479 frameConstRotation.axis.normalize();
480 for (auto &box : boxInitializers)
481 DNDS_assert_info(box.v.size() == nVars, "box initial value dimension incorrect");
482 for (auto &plane : planeInitializers)
483 DNDS_assert_info(plane.v.size() == nVars, "plane initial value dimension incorrect");
484
485 // Compute reference values
488 refUPrim = refU;
489 Gas::IdealGasThermalConservative2Primitive<dim>(refU, refUPrim, idealGasProperty.gamma);
490 DNDS_assert(refUPrim(I4) > 0 && refUPrim(0) > 0);
491 real a = std::sqrt(idealGasProperty.gamma * refUPrim(I4) / (refUPrim(0) + verySmallReal));
492 refU(Seq123).setConstant(refU(Seq123).norm() + a);
493 refUPrim(Seq123).setConstant(refUPrim(Seq123).norm());
494 }
495 };
496}
Iterative angle-of-attack driver for targeting a desired lift coefficient.
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.
#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
Core type definitions, enumerations, and distributed array wrappers for compressible Navier-Stokes / ...
#define DNDS_FV_EULEREVALUATOR_GET_FIXED_EIGEN_SEQS
Declares compile-time Eigen index sequences for sub-vector access into conservative state vectors.
Definition Euler.hpp:56
Ideal-gas Riemann solvers, flux functions, and thermodynamic utilities for the compressible Euler / N...
JSON-to-Eigen conversion utilities and nlohmann_json helper macros.
RangeTag range(double min)
Create a minimum-only range constraint.
RiemannSolverType
Selects the approximate Riemann solver and its entropy-fix variant.
Definition Gas.hpp:62
RANSModel
Enumerates the available RANS turbulence closure models.
Definition Euler.hpp:898
@ RANS_Unknown
Sentinel / uninitialized value.
Definition Euler.hpp:899
@ RANS_SA
Spalart-Allmaras one-equation model.
Definition Euler.hpp:901
@ RANS_None
No turbulence model (laminar or DNS).
Definition Euler.hpp:900
@ RANS_KOWilcox
Wilcox k-omega two-equation model.
Definition Euler.hpp:902
Eigen::Vector3d tPoint
Definition Geometric.hpp:9
Eigen::Matrix3d tGPoint
Definition Geometric.hpp:11
DNDS_CONSTANT const real verySmallReal
Catch-all lower bound ("effectively zero").
Definition Defines.hpp:199
DNDS_CONSTANT const real pi
π in double precision (matches DNDS_E_PI macro).
Definition Defines.hpp:204
double real
Canonical floating-point scalar used throughout DNDSR (double precision).
Definition Defines.hpp:110
DNDS_CONSTANT const real veryLargeReal
Catch-all upper bound ("practically infinity") for physical scalars.
Definition Defines.hpp:195
JSON-configurable settings for the CL (lift coefficient) driver.
Definition CLDriver.hpp:31
Axis-aligned box region for initial condition specification.
Eigen::Vector< real, -1 > v
Initial state vector (size = nVars).
Expression-based initial condition using the ExprTk library.
std::vector< std::string > exprs
ExprTk expression lines (concatenated with newlines).
std::string GetExpr() const
Concatenate all expression lines into a single ExprTk program string.
real Omega()
Compute angular velocity magnitude (rad/s) from RPM.
real rpm
Rotational speed in revolutions per minute.
Geom::tPoint vOmega()
Return the angular velocity vector (axis * Omega).
Geom::tPoint axis
Rotation axis (unit vector; normalized in finalize()).
Geom::tGPoint rtzFrame(const Geom::tPoint &r)
Build the local cylindrical (r, θ, z) coordinate frame at position r.
Geom::tPoint rVec(const Geom::tPoint &r)
Project a position vector onto the plane perpendicular to the axis.
int muModel
Viscosity model: 0 = constant, 1 = Sutherland, 2 = constant_nu.
real muGas
Dynamic viscosity (or reference viscosity for Sutherland).
real CpGas
Heat capacity at constant pressure (derived, not serialized).
real Rgas
Specific gas constant (J/(kg·K) in dimensional runs).
real TRef
Reference temperature (K) for Sutherland's law.
void recomputeDerived()
Recompute derived quantities (CpGas) from gamma and Rgas after deserialization.
Half-space region for initial condition specification.
Eigen::Vector< real, -1 > v
Initial state vector (size = nVars).
real h
Plane equation coefficients: a*x + b*y + c*z = h.
Master configuration struct for the compressible Euler/Navier-Stokes evaluator.
int ransForce2nd
Force 2nd-order accuracy for RANS variables.
real wallDistResTol
Residual tolerance for wall-distance convergence.
bool ignoreSourceTerm
Completely ignore source terms (must be false when RANS or body forces are active).
real centralSmoothEps
Epsilon for central smoothing.
bool forceVolURecBeta
Force volume-based beta in the reconstruction.
real uRecAlphaCompressPower
Alpha compression power for the reconstruction limiter.
int ransSource2nd
Enable 2nd-order RANS source term discretization.
int wallDistExection
Execution mode: 0 = parallel, 1 = serial.
int wallDistNJacobiSweep
Number of Jacobi sweeps per wall-distance iteration.
int rsRotateScheme
Riemann solver rotation scheme selector.
int ransEigScheme
Eigenvalue computation scheme for RANS.
real rsIncFScale
Incremental flux scaling factor.
int rsMeanValueEig
Mean-value eigenvalue computation mode.
int ransSARotCorrection
SA rotation/curvature correction mode.
Gas::RiemannSolverType rsTypeAux
Auxiliary Riemann solver type (UnknownRS = same as primary).
static const int gDim
Geometric dimension (may differ for axi-symmetric).
int wallDistIterStart
Starting iteration count for the wall-distance solver.
std::vector< std::string > cLDriverBCNames
Boundary zone names for CL driver force integration.
real wallDistRefineMax
Maximum wall-distance refinement factor.
int _nVars
Runtime nVars, not serialized. Set by ctor, preserved across from_json.
void finalize()
Post-deserialization finalization: cross-field validation and derived quantities.
int wallDistPoissonP
Poisson equation power in the wall-distance PDE.
real SADESScale
SA-DES length scale (veryLargeReal effectively disables DES).
int SADESMode
SA-DES mode selector (1 = DDES, etc.).
real wallDistDTauScale
Pseudo-time step scaling for wall-distance solver.
int wallDistLinSolver
Linear solver: 0 = Jacobi, 1 = GMRES.
Gas::RiemannSolverType rsTypeWall
Wall-face Riemann solver type (UnknownRS = same as primary).
Eigen::Vector< real, -1 > refUPrim
Reference primitive state (derived from farFieldStaticValue).
static const auto I4
Index of the energy equation in the state vector.
Eigen::Vector< real, -1 > refU
Reference conservative state (derived from farFieldStaticValue).
std::vector< BoxInitializer > boxInitializers
List of box-region initial condition specifiers.
int direct2ndRecMethod
Direct 2nd-order reconstruction method selector.
Gas::RiemannSolverType rsType
Primary Riemann solver type.
int wallDistIter
Maximum iterations for the wall-distance solver.
int source2nd
Enable 2nd-order source term discretization.
int usePrimGradInVisFlux
Use primitive-variable gradients in viscous flux.
real minWallDist
Minimum wall distance clamp (avoids singularities).
bool useRoeJacobian
Use Roe-linearization-based Jacobian.
int wallDistScheme
Wall-distance computation scheme selector.
real uRecBetaCompressPower
Beta compression power for the reconstruction limiter.
real RANSBottomLimit
Lower clamp for RANS turbulence variables.
int nCentralSmoothStep
Number of central-difference smoothing steps.
CLDriverSettings cLDriverSettings
Lift-coefficient (CL) driver settings.
struct DNDS::Euler::EulerEvaluatorSettings::FrameConstRotation frameConstRotation
Rotating reference frame configuration.
EulerEvaluatorSettings(int nVars)
Construct with a known variable count and set model-appropriate defaults.
Eigen::Vector< real, 3 > constMassForce
Constant body force vector [fx, fy, fz].
real rsFixScale
Entropy-fix scaling factor for the Riemann solver.
std::vector< PlaneInitializer > planeInitializers
List of plane-region initial condition specifiers.
int specialBuiltinInitializer
Index of a built-in special initializer (0 = none).
EulerEvaluatorSettings()=default
Default constructor (used for schema emission; _nVars remains 0).
bool ppEpsIsRelaxed
Use relaxed positivity-preserving epsilon.
int ransUseQCR
Enable QCR (Quadratic Constitutive Relation) correction.
bool noRsOnWall
Disable the Riemann solver on wall boundary faces.
int wallDistCellLoadSize
Cell batch size for wall-distance computation.
RANSModel ransModel
RANS turbulence model (RANS_None, RANS_SA, RANS_KOWilcox, etc.).
bool noGRPOnWall
Disable the Generalized Riemann Problem (GRP) on wall faces.
std::vector< ExprtkInitializer > exprtkInitializers
List of ExprTk-based initial condition specifiers.
Eigen::Vector< real, -1 > farFieldStaticValue
Far-field reference state vector (size = nVars).
bool useScalarJacobian
Use scalar (diagonal) Jacobian approximation instead of block.
static const int dim
Physical dimension (2 or 3).
int useSourceGradFixGG
Apply Green-Gauss gradient fix for source terms.
struct DNDS::Euler::EulerEvaluatorSettings::IdealGasProperty idealGasProperty
Ideal gas thermodynamic property configuration.
static const int nVarsFixed
Compile-time variable count.
Compile-time traits for EulerModel variants.
Definition Euler.hpp:1086
static constexpr bool hasSA
True for Spalart-Allmaras models (NS_SA, NS_SA_3D).
Definition Euler.hpp:1095
static constexpr bool has2EQ
True for 2-equation RANS models (NS_2EQ, NS_2EQ_3D).
Definition Euler.hpp:1097
tVec r(NCells)