|
DNDSR 0.1.0.dev1+gcd065ad
Distributed Numeric Data Structure for CFV
|
Core finite-volume evaluator for compressible Navier-Stokes / Euler equations. More...
#include <EulerEvaluator.hpp>
Classes | |
| struct | OutputOverlapDataRefs |
| References to arrays needed by the output data picker. More... | |
Public Types | |
| using | Traits = EulerModelTraits< model > |
| Compile-time traits: hasSA, has2EQ, etc. | |
| typedef Eigen::VectorFMTSafe< real, dim > | TVec |
| Spatial vector (dim components). | |
| typedef Eigen::MatrixFMTSafe< real, dim, Eigen::Dynamic, Eigen::ColMajor, dim, MaxBatch > | TVec_Batch |
| Batch of spatial vectors. | |
| typedef Eigen::MatrixFMTSafe< real, dim, dim > | TMat |
| Spatial matrix (dim x dim). | |
| typedef Eigen::MatrixFMTSafe< real, dim, Eigen::Dynamic, Eigen::ColMajor, dim, MaxBatchMult(3)> | TMat_Batch |
| Batch of spatial matrices. | |
| typedef Eigen::VectorFMTSafe< real, nVarsFixed > | TU |
| Conservative variable vector (nVarsFixed components). | |
| typedef Eigen::MatrixFMTSafe< real, nVarsFixed, Eigen::Dynamic, Eigen::ColMajor, nVarsFixed, MaxBatch > | TU_Batch |
| Batch of conservative variable vectors. | |
| typedef Eigen::MatrixFMTSafe< real, 1, Eigen::Dynamic, Eigen::RowMajor, 1, MaxBatch > | TReal_Batch |
| Batch of scalar values. | |
| typedef Eigen::MatrixFMTSafe< real, nVarsFixed, nVarsFixed > | TJacobianU |
| Jacobian matrix (nVars x nVars) of the flux w.r.t. conserved variables. | |
| typedef Eigen::MatrixFMTSafe< real, dim, nVarsFixed > | TDiffU |
| Gradient of conserved variables (dim x nVars). | |
| typedef Eigen::MatrixFMTSafe< real, Eigen::Dynamic, nVarsFixed, Eigen::ColMajor, MaxBatchMult(3)> | TDiffU_Batch |
| Batch of gradient matrices. | |
| typedef Eigen::MatrixFMTSafe< real, nVarsFixed, dim > | TDiffUTransposed |
| Transposed gradient (nVars x dim). | |
| typedef ArrayDOFV< nVarsFixed > | TDof |
| Cell-centered DOF array (mean values). | |
| typedef ArrayRECV< nVarsFixed > | TRec |
| Reconstruction coefficient array (per-cell polynomial coefficients). | |
| typedef ArrayRECV< 1 > | TScalar |
| Scalar reconstruction coefficient array. | |
| typedef CFV::VariationalReconstruction< gDim > | TVFV |
| Variational reconstruction type for this geometric dimension. | |
| typedef ssp< CFV::VariationalReconstruction< gDim > > | TpVFV |
| Shared pointer to the variational reconstruction object. | |
| typedef ssp< BoundaryHandler< model > > | TpBCHandler |
| Shared pointer to the boundary condition handler. | |
| using | tFCompareField = std::function< TU(const Geom::tPoint &, real)> |
| Callback type for analytical comparison field. | |
| using | tFCompareFieldWeight = std::function< real(const Geom::tPoint &, real)> |
| Callback type for comparison weighting function. | |
Public Member Functions | |
| void | setPassiveDiscardSource (bool n) |
| Enable or disable discarding of passive-scalar source terms. | |
| int | GetAxisSymmetric () |
| Return the axisymmetric mode flag. | |
| EulerEvaluator (const decltype(mesh) &Nmesh, const decltype(vfv) &Nvfv, const decltype(pBCHandler) &npBCHandler, const EulerEvaluatorSettings< model > &nSettings, int n_nVars=getNVars(model)) | |
| Construct an EulerEvaluator and initialize all internal buffers. | |
| void | GetWallDist () |
| Compute wall distance for all cells (dispatcher). | |
| void | EvaluateDt (ArrayDOFV< 1 > &dt, ArrayDOFV< nVarsFixed > &u, ArrayRECV< nVarsFixed > &uRec, real CFL, real &dtMinall, real MaxDt, bool UseLocaldt, real t, uint64_t flags=DT_No_Flags) |
| Estimate the local or global time step for each cell. | |
| void | EvaluateRHS (ArrayDOFV< nVarsFixed > &rhs, JacobianDiagBlock< nVarsFixed > &JSource, ArrayDOFV< nVarsFixed > &u, ArrayRECV< nVarsFixed > &uRecUnlim, ArrayRECV< nVarsFixed > &uRec, ArrayDOFV< 1 > &uRecBeta, ArrayDOFV< 1 > &cellRHSAlpha, bool onlyOnHalfAlpha, real t, uint64_t flags=RHS_No_Flags) |
| Evaluate the spatial right-hand side of the semi-discrete system. | |
| void | LUSGSMatrixInit (JacobianDiagBlock< nVarsFixed > &JDiag, JacobianDiagBlock< nVarsFixed > &JSource, ArrayDOFV< 1 > &dTau, real dt, real alphaDiag, ArrayDOFV< nVarsFixed > &u, ArrayRECV< nVarsFixed > &uRec, int jacobianCode, real t) |
| Assemble the diagonal blocks of the implicit Jacobian for LU-SGS / SGS. | |
| void | LUSGSMatrixVec (real alphaDiag, real t, ArrayDOFV< nVarsFixed > &u, ArrayDOFV< nVarsFixed > &uInc, JacobianDiagBlock< nVarsFixed > &JDiag, ArrayDOFV< nVarsFixed > &AuInc) |
| Compute the matrix-vector product A * uInc for the implicit system. | |
| void | LUSGSMatrixToJacobianLU (real alphaDiag, real t, ArrayDOFV< nVarsFixed > &u, JacobianDiagBlock< nVarsFixed > &JDiag, JacobianLocalLU< nVarsFixed > &jacLU) |
| Build the local LU factorization of the Jacobian for direct solve. | |
| void | UpdateLUSGSForward (real alphaDiag, real t, ArrayDOFV< nVarsFixed > &rhs, ArrayDOFV< nVarsFixed > &u, ArrayDOFV< nVarsFixed > &uInc, JacobianDiagBlock< nVarsFixed > &JDiag, ArrayDOFV< nVarsFixed > &uIncNew) |
| Deprecated: use UpdateSGS with uIncIsZero=true instead. Kept for backward compatibility. Delegates to UpdateSGS. | |
| void | UpdateLUSGSBackward (real alphaDiag, real t, ArrayDOFV< nVarsFixed > &rhs, ArrayDOFV< nVarsFixed > &u, ArrayDOFV< nVarsFixed > &uInc, JacobianDiagBlock< nVarsFixed > &JDiag, ArrayDOFV< nVarsFixed > &uIncNew) |
| Deprecated: use UpdateSGS with uIncIsZero=true instead. Kept for backward compatibility. Delegates to UpdateSGS. | |
| void | UpdateSGS (real alphaDiag, real t, ArrayDOFV< nVarsFixed > &rhs, ArrayDOFV< nVarsFixed > &u, ArrayDOFV< nVarsFixed > &uInc, ArrayDOFV< nVarsFixed > &uIncNew, JacobianDiagBlock< nVarsFixed > &JDiag, bool forward, bool gsUpdate, TU &sumInc, bool uIncIsZero=false) |
| Symmetric Gauss-Seidel update for the implicit linear system. | |
| void | LUSGSMatrixSolveJacobianLU (real alphaDiag, real t, ArrayDOFV< nVarsFixed > &rhs, ArrayDOFV< nVarsFixed > &u, ArrayDOFV< nVarsFixed > &uInc, ArrayDOFV< nVarsFixed > &uIncNew, ArrayDOFV< nVarsFixed > &bBuf, JacobianDiagBlock< nVarsFixed > &JDiag, JacobianLocalLU< nVarsFixed > &jacLU, bool uIncIsZero, TU &sumInc) |
| Solve the implicit linear system using the pre-factored local LU. | |
| void | UpdateSGSWithRec (real alphaDiag, real t, ArrayDOFV< nVarsFixed > &rhs, ArrayDOFV< nVarsFixed > &u, ArrayRECV< nVarsFixed > &uRec, ArrayDOFV< nVarsFixed > &uInc, ArrayRECV< nVarsFixed > &uRecInc, JacobianDiagBlock< nVarsFixed > &JDiag, bool forward, TU &sumInc) |
| SGS sweep coupled with reconstruction update. | |
| void | FixUMaxFilter (ArrayDOFV< nVarsFixed > &u) |
| Clip extreme conserved-variable values to prevent overflow. | |
| void | TimeAverageAddition (ArrayDOFV< nVarsFixed > &w, ArrayDOFV< nVarsFixed > &wAveraged, real dt, real &tCur) |
| Accumulate time-averaged primitive variables for unsteady statistics. | |
| void | MeanValueCons2Prim (ArrayDOFV< nVarsFixed > &u, ArrayDOFV< nVarsFixed > &w) |
| Convert cell-mean conservative variables to primitive variables. | |
| void | MeanValuePrim2Cons (ArrayDOFV< nVarsFixed > &w, ArrayDOFV< nVarsFixed > &u) |
| Convert cell-mean primitive variables to conservative variables. | |
| void | EvaluateNorm (Eigen::Vector< real, -1 > &res, ArrayDOFV< nVarsFixed > &rhs, index P=1, bool volWise=false, bool average=false) |
| Compute the norm of the RHS residual vector. | |
| void | EvaluateRecNorm (Eigen::Vector< real, -1 > &res, ArrayDOFV< nVarsFixed > &u, ArrayRECV< nVarsFixed > &uRec, index P=1, bool compare=false, const tFCompareField &FCompareField=[](const Geom::tPoint &p, real t) { return TU::Zero();}, const tFCompareFieldWeight &FCompareFieldWeight=[](const Geom::tPoint &p, real t) { return 1.0;}, real t=0) |
| Compute the reconstruction error norm (optionally against an analytical field). | |
| void | LimiterUGrad (ArrayDOFV< nVarsFixed > &u, ArrayGRADV< nVarsFixed, gDim > &uGrad, ArrayGRADV< nVarsFixed, gDim > &uGradNew, uint64_t flags=LIMITER_UGRAD_No_Flags) |
| Apply slope limiter to the gradient field. | |
| void | EvaluateURecBeta (ArrayDOFV< nVarsFixed > &u, ArrayRECV< nVarsFixed > &uRec, ArrayDOFV< 1 > &uRecBeta, index &nLim, real &betaMin, int flag) |
| Evaluate the positivity-preserving reconstruction limiter (beta). | |
| bool | AssertMeanValuePP (ArrayDOFV< nVarsFixed > &u, bool panic) |
| Assert that all cell-mean values are physically realizable. | |
| void | EvaluateCellRHSAlpha (ArrayDOFV< nVarsFixed > &u, ArrayRECV< nVarsFixed > &uRec, ArrayDOFV< 1 > &uRecBeta, ArrayDOFV< nVarsFixed > &res, ArrayDOFV< 1 > &cellRHSAlpha, index &nLim, real &alphaMin, real relax, int compress=1, int flag=0) |
| Compute the positivity-preserving RHS scaling factor (alpha) per cell. | |
| void | EvaluateCellRHSAlphaExpansion (ArrayDOFV< nVarsFixed > &u, ArrayRECV< nVarsFixed > &uRec, ArrayDOFV< 1 > &uRecBeta, ArrayDOFV< nVarsFixed > &res, ArrayDOFV< 1 > &cellRHSAlpha, index &nLim, real alphaMin) |
| Expand a previously computed cellRHSAlpha toward 1 where safe. | |
| void | MinSmoothDTau (ArrayDOFV< 1 > &dTau, ArrayDOFV< 1 > &dTauNew) |
| Smooth the local time step across neighboring cells. | |
| real | muEff (const TU &U, real T) |
| Compute effective molecular viscosity using the configured viscosity model. | |
| real | getMuTur (const TU &uMean, const TDiffU &GradUMeanXY, real muRef, real muf, index iFace) |
| Compute turbulent eddy viscosity at a face. | |
| void | visFluxTurVariable (const TU &UMeanXYC, const TDiffU &DiffUxyPrimC, real muRef, real mufPhy, real muTur, const TVec &uNormC, index iFace, TU &VisFlux) |
| Compute viscous flux contribution from turbulence model variables. | |
| void | UFromCell2Face (TU &u, index iFace, index iCell, rowsize if2c) |
| Transform a conservative state vector from cell frame to face frame for periodic BCs. | |
| void | UFromFace2Cell (TU &u, index iFace, index iCell, rowsize if2c) |
| Inverse of UFromCell2Face: transform from face frame back to cell frame. | |
| void | UFromOtherCell (TU &u, index iFace, index iCell, index iCellOther, rowsize if2c) |
| Transform a state from a neighbor cell across a periodic face. | |
| void | DiffUFromCell2Face (TDiffU &u, index iFace, index iCell, rowsize if2c, bool reverse=false) |
| Transform a gradient tensor from cell frame to face frame for periodic BCs. | |
| void | fluxFace (const TU_Batch &ULxy, const TU_Batch &URxy, const TU &ULMeanXy, const TU &URMeanXy, const TDiffU_Batch &DiffUxy, const TDiffU_Batch &DiffUxyPrim, const TVec_Batch &unitNorm, const TVec_Batch &vg, const TVec &unitNormC, const TVec &vgC, TU_Batch &FLfix, TU_Batch &FRfix, TU_Batch &finc, TReal_Batch &lam0V, TReal_Batch &lam123V, TReal_Batch &lam4V, Geom::t_index btype, typename Gas::RiemannSolverType rsType, index iFace, bool ignoreVis) |
| Compute the numerical flux at a face (batched over quadrature points). | |
| TU | source (const TU &UMeanXy, const TDiffU &DiffUxy, const Geom::tPoint &pPhy, TJacobianU &jacobian, index iCell, index ig, int Mode) |
| Compute the source term at a cell quadrature point. | |
| auto | fluxJacobian0_Right (TU &UR, const TVec &uNorm, Geom::t_index btype) |
| inviscid flux approx jacobian (flux term not reconstructed / no riemann) | |
| TU | fluxJacobian0_Right_Times_du (const TU &U, const TU &UOther, const TVec &n, const TVec &vg, Geom::t_index btype, const TU &dU, real lambdaMain, real lambdaC, real lambdaVis, real lambda0, real lambda123, real lambda4, int useRoeTerm, int incFsign=-1, int omitF=0) |
| inviscid flux approx jacobian (flux term not reconstructed / no riemann) if lambdaMain == veryLargeReal, then use lambda0~4 for roe-flux type jacobian | |
| TJacobianU | fluxJacobian0_Right_Times_du_AsMatrix (const TU &U, const TU &UOther, const TVec &n, const TVec &vg, Geom::t_index btype, real lambdaMain, real lambdaC, real lambdaVis, real lambda0, real lambda123, real lambda4, int useRoeTerm, int incFsign=-1, int omitF=0) |
| TU | fluxJacobianC_Right_Times_du (const TU &U, const TVec &n, const TVec &vg, Geom::t_index btype, const TU &dU) |
| TVec | GetFaceVGrid (index iFace, index iG) |
| Get the grid velocity at a face quadrature point (rotating frame). | |
| TVec | GetFaceVGrid (index iFace, index iG, const Geom::tPoint &pPhy) |
| Get the grid velocity at a face quadrature point (with explicit physical point). | |
| TVec | GetFaceVGridFromCell (index iFace, index iCell, int if2c, index iG) |
| Get the grid velocity at a face quadrature point from a specific cell's perspective. | |
| void | TransformVelocityRotatingFrame (TU &U, const Geom::tPoint &pPhysics, int direction) |
| Transform momentum in-place between inertial and rotating frame (velocity only). | |
| void | TransformURotatingFrame (TU &U, const Geom::tPoint &pPhysics, int direction) |
| Transform full conservative state (momentum + total energy) between frames (relative velocity formulation). | |
| void | TransformURotatingFrame_ABS_VELO (TU &U, const Geom::tPoint &pPhysics, int direction) |
| Transform full conservative state for the absolute-velocity rotating frame formulation. | |
| void | updateBCAnchors (ArrayDOFV< nVarsFixed > &u, ArrayRECV< nVarsFixed > &uRec) |
| Update boundary anchor point recorders from current solution. | |
| void | updateBCProfiles (ArrayDOFV< nVarsFixed > &u, ArrayRECV< nVarsFixed > &uRec) |
| Update boundary 1-D profiles from the current solution. | |
| void | updateBCProfilesPressureRadialEq () |
| Update boundary profiles using radial-equilibrium pressure extrapolation. | |
| TU | generateBoundaryValue (TU &ULxy, const TU &ULMeanXy, index iCell, index iFace, int iG, const TVec &uNorm, const TMat &normBase, const Geom::tPoint &pPhysics, real t, Geom::t_index btype, bool fixUL=false, int geomMode=0, int linMode=0) |
| Generate the ghost (boundary) state for a boundary face. | |
| void | PrintBCProfiles (const std::string &name, ArrayDOFV< nVarsFixed > &u, ArrayRECV< nVarsFixed > &uRec) |
| Write boundary profile data to CSV files (rank 0 only). | |
| void | ConsoleOutputBndIntegrations () |
| Print boundary integration results (force/flux) to console on rank 0. | |
| void | BndIntegrationLogWriteLine (const std::string &name, index step, index stage, index iter) |
| Append a line to the per-BC boundary integration CSV log files. | |
| std::tuple< real, real, real > | CLDriverGetIntegrationUpdate (index iter) |
| Query the CL driver for current lift/drag coefficients and update AoA. | |
| TU | CompressRecPart (const TU &umean, const TU &uRecInc, bool &compressed) |
| Compress a reconstruction increment to maintain positivity. | |
| TU | CompressInc (const TU &u, const TU &uInc) |
| Compress a solution increment to maintain positivity. | |
| void | FixIncrement (ArrayDOFV< nVarsFixed > &cx, ArrayDOFV< nVarsFixed > &cxInc, real alpha=1.0) |
| Apply CompressInc to every cell, modifying cxInc in-place. | |
| void | AddFixedIncrement (ArrayDOFV< nVarsFixed > &cx, ArrayDOFV< nVarsFixed > &cxInc, real alpha=1.0) |
| Add a positivity-compressed increment to the solution. | |
| void | CentralSmoothResidual (ArrayDOFV< nVarsFixed > &r, ArrayDOFV< nVarsFixed > &rs, ArrayDOFV< nVarsFixed > &rtemp, int nStep=0) |
| Apply Laplacian smoothing to a residual field. | |
| void | InitializeUDOF (ArrayDOFV< nVarsFixed > &u) |
| Set initial conservative DOF values for all cells. | |
| void | InitializeOutputPicker (OutputPicker &op, OutputOverlapDataRefs dataRefs) |
| Initialize an OutputPicker with field callbacks for VTK/HDF5 output. | |
Static Public Member Functions | |
| static constexpr int | MaxBatchMult (int n) |
| Compute the maximum batch-multiplied column count for batched Eigen matrices. | |
| static void | InitializeFV (ssp< Geom::UnstructuredMesh > &mesh, TpVFV vfv, TpBCHandler pBCHandler) |
| Initialize the finite-volume infrastructure on the mesh. | |
Public Attributes | |
| ssp< Geom::UnstructuredMesh > | mesh |
| Shared pointer to the unstructured mesh. | |
| ssp< CFV::VariationalReconstruction< gDim > > | vfv |
| ssp< BoundaryHandler< model > > | pBCHandler |
| gDim -> 3 for intellisense //!tmptmp ///< Variational reconstruction object. | |
| int | kAv = 0 |
| Artificial viscosity polynomial order (maxOrder + 1). | |
| std::vector< real > | lambdaFace |
| Per-face spectral radius (inviscid + viscous combined). | |
| std::vector< real > | lambdaFaceC |
| Per-face convective spectral radius. | |
| std::vector< real > | lambdaFaceVis |
| Per-face viscous spectral radius. | |
| std::vector< real > | lambdaFace0 |
| Per-face eigenvalue |u·n| (contact wave). | |
| std::vector< real > | lambdaFace123 |
| Per-face eigenvalue |u·n + a| (acoustic wave). | |
| std::vector< real > | lambdaFace4 |
| Per-face eigenvalue |u·n - a| (acoustic wave). | |
| std::vector< real > | deltaLambdaFace |
| Per-face spectral radius difference for implicit diagonal. | |
| ArrayDOFV< 1 > | deltaLambdaCell |
| Per-cell accumulated spectral radius difference. | |
| std::vector< TDiffU > | gradUFix |
| Green-Gauss gradient correction buffer for source term stabilization. | |
| std::vector< Eigen::Vector< real, Eigen::Dynamic > > | dWall |
| Per-cell wall distance (one value per cell node/quadrature point). | |
| std::vector< real > | dWallFace |
| Per-face wall distance (interpolated from cell values). | |
| std::map< Geom::t_index, AnchorPointRecorder< nVarsFixed > > | anchorRecorders |
| Per-BC anchor point recorders for profile extraction. | |
| std::map< Geom::t_index, OneDimProfile< nVarsFixed > > | profileRecorders |
| Per-BC 1-D profile recorders. | |
| std::map< Geom::t_index, IntegrationRecorder > | bndIntegrations |
| Per-BC boundary flux/force integration accumulators. | |
| std::map< Geom::t_index, std::ofstream > | bndIntegrationLogs |
| Per-BC log file streams for integration output. | |
| std::set< Geom::t_index > | cLDriverBndIDs |
| Boundary IDs driven by the CL (lift-coefficient) driver. | |
| std::unique_ptr< CLDriver > | pCLDriver |
| Lift-coefficient driver for AoA adaptation (null if unused). | |
| ArrayGRADV< nVarsFixed, gDim > | uGradBuf |
| ArrayGRADV< nVarsFixed, gDim > | uGradBufNoLim |
| Gradient buffers: limited and unlimited. | |
| Eigen::Vector< real, -1 > | fluxWallSum |
| Accumulated wall flux integral (force on wall boundaries). | |
| std::vector< TU > | fluxBnd |
| Per-boundary-face flux values. | |
| std::vector< TVec > | fluxBndForceT |
| Per-boundary-face tangential force. | |
| index | nFaceReducedOrder = 0 |
| Count of faces where reconstruction order was reduced. | |
| ssp< Direct::SerialSymLUStructure > | symLU |
| Symmetric LU structure for direct preconditioner. | |
| EulerEvaluatorSettings< model > | settings |
| Physics and numerics settings for this evaluator. | |
Static Public Attributes | |
| static const int | nVarsFixed = getnVarsFixed(model) |
| Number of conserved variables (compile-time fixed). | |
| static const int | dim = getDim_Fixed(model) |
| Spatial dimension (2 or 3). | |
| static const int | gDim = getGeomDim_Fixed(model) |
| Geometric dimension of the mesh. | |
| static const auto | I4 = dim + 1 |
| Index of the energy equation (= dim+1). | |
| static const int | MaxBatch = 16 |
| static const uint64_t | DT_No_Flags = 0x0ull |
| No flags for EvaluateDt. | |
| static const uint64_t | DT_Dont_update_lambda01234 = 0x1ull << 0 |
| Skip recomputation of per-face eigenvalues lambda0/123/4. | |
| static const int | EvaluateURecBeta_DEFAULT = 0x00 |
| Default: evaluate beta without compression. | |
| static const int | EvaluateURecBeta_COMPRESS_TO_MEAN = 0x01 |
| Compress reconstruction toward cell mean to enforce positivity. | |
| static const int | EvaluateCellRHSAlpha_DEFAULT = 0x00 |
| Default alpha evaluation mode. | |
| static const int | EvaluateCellRHSAlpha_MIN_IF_NOT_ONE = 0x01 |
| Take min(alpha, prev) only if prev != 1. | |
| static const int | EvaluateCellRHSAlpha_MIN_ALL = 0x02 |
| Always take min(alpha, prev). | |
RHS evaluation flags (bitwise OR combinable) | |
| static const uint64_t | RHS_No_Flags = 0x0ull |
| Default: full RHS evaluation. | |
| static const uint64_t | RHS_Ignore_Viscosity = 0x1ull << 0 |
| Skip viscous flux contribution. | |
| static const uint64_t | RHS_Dont_Update_Integration = 0x1ull << 1 |
| Skip boundary integration accumulation. | |
| static const uint64_t | RHS_Dont_Record_Bud_Flux = 0x1ull << 2 |
| Skip recording per-boundary flux. | |
| static const uint64_t | RHS_Direct_2nd_Rec = 0x1ull << 8 |
| Use direct 2nd-order reconstruction. | |
| static const uint64_t | RHS_Direct_2nd_Rec_1st_Conv = 0x1ull << 9 |
| 2nd-order rec with 1st-order convection. | |
| static const uint64_t | RHS_Direct_2nd_Rec_use_limiter = 0x1ull << 10 |
| Apply limiter when using direct 2nd rec. | |
| static const uint64_t | RHS_Direct_2nd_Rec_already_have_uGradBufNoLim = 0x1ull << 11 |
| uGradBufNoLim is already computed. | |
| static const uint64_t | RHS_Recover_IncFScale = 0x1ull << 12 |
| Recover incremental face scaling. | |
Limiter flags | |
| static const uint64_t | LIMITER_UGRAD_No_Flags = 0x0ull |
| Default limiter flags. | |
| static const uint64_t | LIMITER_UGRAD_Disable_Shock_Limiter = 0x1ull << 0 |
| Disable shock-detecting component of the limiter. | |
Core finite-volume evaluator for compressible Navier-Stokes / Euler equations.
Implements the spatial discretization of the compressible N-S equations on unstructured meshes using Compact Finite Volume (CFV) with variational reconstruction. Handles inviscid and viscous flux evaluation, source terms (including RANS turbulence models), implicit Jacobian assembly, boundary condition ghost-value generation, and positivity-preserving fixes.
| model | EulerModel enum selecting the equation set and spatial dimension (NS, NS_2D, NS_SA, NS_SA_3D, NS_2EQ, NS_2EQ_3D, NS_3D). |
Definition at line 77 of file EulerEvaluator.hpp.
| typedef Eigen::MatrixFMTSafe<real, dim, nVarsFixed> DNDS::Euler::EulerEvaluator< model >::TDiffU |
Gradient of conserved variables (dim x nVars).
Definition at line 98 of file EulerEvaluator.hpp.
| typedef Eigen::MatrixFMTSafe<real, Eigen::Dynamic, nVarsFixed, Eigen::ColMajor, MaxBatchMult(3)> DNDS::Euler::EulerEvaluator< model >::TDiffU_Batch |
Batch of gradient matrices.
Definition at line 99 of file EulerEvaluator.hpp.
| typedef Eigen::MatrixFMTSafe<real, nVarsFixed, dim> DNDS::Euler::EulerEvaluator< model >::TDiffUTransposed |
Transposed gradient (nVars x dim).
Definition at line 100 of file EulerEvaluator.hpp.
| typedef ArrayDOFV<nVarsFixed> DNDS::Euler::EulerEvaluator< model >::TDof |
Cell-centered DOF array (mean values).
Definition at line 101 of file EulerEvaluator.hpp.
| using DNDS::Euler::EulerEvaluator< model >::tFCompareField = std::function<TU(const Geom::tPoint &, real)> |
Callback type for analytical comparison field.
Definition at line 601 of file EulerEvaluator.hpp.
| using DNDS::Euler::EulerEvaluator< model >::tFCompareFieldWeight = std::function<real(const Geom::tPoint &, real)> |
Callback type for comparison weighting function.
Definition at line 602 of file EulerEvaluator.hpp.
| typedef Eigen::MatrixFMTSafe<real, nVarsFixed, nVarsFixed> DNDS::Euler::EulerEvaluator< model >::TJacobianU |
Jacobian matrix (nVars x nVars) of the flux w.r.t. conserved variables.
Definition at line 97 of file EulerEvaluator.hpp.
| typedef Eigen::MatrixFMTSafe<real, dim, dim> DNDS::Euler::EulerEvaluator< model >::TMat |
Spatial matrix (dim x dim).
Definition at line 92 of file EulerEvaluator.hpp.
| typedef Eigen::MatrixFMTSafe<real, dim, Eigen::Dynamic, Eigen::ColMajor, dim, MaxBatchMult(3)> DNDS::Euler::EulerEvaluator< model >::TMat_Batch |
Batch of spatial matrices.
Definition at line 93 of file EulerEvaluator.hpp.
| typedef ssp<BoundaryHandler<model> > DNDS::Euler::EulerEvaluator< model >::TpBCHandler |
Shared pointer to the boundary condition handler.
Definition at line 107 of file EulerEvaluator.hpp.
| typedef ssp<CFV::VariationalReconstruction<gDim> > DNDS::Euler::EulerEvaluator< model >::TpVFV |
Shared pointer to the variational reconstruction object.
Definition at line 106 of file EulerEvaluator.hpp.
| using DNDS::Euler::EulerEvaluator< model >::Traits = EulerModelTraits<model> |
Compile-time traits: hasSA, has2EQ, etc.
Definition at line 80 of file EulerEvaluator.hpp.
| typedef Eigen::MatrixFMTSafe<real, 1, Eigen::Dynamic, Eigen::RowMajor, 1, MaxBatch> DNDS::Euler::EulerEvaluator< model >::TReal_Batch |
Batch of scalar values.
Definition at line 96 of file EulerEvaluator.hpp.
| typedef ArrayRECV<nVarsFixed> DNDS::Euler::EulerEvaluator< model >::TRec |
Reconstruction coefficient array (per-cell polynomial coefficients).
Definition at line 102 of file EulerEvaluator.hpp.
| typedef ArrayRECV<1> DNDS::Euler::EulerEvaluator< model >::TScalar |
Scalar reconstruction coefficient array.
Definition at line 103 of file EulerEvaluator.hpp.
| typedef Eigen::VectorFMTSafe<real, nVarsFixed> DNDS::Euler::EulerEvaluator< model >::TU |
Conservative variable vector (nVarsFixed components).
Definition at line 94 of file EulerEvaluator.hpp.
| typedef Eigen::MatrixFMTSafe<real, nVarsFixed, Eigen::Dynamic, Eigen::ColMajor, nVarsFixed, MaxBatch> DNDS::Euler::EulerEvaluator< model >::TU_Batch |
Batch of conservative variable vectors.
Definition at line 95 of file EulerEvaluator.hpp.
| typedef Eigen::VectorFMTSafe<real, dim> DNDS::Euler::EulerEvaluator< model >::TVec |
Spatial vector (dim components).
Definition at line 90 of file EulerEvaluator.hpp.
| typedef Eigen::MatrixFMTSafe<real, dim, Eigen::Dynamic, Eigen::ColMajor, dim, MaxBatch> DNDS::Euler::EulerEvaluator< model >::TVec_Batch |
Batch of spatial vectors.
Definition at line 91 of file EulerEvaluator.hpp.
| typedef CFV::VariationalReconstruction<gDim> DNDS::Euler::EulerEvaluator< model >::TVFV |
Variational reconstruction type for this geometric dimension.
Definition at line 105 of file EulerEvaluator.hpp.
|
inline |
Construct an EulerEvaluator and initialize all internal buffers.
Allocates per-face spectral-radius buffers, per-boundary flux arrays, gradient correction arrays (if enabled), wall distance fields (for RANS), and the symmetric LU structure for direct preconditioning. Also validates CL-driver boundary configuration.
| Nmesh | Shared pointer to the unstructured mesh. |
| Nvfv | Shared pointer to the variational reconstruction object. |
| npBCHandler | Shared pointer to the boundary condition handler. |
| nSettings | Physics and numerics settings. |
| n_nVars | Runtime number of conserved variables (default from model). |
Definition at line 236 of file EulerEvaluator.hpp.
|
inline |
Add a positivity-compressed increment to the solution.
For each cell, compresses the increment via CompressInc, then adds it to cx. Reports the global minimum compression factor via MPI reduction.
| [in,out] | cx | Solution array (updated in-place). |
| [in] | cxInc | Increment array. |
| [in] | alpha | Scaling factor applied to the increment before compression. |
old inc could be zero, so compresion alpha is always 1
Definition at line 1795 of file EulerEvaluator.hpp.
| bool DNDS::Euler::EulerEvaluator< model >::AssertMeanValuePP | ( | ArrayDOFV< nVarsFixed > & | u, |
| bool | panic | ||
| ) |
Assert that all cell-mean values are physically realizable.
Checks that density > 0 and internal energy > 0 for all cells.
| [in] | u | Cell-centered conservative DOFs. |
| [in] | panic | If true, abort on first violation; otherwise just report. |
Definition at line 1866 of file EulerEvaluator.hxx.
|
inline |
Append a line to the per-BC boundary integration CSV log files.
only 0 needs to write
Definition at line 1466 of file EulerEvaluator.hpp.
|
inline |
Apply Laplacian smoothing to a residual field.
Iteratively smooths r into rs using weighted neighbor averaging. Used to stabilize central-difference residuals on coarse grids.
| [in] | r | Input residual. |
| [in,out] | rs | Smoothed residual (output; also used as work buffer). |
| [out] | rtemp | Temporary buffer (same size as r). |
| [in] | nStep | Number of smoothing passes (0 = use settings.nCentralSmoothStep). |
Definition at line 1923 of file EulerEvaluator.hpp.
|
inline |
Query the CL driver for current lift/drag coefficients and update AoA.
Collects boundary force integrals from CL-driven boundaries, projects them onto lift/drag directions, and feeds the result to the CLDriver for AoA adaptation.
| iter | Current iteration count. |
Definition at line 1502 of file EulerEvaluator.hpp.
|
inline |
Compress a solution increment to maintain positivity.
Given the current state u and an update increment uInc, returns a modified increment that ensures u + result has positive density, positive internal energy, and (for RANS) non-negative turbulent variables. Uses exponential decay clamping for density and a quadratic solve for internal energy.
| u | Current conservative state. |
| uInc | Proposed increment. |
A intuitive fix
need positive perserving technique!
A intuitive fix
might be bad using gmres, add this to gmres inc!
might be bad using gmres, add this to gmres inc!
might be bad using gmres, add this to gmres inc!
Definition at line 1645 of file EulerEvaluator.hpp.
|
inline |
Compress a reconstruction increment to maintain positivity.
Given a cell mean and reconstruction increment, returns umean + uRecInc with the increment clamped so that density, internal energy, and (for RANS) turbulent variables remain non-negative.
| umean | Cell-mean conservative state. |
| uRecInc | Reconstruction polynomial increment at an evaluation point. |
| compressed | Set to true if compression was applied. |
Definition at line 1541 of file EulerEvaluator.hpp.
|
inline |
Print boundary integration results (force/flux) to console on rank 0.
Definition at line 1447 of file EulerEvaluator.hpp.
|
inline |
Transform a gradient tensor from cell frame to face frame for periodic BCs.
Definition at line 928 of file EulerEvaluator.hpp.
| void DNDS::Euler::EulerEvaluator< model >::EvaluateCellRHSAlpha | ( | ArrayDOFV< nVarsFixed > & | u, |
| ArrayRECV< nVarsFixed > & | uRec, | ||
| ArrayDOFV< 1 > & | uRecBeta, | ||
| ArrayDOFV< nVarsFixed > & | res, | ||
| ArrayDOFV< 1 > & | cellRHSAlpha, | ||
| index & | nLim, | ||
| real & | alphaMin, | ||
| real | relax, | ||
| int | compress = 1, |
||
| int | flag = 0 |
||
| ) |
Compute the positivity-preserving RHS scaling factor (alpha) per cell.
Evaluate per-cell RHS scaling factor (alpha) for positivity-preserving time stepping.
Determines the maximum safe scaling alpha in [0,1] such that u + alpha * res remains physically realizable.
| [in] | u | Cell-centered conservative DOFs. |
| [in] | uRec | Reconstruction coefficients. |
| [in] | uRecBeta | Per-cell reconstruction beta from PP limiter. |
| res | Incremental residual (the RHS increment to scale). | |
| [out] | cellRHSAlpha | Per-cell alpha factor. |
| [out] | nLim | Number of cells where alpha < 1. |
| [out] | alphaMin | Global minimum alpha. |
| [in] | relax | Relaxation factor applied to alpha. |
| [in] | compress | Compression mode (1=compress, 0=clamp). |
| [in] | flag | EvaluateCellRHSAlpha_* mode flag. |
Determines the largest alpha in [0,1] such that u + alpha*res remains physically realizable (positive density and pressure) at all quadrature points. Used for under-relaxation to preserve positivity during explicit or implicit updates.
| u | Conservative variable DOF array. |
| uRec | Reconstruction coefficients. |
| uRecBeta | Per-cell reconstruction limiter coefficient. |
| res | RHS residual to be scaled. |
| cellRHSAlpha | Per-cell alpha scaling factor (output). |
| nLim | Total number of limited cells (output, MPI-reduced). |
| alphaMin | Global minimum alpha (output, MPI-reduced). |
| relax | Relaxation factor applied to alpha. |
| compress | Compression mode for reconstruction part. |
| flag | Evaluation mode flags (e.g., EvaluateCellRHSAlpha_MIN_ALL). |
for safety
Definition at line 1938 of file EulerEvaluator.hxx.
| void DNDS::Euler::EulerEvaluator< model >::EvaluateCellRHSAlphaExpansion | ( | ArrayDOFV< nVarsFixed > & | u, |
| ArrayRECV< nVarsFixed > & | uRec, | ||
| ArrayDOFV< 1 > & | uRecBeta, | ||
| ArrayDOFV< nVarsFixed > & | res, | ||
| ArrayDOFV< 1 > & | cellRHSAlpha, | ||
| index & | nLim, | ||
| real | alphaMin | ||
| ) |
Expand a previously computed cellRHSAlpha toward 1 where safe.
Expand the per-cell alpha limiter to neighboring cells for smoothness.
| res | Incremental residual fixed previously. | |
| cellRHSAlpha | Limiting factor evaluated previously (expanded in-place). | |
| [out] | nLim | Number of cells still limited after expansion. |
| [out] | alphaMin | Global minimum alpha after expansion. |
Cells adjacent to limited cells (alpha < 1) are also given reduced alpha values to create a smooth transition region, preventing isolated sharp limiting boundaries.
| u | Conservative variable DOF array. |
| uRec | Reconstruction coefficients. |
| uRecBeta | Per-cell reconstruction limiter coefficient. |
| res | RHS residual. |
| cellRHSAlpha | Per-cell alpha (input/output, expanded to neighbors). |
| nLim | Additional limited cells from expansion (output, accumulated). |
| alphaMin | Global minimum alpha (input). |
Definition at line 2064 of file EulerEvaluator.hxx.
| void DNDS::Euler::EulerEvaluator< model >::EvaluateDt | ( | ArrayDOFV< 1 > & | dt, |
| ArrayDOFV< nVarsFixed > & | u, | ||
| ArrayRECV< nVarsFixed > & | uRec, | ||
| real | CFL, | ||
| real & | dtMinall, | ||
| real | MaxDt, | ||
| bool | UseLocaldt, | ||
| real | t, | ||
| uint64_t | flags = DT_No_Flags |
||
| ) |
Estimate the local or global time step for each cell.
Computes the local pseudo-time step dTau based on CFL number, spectral radii of inviscid and viscous fluxes, and optionally enforces a global minimum. Updates per-face eigenvalue arrays (lambda0, lambda123, lambda4) unless DT_Dont_update_lambda01234 is set.
| [out] | dt | Per-cell time step array (overwritten). |
| [in] | u | Cell-centered conservative variable DOFs. |
| [in] | uRec | Reconstruction coefficients. |
| [in] | CFL | CFL number. |
| [out] | dtMinall | Global minimum time step (MPI-reduced). |
| [in] | MaxDt | Upper bound on the time step. |
| [in] | UseLocaldt | If true, use local (per-cell) time stepping. |
| [in] | t | Current simulation time. |
| [in] | flags | Bitwise combination of DT_* flags. |
| void DNDS::Euler::EulerEvaluator< model >::EvaluateNorm | ( | Eigen::Vector< real, -1 > & | res, |
| ArrayDOFV< nVarsFixed > & | rhs, | ||
| index | P = 1, |
||
| bool | volWise = false, |
||
| bool | average = false |
||
| ) |
Compute the norm of the RHS residual vector.
Evaluate the LP norm of the RHS residual (MPI-global).
| [out] | res | Norm result per variable (resized to nVars). |
| [in] | rhs | Cell residual array. |
| [in] | P | Norm order (1 = L1, 2 = L2, etc.). |
| [in] | volWise | If true, weight by cell volume. |
| [in] | average | If true, divide by total volume/count. |
For P < 3, computes component-wise LP norm; for P >= 3, computes L-infinity norm. Optionally volume-weighted and/or averaged.
| res | Result vector of nVars component norms (output). |
| rhs | Right-hand side residual. |
| P | Norm exponent (1 for L1, 2 for L2, >=3 for L-infinity). |
| volWise | If true, weight by cell volume. |
| average | If true, divide by total weight. |
Definition at line 1394 of file EulerEvaluator.hxx.
| void DNDS::Euler::EulerEvaluator< model >::EvaluateRecNorm | ( | Eigen::Vector< real, -1 > & | res, |
| ArrayDOFV< nVarsFixed > & | u, | ||
| ArrayRECV< nVarsFixed > & | uRec, | ||
| index | P = 1, |
||
| bool | compare = false, |
||
| const tFCompareField & | FCompareField = [](const Geom::tPoint &p, real t) { return TU::Zero(); }, |
||
| const tFCompareFieldWeight & | FCompareFieldWeight = [](const Geom::tPoint &p, real t) { return 1.0; }, |
||
| real | t = 0 |
||
| ) |
Compute the reconstruction error norm (optionally against an analytical field).
Evaluate the LP norm of the reconstructed solution, optionally against a reference field.
| [out] | res | Norm result per variable. |
| [in] | u | Cell-centered conservative DOFs. |
| [in] | uRec | Reconstruction coefficients. |
| [in] | P | Norm order. |
| [in] | compare | If true, compute error against FCompareField. |
| [in] | FCompareField | Analytical field callback. |
| [in] | FCompareFieldWeight | Weighting callback. |
| [in] | t | Current simulation time. |
Integrates the reconstructed solution (or its error against FCompareField) over quadrature points in each cell and reduces globally.
| res | Result vector of nVars component norms (output). |
| u | Conservative variable DOF array. |
| uRec | Reconstruction coefficients. |
| P | Norm exponent (1 for L1, 2 for L2, >=3 for L-infinity). |
| compare | If true, compute error against FCompareField. |
| FCompareField | Reference field function(point, time) -> TU. |
| FCompareFieldWeight | Weight function(point, time) -> TU for error weighting. |
| t | Current simulation time. |
Definition at line 1474 of file EulerEvaluator.hxx.
| void DNDS::Euler::EulerEvaluator< model >::EvaluateRHS | ( | ArrayDOFV< nVarsFixed > & | rhs, |
| JacobianDiagBlock< nVarsFixed > & | JSource, | ||
| ArrayDOFV< nVarsFixed > & | u, | ||
| ArrayRECV< nVarsFixed > & | uRecUnlim, | ||
| ArrayRECV< nVarsFixed > & | uRec, | ||
| ArrayDOFV< 1 > & | uRecBeta, | ||
| ArrayDOFV< 1 > & | cellRHSAlpha, | ||
| bool | onlyOnHalfAlpha, | ||
| real | t, | ||
| uint64_t | flags = RHS_No_Flags |
||
| ) |
Evaluate the spatial right-hand side of the semi-discrete system.
Computes inviscid flux (Riemann solver), viscous flux, and source terms over all internal and boundary faces, accumulating cell residuals. Also records boundary force/flux integrations and updates reduced-order face counts.
| [out] | rhs | Cell residual array (overwritten). |
| [out] | JSource | Diagonal Jacobian block from source terms. |
| [in] | u | Cell-centered conservative DOFs. |
| [in] | uRecUnlim | Unlimited reconstruction coefficients. |
| [in] | uRec | Limited reconstruction coefficients. |
| [in] | uRecBeta | Per-cell reconstruction compression factor (PP limiter). |
| [in] | cellRHSAlpha | Per-cell RHS scaling factor (PP limiter). |
| [in] | onlyOnHalfAlpha | If true, evaluate only cells with alpha < 1. |
| [in] | t | Current simulation time. |
| [in] | flags | Bitwise combination of RHS_* flags. |
| void DNDS::Euler::EulerEvaluator< model >::EvaluateURecBeta | ( | ArrayDOFV< nVarsFixed > & | u, |
| ArrayRECV< nVarsFixed > & | uRec, | ||
| ArrayDOFV< 1 > & | uRecBeta, | ||
| index & | nLim, | ||
| real & | betaMin, | ||
| int | flag | ||
| ) |
Evaluate the positivity-preserving reconstruction limiter (beta).
For each cell, computes the maximum compression factor beta such that u_mean + beta * uRec remains physically realizable (positive density and pressure).
| [in] | u | Cell-centered conservative DOFs. |
| [in] | uRec | Reconstruction coefficients. |
| [out] | uRecBeta | Per-cell compression factor in [0,1]. |
| [out] | nLim | Number of cells where beta < 1. |
| [out] | betaMin | Global minimum beta value. |
| [in] | flag | EvaluateURecBeta_DEFAULT or EvaluateURecBeta_COMPRESS_TO_MEAN. |
|
inline |
Apply CompressInc to every cell, modifying cxInc in-place.
Definition at line 1777 of file EulerEvaluator.hpp.
| void DNDS::Euler::EulerEvaluator< model >::FixUMaxFilter | ( | ArrayDOFV< nVarsFixed > & | u | ) |
Clip extreme conserved-variable values to prevent overflow.
Apply maximum-value spatial filter to conservative variables.
Placeholder for a spatial filter Jacobian approach. Currently returns immediately without modification.
| u | Conservative variable DOF array (unmodified). |
Definition at line 1318 of file EulerEvaluator.hxx.
| void DNDS::Euler::EulerEvaluator< model >::fluxFace | ( | const TU_Batch & | ULxy, |
| const TU_Batch & | URxy, | ||
| const TU & | ULMeanXy, | ||
| const TU & | URMeanXy, | ||
| const TDiffU_Batch & | DiffUxy, | ||
| const TDiffU_Batch & | DiffUxyPrim, | ||
| const TVec_Batch & | unitNorm, | ||
| const TVec_Batch & | vg, | ||
| const TVec & | unitNormC, | ||
| const TVec & | vgC, | ||
| TU_Batch & | FLfix, | ||
| TU_Batch & | FRfix, | ||
| TU_Batch & | finc, | ||
| TReal_Batch & | lam0V, | ||
| TReal_Batch & | lam123V, | ||
| TReal_Batch & | lam4V, | ||
| Geom::t_index | btype, | ||
| typename Gas::RiemannSolverType | rsType, | ||
| index | iFace, | ||
| bool | ignoreVis | ||
| ) |
Compute the numerical flux at a face (batched over quadrature points).
Evaluates the Riemann-solver-based inviscid flux and (optionally) the viscous flux at all quadrature points of a face simultaneously. Supports Roe, HLLC, Lax-Friedrichs, and other Riemann solvers selected by rsType.
| [in] | ULxy | Left state at quadrature points (batched). |
| [in] | URxy | Right state at quadrature points (batched). |
| [in] | ULMeanXy | Left cell mean state. |
| [in] | URMeanXy | Right cell mean state. |
| [in] | DiffUxy | Gradient of conservative variables at quad points. |
| [in] | DiffUxyPrim | Gradient of primitive variables at quad points. |
| [in] | unitNorm | Face outward unit normals at quad points. |
| [in] | vg | Grid velocity at quad points (for ALE / rotating frame). |
| [in] | unitNormC | Face-center outward unit normal. |
| [in] | vgC | Grid velocity at face center. |
| [out] | FLfix | Left-biased flux correction (batched). |
| [out] | FRfix | Right-biased flux correction (batched). |
| [out] | finc | Numerical flux increment (batched). |
| [out] | lam0V | Eigenvalue |u·n| at quad points. |
| [out] | lam123V | Eigenvalue |u·n + a| at quad points. |
| [out] | lam4V | Eigenvalue |u·n - a| at quad points. |
| [in] | btype | Boundary type (UnInitIndex for internal faces). |
| [in] | rsType | Riemann solver type (Roe, HLLC, LF, etc.). |
| [in] | iFace | Face index. |
| [in] | ignoreVis | If true, skip viscous flux computation. |
|
inline |
inviscid flux approx jacobian (flux term not reconstructed / no riemann)
Definition at line 1028 of file EulerEvaluator.hpp.
|
inline |
inviscid flux approx jacobian (flux term not reconstructed / no riemann) if lambdaMain == veryLargeReal, then use lambda0~4 for roe-flux type jacobian
now dF(U, dU) (GasInviscidFluxFacialIncrement) part actually treats the SeqI52Last part (as passive scalars) the RANS models and eulerEX use the same transport jacobian form
Definition at line 1100 of file EulerEvaluator.hpp.
|
inline |
|
inline |
| TU DNDS::Euler::EulerEvaluator< model >::generateBoundaryValue | ( | TU & | ULxy, |
| const TU & | ULMeanXy, | ||
| index | iCell, | ||
| index | iFace, | ||
| int | iG, | ||
| const TVec & | uNorm, | ||
| const TMat & | normBase, | ||
| const Geom::tPoint & | pPhysics, | ||
| real | t, | ||
| Geom::t_index | btype, | ||
| bool | fixUL = false, |
||
| int | geomMode = 0, |
||
| int | linMode = 0 |
||
| ) |
Generate the ghost (boundary) state for a boundary face.
Dispatches to type-specific BC handlers (far-field, wall, outflow, inflow, symmetry, etc.) based on btype. Used by the RHS evaluator to obtain the right state at boundary faces for the Riemann solver.
| [in,out] | ULxy | Left (interior) state; may be modified for certain BCs. |
| [in] | ULMeanXy | Left cell-mean state (base state for linearized mode). |
| [in] | iCell | Cell index adjacent to the boundary face. |
| [in] | iFace | Face index. |
| [in] | iG | Quadrature point index (< -1 for arbitrary location). |
| [in] | uNorm | Face outward unit normal. |
| [in] | normBase | Orthonormal basis with uNorm as first column. |
| [in] | pPhysics | Physical coordinates of the evaluation point. |
| [in] | t | Current simulation time. |
| [in] | btype | Boundary type ID. |
| [in] | fixUL | If true, do not modify the left state. |
| [in] | geomMode | Geometry evaluation mode (0=standard). |
| [in] | linMode | Linearization mode (0=nonlinear, nonzero=linearized about ULMeanXy). |
| ULMeanXy | warning, possible that UL is also modified |
|
inline |
Return the axisymmetric mode flag.
Definition at line 173 of file EulerEvaluator.hpp.
|
inline |
Get the grid velocity at a face quadrature point (rotating frame).
Definition at line 1208 of file EulerEvaluator.hpp.
|
inline |
Get the grid velocity at a face quadrature point (with explicit physical point).
Definition at line 1221 of file EulerEvaluator.hpp.
|
inline |
Get the grid velocity at a face quadrature point from a specific cell's perspective.
Definition at line 1235 of file EulerEvaluator.hpp.
|
inline |
Compute turbulent eddy viscosity at a face.
Dispatches to the appropriate RANS model (SA, k-omega SST, k-omega Wilcox, or Realizable k-epsilon) based on settings.ransModel.
| uMean | Cell-mean conservative state. |
| GradUMeanXY | Gradient of conservative variables in physical coordinates. |
| muRef | Reference dynamic viscosity scaling. |
| muf | Molecular (physical) viscosity at the face. |
| iFace | Face index (for wall distance lookup). |
Definition at line 791 of file EulerEvaluator.hpp.
| void DNDS::Euler::EulerEvaluator< model >::GetWallDist | ( | ) |
Compute wall distance for all cells (dispatcher).
Selects the wall distance algorithm based on the wallDistScheme setting (AABB tree, batched AABB, p-Poisson, or combination) and populates dWall and dWallFace arrays used by RANS turbulence models.
|
inlinestatic |
Initialize the finite-volume infrastructure on the mesh.
Sets periodic transformations on the VFV object, constructs cell/face geometric metrics, builds reconstruction base functions and weights (with BC-type-dependent Dirichlet/Neumann weighting), and computes reconstruction coefficients.
| mesh | Shared pointer to the unstructured mesh. |
| vfv | Shared pointer to the variational reconstruction object. |
| pBCHandler | Shared pointer to the boundary condition handler (used for per-face reconstruction weight selection). |
important caveat! using & captures mesh shared_ptr as reference, lost if inbound pointer is destoried!!
treat as real internal
Definition at line 123 of file EulerEvaluator.hpp.
| void DNDS::Euler::EulerEvaluator< model >::InitializeOutputPicker | ( | OutputPicker & | op, |
| OutputOverlapDataRefs | dataRefs | ||
| ) |
Initialize an OutputPicker with field callbacks for VTK/HDF5 output.
Register available output scalar fields in the OutputPicker.
Populates the output picker with named lambda functions that extract per-cell scalar quantities (conservative components, reconstruction coefficients, limiter values, wall distance, cell volume, turbulent viscosity, etc.) for use in post-processing output.
| op | OutputPicker to populate with named field extractors (output). |
| dataRefs | References to the solution arrays (u, uRec, betaPP, alphaPP). |
Definition at line 2468 of file EulerEvaluator_EvaluateDt.hxx.
| void DNDS::Euler::EulerEvaluator< model >::InitializeUDOF | ( | ArrayDOFV< nVarsFixed > & | u | ) |
Set initial conservative DOF values for all cells.
Populates u with the far-field state (or a problem-specific initial condition) based on the evaluator settings. Handles rotating-frame velocity transformations.
| [out] | u | Cell-centered DOF array to initialize. |
| void DNDS::Euler::EulerEvaluator< model >::LimiterUGrad | ( | ArrayDOFV< nVarsFixed > & | u, |
| ArrayGRADV< nVarsFixed, gDim > & | uGrad, | ||
| ArrayGRADV< nVarsFixed, gDim > & | uGradNew, | ||
| uint64_t | flags = LIMITER_UGRAD_No_Flags |
||
| ) |
Apply slope limiter to the gradient field.
Limits the reconstructed gradient (uGrad) to produce a monotonicity-preserving gradient (uGradNew). Supports WBAP and CWBAP limiter variants.
| [in] | u | Cell-centered conservative DOFs. |
| [in] | uGrad | Input (unlimited) gradient array. |
| [out] | uGradNew | Output limited gradient array. |
| [in] | flags | Bitwise combination of LIMITER_UGRAD_* flags. |
| void DNDS::Euler::EulerEvaluator< model >::LUSGSMatrixInit | ( | JacobianDiagBlock< nVarsFixed > & | JDiag, |
| JacobianDiagBlock< nVarsFixed > & | JSource, | ||
| ArrayDOFV< 1 > & | dTau, | ||
| real | dt, | ||
| real | alphaDiag, | ||
| ArrayDOFV< nVarsFixed > & | u, | ||
| ArrayRECV< nVarsFixed > & | uRec, | ||
| int | jacobianCode, | ||
| real | t | ||
| ) |
Assemble the diagonal blocks of the implicit Jacobian for LU-SGS / SGS.
Computes J_diag = (V/dTau + alphaDiag * sum_faces(spectral_radius)) * I + J_source for each cell, where V is cell volume and dTau is the local pseudo-time step.
| [out] | JDiag | Per-cell diagonal Jacobian block (overwritten). |
| [in] | JSource | Source-term Jacobian contribution to diagonal. |
| [in] | dTau | Per-cell local pseudo-time step. |
| [in] | dt | Physical time step (for dual time stepping). |
| [in] | alphaDiag | Diagonal scaling factor for implicit relaxation. |
| [in] | u | Cell-centered conservative DOFs. |
| [in] | uRec | Reconstruction coefficients. |
| [in] | jacobianCode | Controls Jacobian approximation: 0=scalar, 1=analytical flux Jacobian. |
| [in] | t | Current simulation time. |
| void DNDS::Euler::EulerEvaluator< model >::LUSGSMatrixSolveJacobianLU | ( | real | alphaDiag, |
| real | t, | ||
| ArrayDOFV< nVarsFixed > & | rhs, | ||
| ArrayDOFV< nVarsFixed > & | u, | ||
| ArrayDOFV< nVarsFixed > & | uInc, | ||
| ArrayDOFV< nVarsFixed > & | uIncNew, | ||
| ArrayDOFV< nVarsFixed > & | bBuf, | ||
| JacobianDiagBlock< nVarsFixed > & | JDiag, | ||
| JacobianLocalLU< nVarsFixed > & | jacLU, | ||
| bool | uIncIsZero, | ||
| TU & | sumInc | ||
| ) |
Solve the implicit linear system using the pre-factored local LU.
| [in] | alphaDiag | Diagonal scaling factor. |
| [in] | t | Current simulation time. |
| [in] | rhs | Right-hand side residual. |
| [in] | u | Cell-centered conservative DOFs. |
| [in] | uInc | Current increment (input guess). |
| [out] | uIncNew | Updated increment (output). |
| [out] | bBuf | Buffer for intermediate right-hand side assembly. |
| [in] | JDiag | Diagonal Jacobian blocks. |
| [in] | jacLU | Pre-factored local LU decomposition. |
| [in] | uIncIsZero | If true, skip contributions from zero-increment cells. |
| [out] | sumInc | Accumulated absolute increment for convergence monitoring. |
| void DNDS::Euler::EulerEvaluator< model >::LUSGSMatrixToJacobianLU | ( | real | alphaDiag, |
| real | t, | ||
| ArrayDOFV< nVarsFixed > & | u, | ||
| JacobianDiagBlock< nVarsFixed > & | JDiag, | ||
| JacobianLocalLU< nVarsFixed > & | jacLU | ||
| ) |
Build the local LU factorization of the Jacobian for direct solve.
Assembles the full local Jacobian (including off-diagonal face coupling) and factorizes it, storing the result in jacLU for subsequent direct solves.
| [in] | alphaDiag | Diagonal scaling factor. |
| [in] | t | Current simulation time. |
| [in] | u | Cell-centered conservative DOFs. |
| [in] | JDiag | Pre-assembled diagonal Jacobian blocks. |
| [out] | jacLU | Local LU factorization result (overwritten). |
| void DNDS::Euler::EulerEvaluator< model >::LUSGSMatrixVec | ( | real | alphaDiag, |
| real | t, | ||
| ArrayDOFV< nVarsFixed > & | u, | ||
| ArrayDOFV< nVarsFixed > & | uInc, | ||
| JacobianDiagBlock< nVarsFixed > & | JDiag, | ||
| ArrayDOFV< nVarsFixed > & | AuInc | ||
| ) |
Compute the matrix-vector product A * uInc for the implicit system.
Evaluates the action of the approximate Jacobian on an increment vector, used as the matvec operation inside GMRES.
| [in] | alphaDiag | Diagonal scaling factor. |
| [in] | t | Current simulation time. |
| [in] | u | Cell-centered conservative DOFs. |
| [in] | uInc | Increment vector to multiply. |
| [in] | JDiag | Pre-assembled diagonal Jacobian blocks. |
| [out] | AuInc | Result of A * uInc (overwritten). |
|
inlinestaticconstexpr |
Compute the maximum batch-multiplied column count for batched Eigen matrices.
Definition at line 88 of file EulerEvaluator.hpp.
| void DNDS::Euler::EulerEvaluator< model >::MeanValueCons2Prim | ( | ArrayDOFV< nVarsFixed > & | u, |
| ArrayDOFV< nVarsFixed > & | w | ||
| ) |
Convert cell-mean conservative variables to primitive variables.
Convert cell mean values from conservative to primitive variables.
| u | Conservative variable DOF array (input). |
| w | Primitive variable DOF array (output). |
Definition at line 1348 of file EulerEvaluator.hxx.
| void DNDS::Euler::EulerEvaluator< model >::MeanValuePrim2Cons | ( | ArrayDOFV< nVarsFixed > & | w, |
| ArrayDOFV< nVarsFixed > & | u | ||
| ) |
Convert cell-mean primitive variables to conservative variables.
Convert cell mean values from primitive to conservative variables.
| w | Primitive variable DOF array (input). |
| u | Conservative variable DOF array (output). |
Definition at line 1368 of file EulerEvaluator.hxx.
| void DNDS::Euler::EulerEvaluator< model >::MinSmoothDTau | ( | ArrayDOFV< 1 > & | dTau, |
| ArrayDOFV< 1 > & | dTauNew | ||
| ) |
Smooth the local time step across neighboring cells.
Smooth the local pseudo-time step by taking the minimum of the cell value and its neighbor average.
Prevents large time-step jumps between adjacent cells by capping each cell's dTau at the weighted average of its neighbors' values.
| dTau | Input local pseudo-time step per cell. |
| dTauNew | Smoothed pseudo-time step (output). |
Definition at line 2158 of file EulerEvaluator.hxx.
|
inline |
Compute effective molecular viscosity using the configured viscosity model.
Supports constant viscosity (muModel=0), Sutherland's law (muModel=1), and density-proportional viscosity (muModel=2).
| U | Conservative state vector. |
| T | Temperature. |
Definition at line 751 of file EulerEvaluator.hpp.
|
inline |
Write boundary profile data to CSV files (rank 0 only).
only 0 needs to write
Definition at line 1430 of file EulerEvaluator.hpp.
|
inline |
Enable or disable discarding of passive-scalar source terms.
Definition at line 166 of file EulerEvaluator.hpp.
| TU DNDS::Euler::EulerEvaluator< model >::source | ( | const TU & | UMeanXy, |
| const TDiffU & | DiffUxy, | ||
| const Geom::tPoint & | pPhy, | ||
| TJacobianU & | jacobian, | ||
| index | iCell, | ||
| index | ig, | ||
| int | Mode | ||
| ) |
Compute the source term at a cell quadrature point.
Evaluates body force, axisymmetric, rotating-frame, and RANS turbulence model source terms. Optionally computes the source Jacobian for implicit methods.
| [in] | UMeanXy | Conservative state at the quadrature point. |
| [in] | DiffUxy | Gradient of conservative variables. |
| [in] | pPhy | Physical coordinates of the quadrature point. |
| [out] | jacobian | Source Jacobian matrix (populated if Mode=1 or 2). |
| [in] | iCell | Cell index. |
| [in] | ig | Quadrature point index within the cell. |
| [in] | Mode | 0=source vector only, 1=diagonal Jacobian, 2=full Jacobian. |
| void DNDS::Euler::EulerEvaluator< model >::TimeAverageAddition | ( | ArrayDOFV< nVarsFixed > & | w, |
| ArrayDOFV< nVarsFixed > & | wAveraged, | ||
| real | dt, | ||
| real & | tCur | ||
| ) |
Accumulate time-averaged primitive variables for unsteady statistics.
Accumulate time-averaged solution field using running weighted average.
Updates wAveraged as a running weighted average: wAveraged = (tCur*wAveraged + dt*w) / (tCur+dt).
| w | Current primitive variable field. |
| wAveraged | Running time-averaged field (input/output). |
| dt | Time step weight for current sample. |
| tCur | Accumulated averaging time (input/output, incremented by dt). |
Definition at line 1335 of file EulerEvaluator.hxx.
|
inline |
Transform full conservative state (momentum + total energy) between frames (relative velocity formulation).
Definition at line 1255 of file EulerEvaluator.hpp.
|
inline |
Transform full conservative state for the absolute-velocity rotating frame formulation.
Definition at line 1266 of file EulerEvaluator.hpp.
|
inline |
Transform momentum in-place between inertial and rotating frame (velocity only).
Definition at line 1248 of file EulerEvaluator.hpp.
|
inline |
Transform a conservative state vector from cell frame to face frame for periodic BCs.
Applies the periodic rotation/translation to the momentum components when the face is a periodic boundary and the cell is on the donor side.
| [in,out] | u | Conservative state vector (modified in-place). |
| [in] | iFace | Face index. |
| [in] | iCell | Cell index. |
| [in] | if2c | Face-to-cell local index (0=back, 1=front, <0=auto-detect). |
Definition at line 879 of file EulerEvaluator.hpp.
|
inline |
Inverse of UFromCell2Face: transform from face frame back to cell frame.
Definition at line 896 of file EulerEvaluator.hpp.
|
inline |
Transform a state from a neighbor cell across a periodic face.
Definition at line 913 of file EulerEvaluator.hpp.
|
inline |
Update boundary anchor point recorders from current solution.
Definition at line 1277 of file EulerEvaluator.hpp.
| void DNDS::Euler::EulerEvaluator< model >::updateBCProfiles | ( | ArrayDOFV< nVarsFixed > & | u, |
| ArrayRECV< nVarsFixed > & | uRec | ||
| ) |
Update boundary 1-D profiles from the current solution.
Update boundary condition profiles from current solution for profile-anchored BCs.
For BC zones with anchorOpt==2, integrates the solution at boundary faces into a 1D radial profile recorder. Used for radial equilibrium pressure BCs in turbomachinery applications.
| u | Conservative variable DOF array. |
| uRec | Reconstruction coefficients. |
Definition at line 2195 of file EulerEvaluator.hxx.
| void DNDS::Euler::EulerEvaluator< model >::updateBCProfilesPressureRadialEq | ( | ) |
Update boundary profiles using radial-equilibrium pressure extrapolation.
Compute radial equilibrium pressure distribution from BC profiles.
Integrates the centrifugal pressure gradient (rho * vt^2 / r) along the radial direction in boundary profiles to obtain the radial equilibrium pressure increment for turbomachinery outlet BCs.
Definition at line 2295 of file EulerEvaluator.hxx.
| void DNDS::Euler::EulerEvaluator< model >::UpdateLUSGSBackward | ( | real | alphaDiag, |
| real | t, | ||
| ArrayDOFV< nVarsFixed > & | rhs, | ||
| ArrayDOFV< nVarsFixed > & | u, | ||
| ArrayDOFV< nVarsFixed > & | uInc, | ||
| JacobianDiagBlock< nVarsFixed > & | JDiag, | ||
| ArrayDOFV< nVarsFixed > & | uIncNew | ||
| ) |
Deprecated: use UpdateSGS with uIncIsZero=true instead. Kept for backward compatibility. Delegates to UpdateSGS.
| void DNDS::Euler::EulerEvaluator< model >::UpdateLUSGSForward | ( | real | alphaDiag, |
| real | t, | ||
| ArrayDOFV< nVarsFixed > & | rhs, | ||
| ArrayDOFV< nVarsFixed > & | u, | ||
| ArrayDOFV< nVarsFixed > & | uInc, | ||
| JacobianDiagBlock< nVarsFixed > & | JDiag, | ||
| ArrayDOFV< nVarsFixed > & | uIncNew | ||
| ) |
Deprecated: use UpdateSGS with uIncIsZero=true instead. Kept for backward compatibility. Delegates to UpdateSGS.
| void DNDS::Euler::EulerEvaluator< model >::UpdateSGS | ( | real | alphaDiag, |
| real | t, | ||
| ArrayDOFV< nVarsFixed > & | rhs, | ||
| ArrayDOFV< nVarsFixed > & | u, | ||
| ArrayDOFV< nVarsFixed > & | uInc, | ||
| ArrayDOFV< nVarsFixed > & | uIncNew, | ||
| JacobianDiagBlock< nVarsFixed > & | JDiag, | ||
| bool | forward, | ||
| bool | gsUpdate, | ||
| TU & | sumInc, | ||
| bool | uIncIsZero = false |
||
| ) |
Symmetric Gauss-Seidel update for the implicit linear system.
Symmetric Gauss-Seidel (SGS) sweep for implicit preconditioning.
| forward | Scan ascending (true) or descending (false). |
| gsUpdate | Use Gauss-Seidel update (read from uIncNew for already-processed cells). |
| uIncIsZero | If true, uInc is assumed zero for not-yet-processed cells, so their flux contribution is skipped (LUSGS optimisation). |
| sumInc | Output: accumulated absolute increment for convergence tracking. |
Performs a single forward or backward sweep. When uIncIsZero is true, skips off-diagonal terms involving the zeroed increment in the first forward pass to save computation. Accumulates the global L1 increment change into sumInc.
| alphaDiag | Diagonal scaling factor. |
| t | Current simulation time. |
| rhs | Right-hand side residual. |
| u | Conservative variable DOF array. |
| uInc | Current increment (input, must not alias uIncNew). |
| uIncNew | Updated increment (output). |
| JDiag | Diagonal Jacobian block. |
| forward | If true, sweep from cell 0 to N-1; otherwise N-1 to 0. |
| gsUpdate | If true, use Gauss-Seidel update (use latest values); otherwise Jacobi. |
| sumInc | Global L1 sum of increment changes (output, MPI-reduced). |
| uIncIsZero | If true, skip off-diagonal terms for zero-initialized increment. |
always inner here
Definition at line 550 of file EulerEvaluator.hxx.
| void DNDS::Euler::EulerEvaluator< model >::UpdateSGSWithRec | ( | real | alphaDiag, |
| real | t, | ||
| ArrayDOFV< nVarsFixed > & | rhs, | ||
| ArrayDOFV< nVarsFixed > & | u, | ||
| ArrayRECV< nVarsFixed > & | uRec, | ||
| ArrayDOFV< nVarsFixed > & | uInc, | ||
| ArrayRECV< nVarsFixed > & | uRecInc, | ||
| JacobianDiagBlock< nVarsFixed > & | JDiag, | ||
| bool | forward, | ||
| TU & | sumInc | ||
| ) |
SGS sweep coupled with reconstruction update.
SGS sweep with reconstruction-level off-diagonal terms.
Performs a single SGS sweep that simultaneously updates the conservative increment (uInc) and the reconstruction increment (uRecInc).
| [in] | alphaDiag | Diagonal scaling factor. |
| [in] | t | Current simulation time. |
| [in] | rhs | Right-hand side residual. |
| [in] | u | Cell-centered conservative DOFs. |
| [in] | uRec | Reconstruction coefficients. |
| [in] | uInc | Current increment for conservative variables. |
| [in] | uRecInc | Current reconstruction increment. |
| [in] | JDiag | Diagonal Jacobian blocks. |
| [in] | forward | Sweep direction: ascending (true) or descending (false). |
| [out] | sumInc | Accumulated absolute increment for convergence monitoring. |
Extends UpdateSGS by including reconstruction increment contributions (uRecInc) in the off-diagonal flux Jacobian product, enabling higher-order implicit coupling.
| alphaDiag | Diagonal scaling factor. |
| t | Current simulation time. |
| rhs | Right-hand side residual. |
| u | Conservative variable DOF array. |
| uRec | Reconstruction coefficients. |
| uInc | Current DOF increment. |
| uRecInc | Current reconstruction increment. |
| JDiag | Diagonal Jacobian block. |
| forward | Sweep direction: true for forward, false for backward. |
| sumInc | Global L1 sum of increment changes (output, MPI-reduced). |
TODO periodic here
always inner here
TODO periodic here
Definition at line 713 of file EulerEvaluator.hxx.
|
inline |
Compute viscous flux contribution from turbulence model variables.
Adds the turbulent diffusion flux for SA (nuTilde) or two-equation (k, omega/epsilon) variables, and the turbulent normal-stress correction to the momentum and energy fluxes.
| UMeanXYC | Cell-mean state in physical coordinates. |
| DiffUxyPrimC | Gradient of primitive variables in physical coordinates. |
| muRef | Reference viscosity scaling. |
| mufPhy | Physical molecular viscosity. |
| muTur | Turbulent eddy viscosity. |
| uNormC | Face outward unit normal. |
| iFace | Face index (for wall distance lookup). |
| VisFlux | Viscous flux vector (accumulated in-place). |
SA's normal stress
Definition at line 835 of file EulerEvaluator.hpp.
| std::map<Geom::t_index, AnchorPointRecorder<nVarsFixed> > DNDS::Euler::EulerEvaluator< model >::anchorRecorders |
Per-BC anchor point recorders for profile extraction.
Definition at line 201 of file EulerEvaluator.hpp.
| std::map<Geom::t_index, std::ofstream> DNDS::Euler::EulerEvaluator< model >::bndIntegrationLogs |
Per-BC log file streams for integration output.
Definition at line 204 of file EulerEvaluator.hpp.
| std::map<Geom::t_index, IntegrationRecorder> DNDS::Euler::EulerEvaluator< model >::bndIntegrations |
Per-BC boundary flux/force integration accumulators.
Definition at line 203 of file EulerEvaluator.hpp.
| std::set<Geom::t_index> DNDS::Euler::EulerEvaluator< model >::cLDriverBndIDs |
Boundary IDs driven by the CL (lift-coefficient) driver.
Definition at line 206 of file EulerEvaluator.hpp.
| ArrayDOFV<1> DNDS::Euler::EulerEvaluator< model >::deltaLambdaCell |
Per-cell accumulated spectral radius difference.
Definition at line 191 of file EulerEvaluator.hpp.
| std::vector<real> DNDS::Euler::EulerEvaluator< model >::deltaLambdaFace |
Per-face spectral radius difference for implicit diagonal.
Definition at line 190 of file EulerEvaluator.hpp.
|
static |
Spatial dimension (2 or 3).
Definition at line 82 of file EulerEvaluator.hpp.
|
static |
Skip recomputation of per-face eigenvalues lambda0/123/4.
Definition at line 340 of file EulerEvaluator.hpp.
|
static |
No flags for EvaluateDt.
Definition at line 339 of file EulerEvaluator.hpp.
| std::vector<Eigen::Vector<real, Eigen::Dynamic> > DNDS::Euler::EulerEvaluator< model >::dWall |
Per-cell wall distance (one value per cell node/quadrature point).
Definition at line 197 of file EulerEvaluator.hpp.
| std::vector<real> DNDS::Euler::EulerEvaluator< model >::dWallFace |
Per-face wall distance (interpolated from cell values).
Definition at line 198 of file EulerEvaluator.hpp.
|
static |
Default alpha evaluation mode.
Definition at line 691 of file EulerEvaluator.hpp.
|
static |
Always take min(alpha, prev).
Definition at line 693 of file EulerEvaluator.hpp.
|
static |
Take min(alpha, prev) only if prev != 1.
Definition at line 692 of file EulerEvaluator.hpp.
|
static |
Compress reconstruction toward cell mean to enforce positivity.
Definition at line 660 of file EulerEvaluator.hpp.
|
static |
Default: evaluate beta without compression.
Definition at line 659 of file EulerEvaluator.hpp.
| std::vector<TU> DNDS::Euler::EulerEvaluator< model >::fluxBnd |
Per-boundary-face flux values.
Definition at line 214 of file EulerEvaluator.hpp.
| std::vector<TVec> DNDS::Euler::EulerEvaluator< model >::fluxBndForceT |
Per-boundary-face tangential force.
Definition at line 215 of file EulerEvaluator.hpp.
| Eigen::Vector<real, -1> DNDS::Euler::EulerEvaluator< model >::fluxWallSum |
Accumulated wall flux integral (force on wall boundaries).
Definition at line 213 of file EulerEvaluator.hpp.
|
static |
Geometric dimension of the mesh.
Definition at line 83 of file EulerEvaluator.hpp.
| std::vector<TDiffU> DNDS::Euler::EulerEvaluator< model >::gradUFix |
Green-Gauss gradient correction buffer for source term stabilization.
Definition at line 194 of file EulerEvaluator.hpp.
|
static |
Index of the energy equation (= dim+1).
Definition at line 84 of file EulerEvaluator.hpp.
| int DNDS::Euler::EulerEvaluator< model >::kAv = 0 |
Artificial viscosity polynomial order (maxOrder + 1).
Definition at line 180 of file EulerEvaluator.hpp.
| std::vector<real> DNDS::Euler::EulerEvaluator< model >::lambdaFace |
Per-face spectral radius (inviscid + viscous combined).
Definition at line 184 of file EulerEvaluator.hpp.
| std::vector<real> DNDS::Euler::EulerEvaluator< model >::lambdaFace0 |
Per-face eigenvalue |u·n| (contact wave).
Definition at line 187 of file EulerEvaluator.hpp.
| std::vector<real> DNDS::Euler::EulerEvaluator< model >::lambdaFace123 |
Per-face eigenvalue |u·n + a| (acoustic wave).
Definition at line 188 of file EulerEvaluator.hpp.
| std::vector<real> DNDS::Euler::EulerEvaluator< model >::lambdaFace4 |
Per-face eigenvalue |u·n - a| (acoustic wave).
Definition at line 189 of file EulerEvaluator.hpp.
| std::vector<real> DNDS::Euler::EulerEvaluator< model >::lambdaFaceC |
Per-face convective spectral radius.
Definition at line 185 of file EulerEvaluator.hpp.
| std::vector<real> DNDS::Euler::EulerEvaluator< model >::lambdaFaceVis |
Per-face viscous spectral radius.
Definition at line 186 of file EulerEvaluator.hpp.
|
static |
Disable shock-detecting component of the limiter.
Definition at line 642 of file EulerEvaluator.hpp.
|
static |
Default limiter flags.
Definition at line 641 of file EulerEvaluator.hpp.
|
static |
Maximum batch size for vectorized quadrature-point evaluation.
Definition at line 86 of file EulerEvaluator.hpp.
| ssp<Geom::UnstructuredMesh> DNDS::Euler::EulerEvaluator< model >::mesh |
Shared pointer to the unstructured mesh.
Definition at line 177 of file EulerEvaluator.hpp.
| index DNDS::Euler::EulerEvaluator< model >::nFaceReducedOrder = 0 |
Count of faces where reconstruction order was reduced.
Definition at line 216 of file EulerEvaluator.hpp.
|
static |
Number of conserved variables (compile-time fixed).
Definition at line 81 of file EulerEvaluator.hpp.
| ssp<BoundaryHandler<model> > DNDS::Euler::EulerEvaluator< model >::pBCHandler |
gDim -> 3 for intellisense //!tmptmp ///< Variational reconstruction object.
Boundary condition handler.
Definition at line 179 of file EulerEvaluator.hpp.
| std::unique_ptr<CLDriver> DNDS::Euler::EulerEvaluator< model >::pCLDriver |
Lift-coefficient driver for AoA adaptation (null if unused).
Definition at line 207 of file EulerEvaluator.hpp.
| std::map<Geom::t_index, OneDimProfile<nVarsFixed> > DNDS::Euler::EulerEvaluator< model >::profileRecorders |
Per-BC 1-D profile recorders.
Definition at line 202 of file EulerEvaluator.hpp.
|
static |
Use direct 2nd-order reconstruction.
Definition at line 374 of file EulerEvaluator.hpp.
|
static |
2nd-order rec with 1st-order convection.
Definition at line 375 of file EulerEvaluator.hpp.
|
static |
uGradBufNoLim is already computed.
Definition at line 377 of file EulerEvaluator.hpp.
|
static |
Apply limiter when using direct 2nd rec.
Definition at line 376 of file EulerEvaluator.hpp.
|
static |
Skip recording per-boundary flux.
Definition at line 373 of file EulerEvaluator.hpp.
|
static |
Skip boundary integration accumulation.
Definition at line 372 of file EulerEvaluator.hpp.
|
static |
Skip viscous flux contribution.
Definition at line 371 of file EulerEvaluator.hpp.
|
static |
Default: full RHS evaluation.
Definition at line 370 of file EulerEvaluator.hpp.
|
static |
Recover incremental face scaling.
Definition at line 378 of file EulerEvaluator.hpp.
| EulerEvaluatorSettings<model> DNDS::Euler::EulerEvaluator< model >::settings |
Physics and numerics settings for this evaluator.
Definition at line 220 of file EulerEvaluator.hpp.
| ssp<Direct::SerialSymLUStructure> DNDS::Euler::EulerEvaluator< model >::symLU |
Symmetric LU structure for direct preconditioner.
Definition at line 218 of file EulerEvaluator.hpp.
| ArrayGRADV<nVarsFixed, gDim> DNDS::Euler::EulerEvaluator< model >::uGradBuf |
Definition at line 211 of file EulerEvaluator.hpp.
| ArrayGRADV<nVarsFixed, gDim> DNDS::Euler::EulerEvaluator< model >::uGradBufNoLim |
Gradient buffers: limited and unlimited.
Definition at line 211 of file EulerEvaluator.hpp.
| ssp<CFV::VariationalReconstruction<gDim> > DNDS::Euler::EulerEvaluator< model >::vfv |
Definition at line 178 of file EulerEvaluator.hpp.