DNDSR 0.1.0.dev1+gcd065ad
Distributed Numeric Data Structure for CFV
Loading...
Searching...
No Matches
DNDS::Euler Namespace Reference

Namespaces

namespace  Gas
 
namespace  RANS
 RANS turbulence model functions (eddy viscosity, viscous flux, source terms).
 
namespace  SpecialFields
 

Classes

struct  AnchorPointRecorder
 Finds the cell closest to a specified anchor point across all MPI ranks. More...
 
class  ArrayDOFV
 MPI-parallel distributed array of per-cell Degrees-of-Freedom (conservative variable vectors). More...
 
class  ArrayGRADV
 MPI-parallel distributed array of per-cell gradient matrices. More...
 
class  ArrayRECV
 MPI-parallel distributed array of per-cell reconstruction coefficient matrices. More...
 
class  BoundaryHandler
 Per-zone boundary condition handler for Euler/Navier-Stokes solvers. More...
 
class  CLDriver
 Lift-coefficient driver controller. More...
 
struct  CLDriverSettings
 JSON-configurable settings for the CL (lift coefficient) driver. More...
 
class  EulerEvaluator
 Core finite-volume evaluator for compressible Navier-Stokes / Euler equations. More...
 
struct  EulerEvaluatorSettings
 Master configuration struct for the compressible Euler/Navier-Stokes evaluator. More...
 
struct  EulerModelTraits
 Compile-time traits for EulerModel variants. More...
 
class  EulerSolver
 Top-level solver orchestrator for compressible Navier-Stokes / Euler equations. More...
 
struct  IntegrationRecorder
 Accumulator for MPI-reduced boundary-face integrals. More...
 
class  JacobianDiagBlock
 Per-cell diagonal or block-diagonal Jacobian storage for implicit time stepping. More...
 
struct  JacobianLocalLDLT
 Symmetric LDLT factorization of the cell-local Jacobian. More...
 
struct  JacobianLocalLU
 Complete LU factorization of the cell-local Jacobian for GMRES preconditioning. More...
 
class  JacobianValue
 Placeholder container for implicit-solver Jacobian matrix storage. More...
 
struct  OneDimProfile
 One-dimensional profile for inlet/outlet boundary data. More...
 
class  OutputPicker
 Registry that maps named cell-scalar output fields to getter lambdas. More...
 
struct  TU_P_Reduction
 

Typedefs

using tCellScalarFGet = std::function< real(index)>
 Function type returning a scalar value for a given cell index [0, NumCell).
 
using tCellScalarList = std::vector< std::tuple< std::string, const tCellScalarFGet > >
 List of (field_name, getter_function) pairs for cell-scalar output.
 

Enumerations

enum  EulerModel {
  NS = 0 , NS_SA = 1 , NS_2D = 2 , NS_3D = 3 ,
  NS_SA_3D = 4 , NS_2EQ = 5 , NS_2EQ_3D = 6 , NS_EX = 101 ,
  NS_EX_3D = 102
}
 Enumerates the available compressible flow solver model configurations. More...
 
enum  RANSModel {
  RANS_Unknown = 0 , RANS_None , RANS_SA , RANS_KOWilcox ,
  RANS_KOSST , RANS_RKE
}
 Enumerates the available RANS turbulence closure models. More...
 
enum  EulerBCType {
  BCUnknown = 0 , BCFar , BCWall , BCWallInvis ,
  BCWallIsothermal , BCOut , BCOutP , BCIn ,
  BCInPsTs , BCSym , BCSpecial
}
 Boundary condition type identifiers for compressible flow solvers. More...
 

Functions

 DNDS_DEFINE_ENUM_JSON (RANSModel, { {RANS_Unknown, nullptr}, {RANS_None, "RANS_None"}, {RANS_SA, "RANS_SA"}, {RANS_KOWilcox, "RANS_KOWilcox"}, {RANS_KOSST, "RANS_KOSST"}, {RANS_RKE, "RANS_RKE"}, }) const expr static inline int getnVarsFixed(const EulerModel model)
 Returns the compile-time (fixed) number of conservative variables for a given model.
 
 DNDS_DEFINE_ENUM_JSON (EulerBCType, { {BCUnknown, nullptr}, {BCFar, "BCFar"}, {BCWall, "BCWall"}, {BCWallIsothermal, "BCWallIsothermal"}, {BCWallInvis, "BCWallInvis"}, {BCOut, "BCOut"}, {BCOutP, "BCOutP"}, {BCIn, "BCIn"}, {BCInPsTs, "BCInPsTs"}, {BCSym, "BCSym"}, {BCSpecial, "BCSpecial"}, }) inline std
 Convert an EulerBCType enumerator to its JSON string representation.
 
