|
DNDSR 0.1.0.dev1+gcd065ad
Distributed Numeric Data Structure for CFV
|
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. | |
| 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.
| 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.
Boundary condition type identifiers for compressible flow solvers.
Each enumerator selects the physics applied at a boundary zone:
Definition at line 35 of file EulerBC.hpp.
Enumerates the available compressible flow solver model configurations.
Each enumerant selects a combination of:
dim): number of velocity components in the equations. All models use dim=3 (three velocity components) except NS_2D which uses dim=2.gDim): mesh spatial dimension (2D or 3D). Models without _3D suffix use gDim=2; those with _3D use gDim=3.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) |
| 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). |
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.
|
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:
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.
| 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.
| type | The boundary condition type to convert. |
Definition at line 50 of file EulerBC.hpp.
| 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:
Eigen::Dynamic (determined at runtime)| model | The Euler model enumerant. |
Eigen::Dynamic (-1) for extended models.| 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.
| alphaDiag | Diagonal scaling factor. |
| t | Current simulation time. |
| rhs | Right-hand side residual. |
| u | Conservative variable DOF array. |
| uInc | Current increment (input, may alias uIncNew). |
| JDiag | Diagonal Jacobian block. |
| uIncNew | Updated 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.
| alphaDiag | Diagonal scaling factor. |
| t | Current simulation time. |
| rhs | Right-hand side residual. |
| u | Conservative variable DOF array. |
| uInc | Current increment (input). |
| JDiag | Diagonal Jacobian block. |
| uIncNew | Updated 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.
| alphaDiag | Diagonal scaling factor. |
| t | Current simulation time. |
| rhs | Right-hand side residual. |
| u | Conservative variable DOF array. |
| uInc | Current increment (input). |
| uIncNew | Updated increment (output, must not alias uInc). |
| bBuf | Temporary buffer for RHS assembly. |
| JDiag | Diagonal Jacobian block. |
| jacLU | Pre-factorized local LU structure. |
| uIncIsZero | If true, uInc is zero and can be skipped in accumulation. |
| sumInc | Component-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.
| u | Conservative variable DOF array. |
| uGrad | Input gradient array. |
| uGradNew | Limited gradient array (output). |
| flags | Bitfield 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.
| u | Conservative variable DOF array. |
| uRec | Reconstruction coefficients (may be modified for order reduction). |
| uRecBeta | Per-cell limiter coefficient (output). |
| nLim | Total number of limited cells across all MPI ranks (output). |
| betaMin | Global minimum beta value (output). |
| flag | Evaluation 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.
| useQuadPatches | If true, create triangle patches from quadrature points. |
| trianglesOut | Collected 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).
| dt | Local time step per cell (output). |
| u | Conservative variable DOF array. |
| uRec | Reconstruction coefficients. |
| CFL | CFL number for time step scaling. |
| dtMinall | Global minimum time step (output, MPI-reduced). |
| MaxDt | Maximum allowed time step (upper clamp). |
| UseLocaldt | If true, use local time stepping; otherwise set all cells to dtMinall. |
| t | Current simulation time. |
| flags | Bitfield 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).
| ULxy | Batched left states at quadrature points. |
| URxy | Batched right states at quadrature points. |
| ULMeanXy | Left cell mean state. |
| URMeanXy | Right cell mean state. |
| DiffUxy | Batched conservative gradient at face. |
| DiffUxyPrim | Batched primitive gradient at face. |
| unitNorm | Batched outward unit normals at quadrature points. |
| vgXY | Batched grid velocity at quadrature points. |
| unitNormC | Face-center unit normal. |
| vgC | Face-center grid velocity. |
| FLfix | Left viscous correction (output). |
| FRfix | Right viscous correction (output). |
| finc | Total flux increment (output). |
| lam0V | Eigenvalue estimate for wave 0 (output). |
| lam123V | Eigenvalue estimate for waves 1-3 (output). |
| lam4V | Eigenvalue estimate for wave 4 (output). |
| btype | Boundary type index of the face. |
| rsType | Riemann solver type to use. |
| iFace | Face index. |
| ignoreVis | If 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.
| UMeanXy | Cell mean conservative state. |
| DiffUxy | Gradient of conservative variables at the quadrature point. |
| pPhy | Physical coordinates of the quadrature point. |
| jacobian | Source-term Jacobian matrix (output when Mode==2). |
| iCell | Cell index. |
| ig | Quadrature point index within the cell (-1 for cell center). |
| Mode | 0: return source vector; 1: return diagonal Jacobian; 2: return full Jacobian. |
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.
| ULxy | Left (interior) state; may be modified by some handlers (e.g., fixUL). |
| ULMeanXy | Left cell mean state. |
| iCell | Owner cell index. |
| iFace | Face index. |
| iG | Quadrature point index (-1 for face center). |
| uNorm | Outward unit normal vector. |
| normBase | Local coordinate frame built from the normal. |
| pPhysics | Physical coordinates of the quadrature point. |
| t | Current simulation time. |
| btype | Boundary zone identifier. |
| fixUL | If true, handler may modify ULxy (e.g., enforce zero SA at wall). |
| geomMode | Geometry mode for normal evaluation. |
| linMode | If 1, return linearized boundary perturbation (for Jacobian). |
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.
| ULxy | Left (interior) state. |
| ULMeanXy | Left cell mean state. |
| iCell | Owner cell index. |
| iFace | Face index. |
| iG | Quadrature point index. |
| uNorm | Outward unit normal. |
| normBase | Local coordinate frame. |
| pPhysics | Physical coordinates. |
| t | Current simulation time. |
| btype | Boundary zone identifier. |
| fixUL | If true, may modify ULxy. |
| geomMode | Geometry mode. |
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.
| ULxy | Left (interior) state. |
| ULMeanXy | Left cell mean state. |
| iCell | Owner cell index. |
| iFace | Face index. |
| iG | Quadrature point index. |
| uNorm | Outward unit normal. |
| normBase | Local coordinate frame. |
| pPhysics | Physical coordinates. |
| t | Current simulation time. |
| btype | Boundary zone identifier (special reserved ID). |
Mirrors the velocity component normal to the wall while preserving the tangential component, density, and energy. Handles rotating-frame transformations.
| ULxy | Left (interior) state. |
| ULMeanXy | Left cell mean state. |
| iCell | Owner cell index. |
| iFace | Face index. |
| iG | Quadrature point index. |
| uNorm | Outward unit normal. |
| normBase | Local coordinate frame. |
| pPhysics | Physical coordinates. |
| t | Current simulation time. |
| btype | Boundary zone identifier. |
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.
| ULxy | Left (interior) state (may be modified when fixUL is true for SA). |
| ULMeanXy | Left cell mean state. |
| iCell | Owner cell index. |
| iFace | Face index. |
| iG | Quadrature point index. |
| uNorm | Outward unit normal. |
| normBase | Local coordinate frame. |
| pPhysics | Physical coordinates. |
| t | Current simulation time. |
| btype | Boundary zone identifier. |
| fixUL | If true, may zero SA variables in ULxy at the wall. |
| linMode | If 1, return linearized ghost perturbation for Jacobian. |
Simply returns the interior state as the ghost value, allowing all waves to exit the domain without reflection.
| ULxy | Left (interior) state. |
| ULMeanXy | Left cell mean state. |
| iCell | Owner cell index. |
| iFace | Face index. |
| iG | Quadrature point index. |
| uNorm | Outward unit normal. |
| normBase | Local coordinate frame. |
| pPhysics | Physical coordinates. |
| t | Current simulation time. |
| btype | Boundary zone identifier. |
Sets the ghost state to the BC-prescribed conservative values. Applies CL-driver AoA rotation and rotating-frame transformation if applicable.
| ULxy | Left (interior) state (unused for pure inflow). |
| ULMeanXy | Left cell mean state. |
| iCell | Owner cell index. |
| iFace | Face index. |
| iG | Quadrature point index. |
| uNorm | Outward unit normal. |
| normBase | Local coordinate frame. |
| pPhysics | Physical coordinates. |
| t | Current simulation time. |
| btype | Boundary zone identifier. |
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.
| ULxy | Left (interior) state (used for velocity magnitude). |
| ULMeanXy | Left cell mean state. |
| iCell | Owner cell index. |
| iFace | Face index. |
| iG | Quadrature point index. |
| uNorm | Outward unit normal. |
| normBase | Local coordinate frame. |
| pPhysics | Physical coordinates. |
| t | Current simulation time. |
| btype | Boundary zone identifier. |
This is the core spatial operator. It performs the following steps:
| rhs | Cell RHS residual (output, zeroed then accumulated). |
| JSource | Source-term Jacobian diagonal block (output for implicit). |
| u | Conservative variable DOF array. |
| uRecUnlim | Unlimited reconstruction coefficients. |
| uRec | Limited reconstruction coefficients. |
| uRecBeta | Per-cell reconstruction limiter coefficient. |
| cellRHSAlpha | Per-cell RHS scaling factor for positivity preservation. |
| onlyOnHalfAlpha | If true, evaluate RHS only on cells with alpha < 1. |
| t | Current simulation time. |
| flags | Bitfield flags controlling viscosity, integration recording, direct 2nd-order reconstruction, and other options. |
Performs the following at each time step:
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.
| alphaDiag | Diagonal scaling factor for the implicit operator. |
| t | Current simulation time. |
| cres | RHS residual (right-hand side of the linear system). |
| cx | Current solution DOF. |
| cxInc | Solution increment (input/output). |
| uRecC | Reconstruction coefficients. |
| uRecIncC | Reconstruction increment (for SGS-with-rec mode). |
| JDC | Diagonal Jacobian block. |
| gmres | GMRES solver object. |
| gridLevel | Multigrid 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.
| alphaDiag | Diagonal scaling factor. |
| t | Current simulation time. |
| cres | Right-hand side for preconditioning. |
| cx | Current solution DOF. |
| cxInc | Preconditioned result (output). |
| uTemp | Temporary DOF buffer. |
| JDC | Diagonal Jacobian block. |
| sgsRes | SGS residual for convergence tracking. |
| inputIsZero | Tracks whether cxInc is zero (updated in-place). |
| hasLUDone | Tracks whether LU factorization has been computed (updated). |
| gridLevel | Multigrid 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.
| runningEnvironment | Running 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.
| 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.
| JDiag | Diagonal Jacobian block storage (output, cleared then filled). |
| JSource | Source-term Jacobian diagonal block (must match JDiag block mode). |
| dTau | Local pseudo-time step per cell. |
| dt | Global physical time step. |
| alphaDiag | Diagonal scaling factor for the spectral radius. |
| u | Conservative variable DOF array. |
| uRec | Reconstruction coefficients. |
| jacobianCode | Jacobian assembly mode (must be 0). |
| t | Current 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.
| alphaDiag | Diagonal scaling factor. |
| t | Current simulation time. |
| u | Conservative variable DOF array. |
| uInc | Increment vector to multiply. |
| JDiag | Diagonal Jacobian block. |
| AuInc | Result 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.
| alphaDiag | Diagonal scaling factor for face flux Jacobian. |
| t | Current simulation time. |
| u | Conservative variable DOF array. |
| JDiag | Diagonal Jacobian block. |
| jacLU | Local 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.
| u | Conservative 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.
| 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.
| self | Reference to the EulerSolver instance. |
| u | Target DOF array (output, reordered). |
| uRead | DOF array as read from file (input). |
| serializerP | Serializer used to read cell ordering metadata. |
Definition at line 943 of file EulerSolver_PrintData.hxx.
| 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:
config — path to the user JSON configuration file.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.--emit-schema mode: instantiate a default EulerSolver, emit the schema produced by DNDS_DECLARE_CONFIG, and return.| model | Compile-time EulerModel selecting the equation set / dimension. |
Definition at line 79 of file SingleBlockApp.hpp.