nlohmann::ordered_json bcSettingsSchema ()
 Build a JSON Schema (draft-07) for the bcSettings array.
 
 DNDS_SWITCH_INTELLISENSE (template< EulerModel model >,) void EulerEvaluator< model >
 Assemble diagonal Jacobian block for each cell for LU-SGS preconditioning.
 
 DNDS_SWITCH_INTELLISENSE (template< EulerModel model >, template<>) void EulerEvaluator< model >
 [[deprecated]] Forward LU-SGS sweep. Use UpdateSGS with uIncIsZero=true instead.
 
void paste_read_restart_with_cell_ordering (EulerSolver< model > &self, typename EulerSolver< model >::TDof &u, typename EulerSolver< model >::TDof &uRead, Serializer::SerializerBaseSSP serializerP)
 Reorder restart data to match current cell ordering using cell2cellOrig mapping.
 
template<EulerModel model>
int RunSingleBlockConsoleApp (int argc, char *argv[])
 Main entry point for single-block solver console applications.
 

Typedef Documentation

◆ tCellScalarFGet

using DNDS::Euler::tCellScalarFGet = typedef std::function<real( index )>

Function type returning a scalar value for a given cell index [0, NumCell).

Definition at line 702 of file EulerBC.hpp.

◆ tCellScalarList

using DNDS::Euler::tCellScalarList = typedef std::vector<std::tuple<std::string, const tCellScalarFGet> >

List of (field_name, getter_function) pairs for cell-scalar output.

Definition at line 706 of file EulerBC.hpp.

Enumeration Type Documentation

◆ EulerBCType

Boundary condition type identifiers for compressible flow solvers.

Each enumerator selects the physics applied at a boundary zone:

  • Far-field and inflow/outflow conditions use characteristic or prescribed states.
  • Wall conditions enforce no-penetration (inviscid) or no-slip (viscous) constraints.
  • Special conditions implement test-case-specific setups (Riemann problems, DMR, isentropic vortex, Rayleigh–Taylor).
Enumerator
BCUnknown 

Uninitialized / invalid sentinel.

BCFar 

Far-field (characteristic-based) boundary.

BCWall 

No-slip viscous wall boundary.

BCWallInvis 

Inviscid slip wall boundary.

BCWallIsothermal 

Isothermal wall; requires wall temperature in the BC value vector.

BCOut 

Supersonic outflow (extrapolation) boundary.

BCOutP 

Back-pressure (subsonic) outflow boundary.

BCIn 

Supersonic inflow (fully prescribed state) boundary.

BCInPsTs 

Subsonic inflow with prescribed total pressure and temperature.

BCSym 

Symmetry plane boundary.

BCSpecial 

Test-case-specific boundary (Riemann, DMR, isentropic vortex, RT).

Definition at line 35 of file EulerBC.hpp.

◆ EulerModel

Enumerates the available compressible flow solver model configurations.

Each enumerant selects a combination of:

  • Physics dimension (dim): number of velocity components in the equations. All models use dim=3 (three velocity components) except NS_2D which uses dim=2.
  • Geometry dimension (gDim): mesh spatial dimension (2D or 3D). Models without _3D suffix use gDim=2; those with _3D use gDim=3.
  • Turbulence model: plain NS (no turbulence transport), Spalart-Allmaras (SA, +1 var), two-equation RANS (2EQ, +2 vars), or extended/dynamic (EX, runtime-determined).

The compile-time number of conservative variables (nVarsFixed) follows from the model:

Model nVarsFixed Variables
NS, NS_3D 5 \(\rho,\; \rho u,\; \rho v,\; \rho w,\; E\)
NS_2D 4 \(\rho,\; \rho u,\; \rho v,\; E\)
NS_SA, NS_SA_3D 6 NS + \(\rho\tilde{\nu}\)
NS_2EQ, NS_2EQ_3D 7 NS + \(\rho k,\; \rho\omega\)
NS_EX, NS_EX_3D Dynamic Runtime-determined (Eigen::Dynamic)
See also
getnVarsFixed() for the nVarsFixed mapping.
getDim_Fixed() for the physics dimension mapping.
getGeomDim_Fixed() for the geometry dimension mapping.
EulerModelTraits for a compile-time trait bundle.
Enumerator
NS 

Navier-Stokes, 2D geometry, 3D physics (5 vars).

NS_SA 

NS + Spalart-Allmaras, 2D geometry (6 vars).

NS_2D 

NS with 2D physics (2 velocity components, 4 vars). The only model with dim=2.

NS_3D 

Navier-Stokes, 3D geometry, 3D physics (5 vars).

NS_SA_3D 

NS + Spalart-Allmaras, 3D geometry (6 vars).

NS_2EQ 

NS + 2-equation RANS, 2D geometry (7 vars).

NS_2EQ_3D 

NS + 2-equation RANS, 3D geometry (7 vars).

NS_EX 

Extended NS, 2D geometry, dynamic nVars (Eigen::Dynamic).

NS_EX_3D 

Extended NS, 3D geometry, dynamic nVars (Eigen::Dynamic).

Definition at line 874 of file Euler.hpp.

◆ RANSModel

Enumerates the available RANS turbulence closure models.

Selects which turbulence transport equations are solved at runtime. This is orthogonal to the compile-time EulerModel selection: EulerModel determines the number of transport variables, while RANSModel selects the specific closure (e.g. which two-equation model to use when EulerModel is NS_2EQ or NS_2EQ_3D).

JSON serialization is provided via DNDS_DEFINE_ENUM_JSON.

Enumerator
RANS_Unknown 

Sentinel / uninitialized value.

RANS_None 

No turbulence model (laminar or DNS).

RANS_SA 

Spalart-Allmaras one-equation model.

RANS_KOWilcox 

Wilcox k-omega two-equation model.

RANS_KOSST 

Menter k-omega SST two-equation model.

RANS_RKE 

Realizable k-epsilon two-equation model.

Definition at line 897 of file Euler.hpp.

Function Documentation

◆ bcSettingsSchema()

nlohmann::ordered_json DNDS::Euler::bcSettingsSchema ( )
inline

Build a JSON Schema (draft-07) for the bcSettings array.

Each BC entry is a discriminated union keyed on "type". The schema uses oneOf with "const" on the discriminator so that IDEs can provide type-specific autocomplete.

There are three groups:

  • Flow BCs (BCFar, BCOut, BCOutP, BCIn, BCInPsTs): required: type, name, value optional: frameOption, anchorOption, integrationOption, valueExtra
  • Wall BCs (BCWall, BCWallInvis, BCWallIsothermal, BCSpecial): required: type, name optional: value, frameOption, integrationOption, specialOption, valueExtra (value is required for BCWallIsothermal)
  • Symmetry (BCSym): required: type, name optional: rectifyOption, integrationOption, valueExtra
Returns
An ordered_json object representing the JSON Schema array definition suitable for validating a user-supplied boundary condition configuration.

Definition at line 96 of file EulerBC.hpp.

Here is the caller graph for this function:

◆ DNDS_DEFINE_ENUM_JSON() [1/2]

DNDS::Euler::DNDS_DEFINE_ENUM_JSON ( EulerBCType  ,
{ {BCUnknown, nullptr}, {BCFar, "BCFar"}, {BCWall, "BCWall"}, {BCWallIsothermal, "BCWallIsothermal"}, {BCWallInvis, "BCWallInvis"}, {BCOut, "BCOut"}, {BCOutP, "BCOutP"}, {BCIn, "BCIn"}, {BCInPsTs, "BCInPsTs"}, {BCSym, "BCSym"}, {BCSpecial, "BCSpecial"}, }   
)

Convert an EulerBCType enumerator to its JSON string representation.

Parameters
typeThe boundary condition type to convert.
Returns
A human-readable string such as "BCFar", "BCWall", etc.

Definition at line 50 of file EulerBC.hpp.

◆ DNDS_DEFINE_ENUM_JSON() [2/2]

DNDS::Euler::DNDS_DEFINE_ENUM_JSON ( RANSModel  ,
{ {RANS_Unknown, nullptr}, {RANS_None, "RANS_None"}, {RANS_SA, "RANS_SA"}, {RANS_KOWilcox, "RANS_KOWilcox"}, {RANS_KOSST, "RANS_KOSST"}, {RANS_RKE, "RANS_RKE"}, }   
) const

Returns the compile-time (fixed) number of conservative variables for a given model.

This is a constexpr function suitable for use as a template argument. The mapping is:

  • NS, NS_3D: 5 ( \(\rho, \rho u, \rho v, \rho w, E\))
  • NS_SA, NS_SA_3D: 6 (+ \(\rho\tilde{\nu}\))
  • NS_2D: 4 ( \(\rho, \rho u, \rho v, E\))
  • NS_2EQ, NS_2EQ_3D: 7 (+ \(\rho k, \rho\omega\))
  • NS_EX, NS_EX_3D: Eigen::Dynamic (determined at runtime)
Parameters
modelThe Euler model enumerant.
Returns
Compile-time variable count, or Eigen::Dynamic (-1) for extended models.
See also
getNVars() for a version that always returns a positive runtime value.

Definition at line 907 of file Euler.hpp.

◆ DNDS_SWITCH_INTELLISENSE() [1/2]

DNDS::Euler::DNDS_SWITCH_INTELLISENSE ( template< EulerModel model >  ,
template<>   
)

[[deprecated]] Forward LU-SGS sweep. Use UpdateSGS with uIncIsZero=true instead.

Initialize the running environment for the time-marching loop.

Apply the preconditioner for the GMRES linear solver.

Dispatch the linear solve for implicit time stepping.

Main implicit time-marching loop for the compressible Navier-Stokes solver.

Evaluate the spatial right-hand side (RHS) of the semi-discrete equations.

Generate ghost state for total-condition (stagnation) inflow BC.

Generate ghost state for prescribed inflow boundary condition.

Generate ghost state for supersonic/extrapolation outflow BC.

Generate ghost state for viscous (no-slip) wall boundary conditions.

Generate ghost state for inviscid wall (slip wall) and symmetry boundary conditions.

Generate ghost state for special farfield BCs (DMR, Rayleigh-Taylor, etc.).

Generate ghost state for farfield and back-pressure outlet boundary conditions.

Dispatcher for boundary ghost-value generation, routing to per-BC-type handlers.

Compute source terms for a cell at a quadrature point.

Evaluate the numerical flux at a single face for a batch of quadrature points.

evaluates dt and facial spectral radius

Dispatcher for wall distance computation, selecting method based on settings.

Compute wall distance by solving a Poisson equation with pseudo-time stepping.

Compute geometric wall distance using batched CGAL AABB tree queries.

Compute geometric wall distance using CGAL AABB tree (all-to-all broadcast).

Compute face-averaged wall distances from cell wall distances.

Collect wall boundary face triangulations for AABB-based wall distance computation.

Evaluate the positivity-preserving reconstruction limiter coefficient (beta) per cell.

Apply positivity-preserving gradient limiter to reconstructed solution gradients.

Solve the local LU-factorized Jacobian system for the increment.

[[deprecated]] Backward LU-SGS sweep. Use UpdateSGS instead.

Performs a forward Gauss-Seidel sweep using diagonal inversion and lower-triangular off-diagonal flux Jacobian contributions. Retained for callers that alias uInc==uIncNew.

Parameters
alphaDiagDiagonal scaling factor.
tCurrent simulation time.
rhsRight-hand side residual.
uConservative variable DOF array.
uIncCurrent increment (input, may alias uIncNew).
JDiagDiagonal Jacobian block.
uIncNewUpdated increment (output, may alias uInc).

Performs a backward Gauss-Seidel sweep using diagonal inversion and upper-triangular off-diagonal flux Jacobian contributions. Adds the correction to uIncNew in-place.

Parameters
alphaDiagDiagonal scaling factor.
tCurrent simulation time.
rhsRight-hand side residual.
uConservative variable DOF array.
uIncCurrent increment (input).
JDiagDiagonal Jacobian block.
uIncNewUpdated increment (input/output, correction added).

Applies the pre-factorized local LU solve (from LUSGSMatrixToJacobianLU) to compute the new increment from the RHS residual and off-diagonal contributions. Supports accumulation onto a nonzero initial increment when uIncIsZero is false.

Parameters
alphaDiagDiagonal scaling factor.
tCurrent simulation time.
rhsRight-hand side residual.
uConservative variable DOF array.
uIncCurrent increment (input).
uIncNewUpdated increment (output, must not alias uInc).
bBufTemporary buffer for RHS assembly.
JDiagDiagonal Jacobian block.
jacLUPre-factorized local LU structure.
uIncIsZeroIf true, uInc is zero and can be skipped in accumulation.
sumIncComponent-wise L1 norm of the new increment (output).

Limits the gradient so that density and internal energy remain above threshold values at all face quadrature points (Zhang-Shu style for gradients). Optionally disables shock-based limiting via flags.

Parameters
uConservative variable DOF array.
uGradInput gradient array.
uGradNewLimited gradient array (output).
flagsBitfield flags (e.g., LIMITER_UGRAD_Disable_Shock_Limiter).

Computes a scalar beta in [0,1] for each cell such that the reconstructed solution at all face and (optionally) volume quadrature points maintains positive density and pressure. Uses bisection to find the maximum allowable beta. Reports the total number of limited cells and the global minimum beta.

Parameters
uConservative variable DOF array.
uRecReconstruction coefficients (may be modified for order reduction).
uRecBetaPer-cell limiter coefficient (output).
nLimTotal number of limited cells across all MPI ranks (output).
betaMinGlobal minimum beta value (output).
flagEvaluation mode flags.

Iterates over wall and isothermal-wall boundary faces, decomposes them into triangles (handling Line, Tri, and Quad element types), and appends to the output vector. When useQuadPatches is true, uses quadrature-point-based sub-patches instead of vertex-based triangulation.

Parameters
useQuadPatchesIf true, create triangle patches from quadrature points.
trianglesOutCollected triangles as 3x3 matrices (columns = vertices) (output).

For each face, averages the wall distance of its two adjacent cells (or uses the single owner cell for boundary faces) and stores in dWallFace.

Gathers all wall triangles globally, builds a CGAL AABB tree, and queries the nearest distance for each cell's quadrature points. Uses execution serialization (wallDistExection) to limit concurrent memory usage. Clamps result to minWallDist.

Builds a local AABB tree from wall triangles (using quadrature patches), then performs batched distributed queries where cells are processed in configurable load-sized chunks. Supports refinement mode (wallDistScheme==20) for selective re-querying of cells below wallDistRefineMax.

Solves a p-Poisson equation (grad^2 phi = -1) with Dirichlet BC (phi=0) on wall boundaries and Neumann BC elsewhere, using Jacobi or GMRES iterations. The wall distance is recovered from the solution as d = |grad(phi)| + sqrt(|grad(phi)|^2 + 2*phi). Configurable via settings: wallDistNJacobiSweep, wallDistIter, wallDistPoissonP, etc.

Routes to GetWallDist_AABB (schemes 0,1,20), GetWallDist_BatchedAABB (schemes 2,20), or GetWallDist_Poisson (scheme 3), then computes face-averaged distances.

Estimate the local time step and face spectral radii for CFL-based time stepping.

For each face, evaluates L/R states using 2nd-order reconstruction, computes the convective and viscous eigenvalue estimates, and stores face spectral radii (lambdaFace, lambdaFace0..4, lambdaFaceC, lambdaFaceVis). Then for each cell, sums the face spectral radii to obtain a CFL-limited local time step, optionally clamped by MaxDt or made uniform (global minimum).

Parameters
dtLocal time step per cell (output).
uConservative variable DOF array.
uRecReconstruction coefficients.
CFLCFL number for time step scaling.
dtMinallGlobal minimum time step (output, MPI-reduced).
MaxDtMaximum allowed time step (upper clamp).
UseLocaldtIf true, use local time stepping; otherwise set all cells to dtMinall.
tCurrent simulation time.
flagsBitfield flags (e.g., DT_Dont_update_lambda01234).

Computes the inviscid and viscous numerical flux at a face given batched L/R states, gradients, and face geometry. Dispatches to the configured Riemann solver type. Outputs left/right viscous corrections (FLfix, FRfix), the total flux increment (finc), and per-quadrature-point eigenvalue estimates (lam0V, lam123V, lam4V).

Parameters
ULxyBatched left states at quadrature points.
URxyBatched right states at quadrature points.
ULMeanXyLeft cell mean state.
URMeanXyRight cell mean state.
DiffUxyBatched conservative gradient at face.
DiffUxyPrimBatched primitive gradient at face.
unitNormBatched outward unit normals at quadrature points.
vgXYBatched grid velocity at quadrature points.
unitNormCFace-center unit normal.
vgCFace-center grid velocity.
FLfixLeft viscous correction (output).
FRfixRight viscous correction (output).
fincTotal flux increment (output).
lam0VEigenvalue estimate for wave 0 (output).
lam123VEigenvalue estimate for waves 1-3 (output).
lam4VEigenvalue estimate for wave 4 (output).
btypeBoundary type index of the face.
rsTypeRiemann solver type to use.
iFaceFace index.
ignoreVisIf true, skip viscous flux computation.

Evaluates volume source terms including: constant mass force, rotating-frame fictitious forces (Coriolis and centrifugal), axisymmetric source terms, SA turbulence model source, and k-omega/k-epsilon RANS model sources.

Parameters
UMeanXyCell mean conservative state.
DiffUxyGradient of conservative variables at the quadrature point.
pPhyPhysical coordinates of the quadrature point.
jacobianSource-term Jacobian matrix (output when Mode==2).
iCellCell index.
igQuadrature point index within the cell (-1 for cell center).
Mode0: return source vector; 1: return diagonal Jacobian; 2: return full Jacobian.
Returns
Source term vector.

Selects the appropriate boundary value generator based on the BC type (farfield, wall, inviscid wall, outflow, inflow, total-condition inflow, symmetry, etc.) and returns the ghost (right) state for the boundary face.

Parameters
ULxyLeft (interior) state; may be modified by some handlers (e.g., fixUL).
ULMeanXyLeft cell mean state.
iCellOwner cell index.
iFaceFace index.
iGQuadrature point index (-1 for face center).
uNormOutward unit normal vector.
normBaseLocal coordinate frame built from the normal.
pPhysicsPhysical coordinates of the quadrature point.
tCurrent simulation time.
btypeBoundary zone identifier.
fixULIf true, handler may modify ULxy (e.g., enforce zero SA at wall).
geomModeGeometry mode for normal evaluation.
linModeIf 1, return linearized boundary perturbation (for Jacobian).
Returns
Ghost (right) state for the boundary face.

Uses Riemann-invariant-based switching: supersonic outflow extrapolates the interior state; subsonic outflow imposes farfield pressure; subsonic inflow imposes farfield values but retains interior pressure; supersonic inflow uses the full farfield state. Supports anchor-point and radial-profile pressure corrections, CL-driver AoA rotation, and rotating-frame transformations.

Parameters
ULxyLeft (interior) state.
ULMeanXyLeft cell mean state.
iCellOwner cell index.
iFaceFace index.
iGQuadrature point index.
uNormOutward unit normal.
normBaseLocal coordinate frame.
pPhysicsPhysical coordinates.
tCurrent simulation time.
btypeBoundary zone identifier.
fixULIf true, may modify ULxy.
geomModeGeometry mode.
Returns
Ghost (right) state for the farfield face.

Handles built-in special boundary conditions identified by reserved BC IDs: double Mach reflection (DMR), Rayleigh-Taylor instability, and user-specified special options. Uses Riemann-invariant switching like the standard farfield.

Parameters
ULxyLeft (interior) state.
ULMeanXyLeft cell mean state.
iCellOwner cell index.
iFaceFace index.
iGQuadrature point index.
uNormOutward unit normal.
normBaseLocal coordinate frame.
pPhysicsPhysical coordinates.
tCurrent simulation time.
btypeBoundary zone identifier (special reserved ID).
Returns
Ghost (right) state for the special farfield face.

Mirrors the velocity component normal to the wall while preserving the tangential component, density, and energy. Handles rotating-frame transformations.

Parameters
ULxyLeft (interior) state.
ULMeanXyLeft cell mean state.
iCellOwner cell index.
iFaceFace index.
iGQuadrature point index.
uNormOutward unit normal.
normBaseLocal coordinate frame.
pPhysicsPhysical coordinates.
tCurrent simulation time.
btypeBoundary zone identifier.
Returns
Ghost (right) state with reflected normal velocity.

Reverses the velocity to enforce zero velocity at the wall (ghost mirroring). For isothermal walls, adjusts density to match the prescribed wall temperature. Handles SA turbulence variable (negated or zeroed at wall) and k-omega/k-epsilon wall boundary conditions including the realizable k-epsilon wall epsilon formula.

Parameters
ULxyLeft (interior) state (may be modified when fixUL is true for SA).
ULMeanXyLeft cell mean state.
iCellOwner cell index.
iFaceFace index.
iGQuadrature point index.
uNormOutward unit normal.
normBaseLocal coordinate frame.
pPhysicsPhysical coordinates.
tCurrent simulation time.
btypeBoundary zone identifier.
fixULIf true, may zero SA variables in ULxy at the wall.
linModeIf 1, return linearized ghost perturbation for Jacobian.
Returns
Ghost (right) state with no-slip wall condition.

Simply returns the interior state as the ghost value, allowing all waves to exit the domain without reflection.

Parameters
ULxyLeft (interior) state.
ULMeanXyLeft cell mean state.
iCellOwner cell index.
iFaceFace index.
iGQuadrature point index.
uNormOutward unit normal.
normBaseLocal coordinate frame.
pPhysicsPhysical coordinates.
tCurrent simulation time.
btypeBoundary zone identifier.
Returns
Ghost state equal to the interior state (full extrapolation).

Sets the ghost state to the BC-prescribed conservative values. Applies CL-driver AoA rotation and rotating-frame transformation if applicable.

Parameters
ULxyLeft (interior) state (unused for pure inflow).
ULMeanXyLeft cell mean state.
iCellOwner cell index.
iFaceFace index.
iGQuadrature point index.
uNormOutward unit normal.
normBaseLocal coordinate frame.
pPhysicsPhysical coordinates.
tCurrent simulation time.
btypeBoundary zone identifier.
Returns
Ghost state set to prescribed inflow values.

Given prescribed stagnation pressure and temperature, computes the static conditions from the interior velocity magnitude using isentropic relations. The flow direction is taken from the BC-prescribed direction vector. Applies CL-driver AoA rotation and rotating-frame transformation.

Parameters
ULxyLeft (interior) state (used for velocity magnitude).
ULMeanXyLeft cell mean state.
iCellOwner cell index.
iFaceFace index.
iGQuadrature point index.
uNormOutward unit normal.
normBaseLocal coordinate frame.
pPhysicsPhysical coordinates.
tCurrent simulation time.
btypeBoundary zone identifier.
Returns
Ghost state computed from total conditions and interior velocity.

This is the core spatial operator. It performs the following steps:

  1. Loop over internal faces: reconstruct L/R states, compute inviscid numerical flux via the selected Riemann solver, and accumulate face contributions to cell RHS.
  2. Loop over boundary faces: generate ghost boundary values, compute boundary flux.
  3. Loop over cells: compute viscous flux, source terms (RANS, mass force, rotating frame), and add volume contributions.
  4. Track faces where reduced-order reconstruction is used for robustness.
  5. Optionally use direct 2nd-order reconstruction methods for multigrid.
Parameters
rhsCell RHS residual (output, zeroed then accumulated).
JSourceSource-term Jacobian diagonal block (output for implicit).
uConservative variable DOF array.
uRecUnlimUnlimited reconstruction coefficients.
uRecLimited reconstruction coefficients.
uRecBetaPer-cell reconstruction limiter coefficient.
cellRHSAlphaPer-cell RHS scaling factor for positivity preservation.
onlyOnHalfAlphaIf true, evaluate RHS only on cells with alpha < 1.
tCurrent simulation time.
flagsBitfield flags controlling viscosity, integration recording, direct 2nd-order reconstruction, and other options.

Performs the following at each time step:

  1. CFL ramping with exponential growth and linear fallback on divergence.
  2. Reconstruction update: VR coefficient computation, gradient limiting, positivity-preserving beta evaluation.
  3. Local time step (dTau) estimation from face eigenvalues.
  4. RHS evaluation via EvaluateRHS.
  5. Implicit time stepping using either LU-SGS sweeps or GMRES with LU preconditioner.
  6. Positivity-preserving alpha limiter for the update.
  7. Component-wise L1 residual monitoring with convergence checks.
  8. CL-driver integration for angle-of-attack adjustment.
  9. Checkpoint/restart writing on schedule or SIGUSR1 signal.
  10. Time averaging accumulation.
  11. VTK output scheduling based on time or iteration.

Solves the linearized system using either SGS sweeps (gmresCode=0) or FGMRES with preconditioner (gmresCode=1/2). Supports multi-grid levels with per-level configuration of iteration counts and solver parameters.

Parameters
alphaDiagDiagonal scaling factor for the implicit operator.
tCurrent simulation time.
cresRHS residual (right-hand side of the linear system).
cxCurrent solution DOF.
cxIncSolution increment (input/output).
uRecCReconstruction coefficients.
uRecIncCReconstruction increment (for SGS-with-rec mode).
JDCDiagonal Jacobian block.
gmresGMRES solver object.
gridLevelMultigrid level (0 = finest).

Depending on jacobiCode: applies symmetric SGS sweeps (code 0/1) or local LU factorization solve (code 2). Manages the inputIsZero and hasLUDone state flags across GMRES restarts.

Parameters
alphaDiagDiagonal scaling factor.
tCurrent simulation time.
cresRight-hand side for preconditioning.
cxCurrent solution DOF.
cxIncPreconditioned result (output).
uTempTemporary DOF buffer.
JDCDiagonal Jacobian block.
sgsResSGS residual for convergence tracking.
inputIsZeroTracks whether cxInc is zero (updated in-place).
hasLUDoneTracks whether LU factorization has been computed (updated).
gridLevelMultigrid level.

Sets up the ODE solver, error logging stream, temporary arrays, signal handlers, output data structures, and applies steady-mode overrides when configured. Must be called before entering the main time-stepping loop.

Parameters
runningEnvironmentRunning environment structure to populate (output).

always inner here

early exit, reconstruction is good it self

todo: periodic!!

viscous flux

high fix

using velo instead of velo + a

using velo instead of velo + a

k's normal stress

using velo instead of velo + a

DES mesh lengths should be cached!

missing hWallNormal!

TODO: make really block jacobian

TODO: make really block jacobian

disabled now, as ULxyUnlim may have negative pressure failing in generateBV

USING REAL GEOMETRICAL

don't forget this

using consistent rho U scale

dim balancing here

consecutive pcg is bad in 0012, using separate pcg

uRecIncC now stores uRecIncrement

uGradBufNoLim already existent in fdtau

note: in HM3 IV test for TPMG, this O1 version makes convergence slower compared to the O2 one, why?

this also need to update!

here we borrow PMG's level1 setting into TPMG

uGradBufNoLim already existent in fdtau

uGradBufNoLim already existent in fdtau

warning! dTauC is overwritten

overwrites cxInc, cxTemp and resTemp do not overwrite cres! (as we might use cres for evaluation of convergence)

need dim balancing here

using consistent rho U scale

start as zero

always inner here

early exit, reconstruction is good it self

todo: periodic!!

viscous flux

high fix

using velo instead of velo + a

using velo instead of velo + a

k's normal stress

using velo instead of velo + a

DES mesh lengths should be cached!

missing hWallNormal!

TODO: make really block jacobian

TODO: make really block jacobian

disabled now, as ULxyUnlim may have negative pressure failing in generateBV

USING REAL GEOMETRICAL

don't forget this

using consistent rho U scale

dim balancing here

consecutive pcg is bad in 0012, using separate pcg

uRecIncC now stores uRecIncrement

uGradBufNoLim already existent in fdtau

note: in HM3 IV test for TPMG, this O1 version makes convergence slower compared to the O2 one, why?

this also need to update!

here we borrow PMG's level1 setting into TPMG

uGradBufNoLim already existent in fdtau

uGradBufNoLim already existent in fdtau

warning! dTauC is overwritten

overwrites cxInc, cxTemp and resTemp do not overwrite cres! (as we might use cres for evaluation of convergence)

need dim balancing here

using consistent rho U scale

start as zero

Definition at line 361 of file EulerEvaluator.hxx.

Here is the call graph for this function:

◆ DNDS_SWITCH_INTELLISENSE() [2/2]

DNDS::Euler::DNDS_SWITCH_INTELLISENSE ( template< EulerModel model >  )

Assemble diagonal Jacobian block for each cell for LU-SGS preconditioning.

Initialize conservative variable DOF from farfield and special-field initializers.

Assemble and factorize the local LU Jacobian from diagonal and off-diagonal blocks.

Compute the implicit matrix-vector product A*uInc for GMRES.

Computes the spectral radius contribution from face eigenvalues and accumulates into JDiag. When settings.useRoeJacobian is enabled, also computes Roe flux Jacobian block contributions including boundary face linearization.

Parameters
JDiagDiagonal Jacobian block storage (output, cleared then filled).
JSourceSource-term Jacobian diagonal block (must match JDiag block mode).
dTauLocal pseudo-time step per cell.
dtGlobal physical time step.
alphaDiagDiagonal scaling factor for the spectral radius.
uConservative variable DOF array.
uRecReconstruction coefficients.
jacobianCodeJacobian assembly mode (must be 0).
tCurrent simulation time.

Evaluates the action of the implicit operator (diagonal block + off-diagonal flux Jacobian contributions) on the increment vector, storing the result in AuInc.

Parameters
alphaDiagDiagonal scaling factor.
tCurrent simulation time.
uConservative variable DOF array.
uIncIncrement vector to multiply.
JDiagDiagonal Jacobian block.
AuIncResult of the matrix-vector product (output).

Populates the sparse local LU structure with diagonal blocks from JDiag and off-diagonal Roe flux Jacobian entries for local-partition neighbor cells, then performs in-place LU decomposition.

Parameters
alphaDiagDiagonal scaling factor for face flux Jacobian.
tCurrent simulation time.
uConservative variable DOF array.
JDiagDiagonal Jacobian block.
jacLULocal LU factorization structure (output).

Sets all cells to the farfield static value, then applies model-specific initialization (SA wall-distance scaling), box/plane/exprtk region overrides, and built-in special fields (isentropic vortex, Rayleigh-Taylor, Taylor-Green, double Mach reflection, etc.) based on settings.

Parameters
uConservative variable DOF array (output, fully initialized).

always inner here

always inner here

always inner here

always inner here

Definition at line 19 of file EulerEvaluator.hxx.

◆ paste_read_restart_with_cell_ordering()

void DNDS::Euler::paste_read_restart_with_cell_ordering ( EulerSolver< model > &  self,
typename EulerSolver< model >::TDof &  u,
typename EulerSolver< model >::TDof &  uRead,
Serializer::SerializerBaseSSP  serializerP 
)

Reorder restart data to match current cell ordering using cell2cellOrig mapping.

When restart data was saved with a different partition layout (JSON path), reads the original cell ordering and permutes the read DOF to match the current mesh partition. Falls back to direct copy with a warning if cell2cellOrig is not available.

Parameters
selfReference to the EulerSolver instance.
uTarget DOF array (output, reordered).
uReadDOF array as read from file (input).
serializerPSerializer used to read cell ordering metadata.

Definition at line 943 of file EulerSolver_PrintData.hxx.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ RunSingleBlockConsoleApp()

template<EulerModel model>
int DNDS::Euler::RunSingleBlockConsoleApp ( int  argc,
char *  argv[] 
)

Main entry point for single-block solver console applications.

This function template is instantiated for each EulerModel variant and serves as the complete main() implementation for that solver executable.

Execution flow:

  1. Initialize MPI and build default/user configuration file paths.
  2. Parse command-line arguments via argparse:
    • Positional config — path to the user JSON configuration file.
    • Positional field_n_variables (extended models with DynamicSize only).
    • -k/--overwrite_key and -v/--overwrite_value — repeated pairs for ad-hoc JSON config overrides.
    • --debug — attach-debugger mode (MPI hold).
    • --emit-schema — print the full JSON Schema for this solver's configuration and exit.
  3. In --emit-schema mode: instantiate a default EulerSolver, emit the schema produced by DNDS_DECLARE_CONFIG, and return.
  4. In normal mode: construct an EulerSolver, load the default config then the user config (with key/value overrides), call ReadMeshAndInitialize(), and RunImplicitEuler().
Template Parameters
modelCompile-time EulerModel selecting the equation set / dimension.
Parameters
argcArgument count forwarded from main().
argvArgument vector forwarded from main().
Returns
0 on success (abnormal paths call std::abort()).

Definition at line 79 of file SingleBlockApp.hpp.

Here is the call graph for this function: