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

Namespace for the EulerP alternative evaluator module with GPU support. More...

Namespaces

namespace  detail
 Internal implementation detail namespace for EulerP Evaluator utilities.
 

Classes

class  BC
 Host-side boundary condition object managing parameter values and device transfer. More...
 
class  BC_DeviceView
 Device-callable view of a single boundary condition. More...
 
struct  BCFunc_Impl
 Primary template for boundary condition function implementations. More...
 
struct  BCFunc_Impl< B, BCType::Far >
 Farfield BC implementation (TODO: not yet implemented). More...
 
struct  BCFunc_Impl< B, BCType::In >
 General inflow BC implementation (TODO: not yet implemented). More...
 
struct  BCFunc_Impl< B, BCType::InPsTs >
 Inflow with stagnation pressure/temperature BC implementation (TODO: not yet implemented). More...
 
struct  BCFunc_Impl< B, BCType::Out >
 Outflow BC implementation (TODO: not yet implemented). More...
 
struct  BCFunc_Impl< B, BCType::OutP >
 Outflow with back-pressure BC implementation (TODO: not yet implemented). More...
 
struct  BCFunc_Impl< B, BCType::Special >
 Special-purpose BC implementation for benchmark cases (TODO: not yet implemented). More...
 
struct  BCFunc_Impl< B, BCType::Sym >
 Symmetry plane BC implementation (TODO: not yet implemented). More...
 
struct  BCFunc_Impl< B, BCType::Wall >
 No-slip viscous wall BC implementation (TODO: not yet implemented). More...
 
struct  BCFunc_Impl< B, BCType::WallInvis >
 Inviscid (slip) wall BC implementation (TODO: not yet implemented). More...
 
struct  BCFunc_Impl< B, BCType::WallIsothermal >
 Isothermal viscous wall BC implementation (TODO: not yet implemented). More...
 
class  BCHandler
 Host-side boundary condition handler managing all BC objects for a simulation. More...
 
class  BCHandlerDeviceView
 Device-callable view of the BC handler providing BC lookup by zone ID. More...
 
struct  BCInput
 Simple struct for JSON-deserialized boundary condition input specification. More...
 
class  Evaluator
 Main EulerP evaluator orchestrating CFD kernel dispatch for compressible N-S solvers. More...
 
struct  Evaluator_impl
 Backend-specific implementation of EulerP Evaluator kernels. More...
 
class  EvaluatorArgBase
 CRTP base class for packed kernel argument structs used by the Evaluator. More...
 
class  EvaluatorConfig
 JSON-configurable settings for the Evaluator. More...
 
class  EvaluatorDeviceView
 Device-callable view of the Evaluator, combining finite volume, BC handler, and physics views. More...
 
struct  Physics
 Host-side physics object managing gas parameters and device-transferable reference values. More...
 
struct  PhysicsDeviceView
 Device-callable view of physics parameters providing thermodynamic operations. More...
 
struct  PhysicsParams
 POD struct holding gas physical properties for the ideal gas model. More...
 

Typedefs

using TU = Eigen::Vector< real, nVarsFlow >
 Fixed-size 5-component conservative state vector (rho, rhoU, rhoV, rhoW, E).
 
using TDiffU = Eigen::Matrix< real, 3, nVarsFlow >
 Fixed-size 3x5 spatial gradient of the conservative state (one row per spatial dimension).
 
using TUFull = Eigen::Vector< real, Eigen::Dynamic >
 Dynamic-size state vector for extended variables (flow + turbulence model variables).
 
using TDiffUFull = Eigen::Matrix< real, 3, Eigen::Dynamic >
 Dynamic-size 3xN spatial gradient for extended variables.
 
using TUFullMap = Eigen::Map< TUFull >
 Eigen::Map wrapper for non-owning access to a dynamic-size state vector.
 
using TDiffUFullMap = Eigen::Map< TDiffUFull >
 Eigen::Map wrapper for non-owning access to a dynamic-size gradient matrix.
 
using TUDof = ArrayDof< nVarsFlow, 1 >
 Distributed array of 5-component flow state vectors (one per DOF).
 
using TUGrad = ArrayDof< 3, nVarsFlow >
 Distributed array of 3x5 gradient matrices (one per DOF).
 
using TUScalar = ArrayDof< 1, 1 >
 Distributed array of scalar values (e.g., pressure, temperature).
 
using TUScalarGrad = ArrayDof< 3, 1 >
 Distributed array of 3-component scalar gradients.
 
using TUScalar2 = ArrayDof< 2, 1 >
 Distributed array of 2-component scalar values (e.g., paired quantities).
 
using TUVec = ArrayDof< 3, 1 >
 Distributed array of 3-component vectors (e.g., velocity, coordinates).
 
using TUVecGrad = ArrayDof< 3, 3 >
 Distributed array of 3x3 vector gradient tensors.
 
using TFiniteVolume = CFV::FiniteVolume
 

Enumerations

enum class  BCType : uint8_t {
  Unknown = 0 , Far , Wall , WallInvis ,
  WallIsothermal , Out , OutP , In ,
  InPsTs , Sym , Special
}
 Enumeration of boundary condition types for the EulerP module. More...
 

Functions

template<class TU >
DNDS_DEVICE_CALLABLE DNDS_FORCEINLINE auto U123 (TU &&v)
 Extracts the momentum components (indices 1,2,3) from a state vector as a 3x1 block.
 
template<class TU >
DNDS_DEVICE_CALLABLE DNDS_FORCEINLINE auto U012 (TU &&v)
 Extracts the first 3 components (indices 0,1,2) from a vector as a 3x1 block.
 
DNDS_DEVICE_CALLABLE void RoeEigenValueFixer (real aL, real aR, real vnL, real vnR, real dLambda, real fixScale, real &lam0, real &lam123, real &lam4)
 Applies entropy fix to Roe eigenvalues to prevent expansion shocks.
 
template<DeviceBackend B, class TUL , class TUR , class TULPrim , class TURPrim >
DNDS_DEVICE_CALLABLE void RoeAverageNS (TUL &&UL, TUR &&UR, TULPrim &&ULPrim, TURPrim &&URPrim, int nVars, real pL, real pR, PhysicsDeviceView< B > &phy, Geom::tPoint &veloRoe, real &vsqrRoe, real &HRoe, real &rhoRoe, real &aSqrRoe)
 Computes Roe-averaged quantities from left and right states.
 
template<DeviceBackend B>
DNDS_DEVICE_CALLABLE void RoeFluxFlow (const TU &UL, const TU &UR, real pL, real pR, const Geom::tPoint &veloRoe, real vsqrRoe, real vgn, const Geom::tPoint &n, real asqrRoe, real aRoe, real HRoe, PhysicsDeviceView< B > &phy, real lam0, real lam123, real lam4, TU &F)
 Computes the complete Roe numerical flux for the 5-equation flow system.
 
 DNDS_DEFINE_ENUM_JSON (BCType, { {BCType::Unknown, nullptr}, {BCType::Far, "Far"}, {BCType::Wall, "Wall"}, {BCType::WallInvis, "WallInvis"}, {BCType::WallIsothermal, "WallIsothermal"}, {BCType::Out, "Out"}, {BCType::OutP, "OutP"}, {BCType::In, "In"}, {BCType::InPsTs, "InPsTs"}, {BCType::Sym, "Sym"}, {BCType::Special, "Special"}, }) template< DeviceBackend B > using BCStorageVecDeviceView
 Type alias for the device-resident BC parameter storage vector.
 
void pybind11_BCType_define (py::module_ &m)
 Registers the BCType enum as a Python enum class.
 
void pybind11_BCInput_define (py::module_ &m)
 Registers pybind11 bindings for the BCInput struct.
 
void pybind11_BC_define (py::module_ &m)
 Registers pybind11 bindings for the BC class.
 
void pybind11_BCHandler_define (py::module_ &m)
 Registers pybind11 bindings for the BCHandler class.
 
void pybind11_BC_bind (py::module_ &m)
 Top-level binding function for all EulerP boundary condition Python API.
 
void pybind11_Evaluator_define (py::module_ &m)
 Defines pybind11 bindings for the Evaluator class and its kernel argument structs.
 
void pybind11_Evaluator_bind (py::module_ &m)
 Top-level binding function for the EulerP Evaluator Python API.
 
template<DeviceBackend B = DeviceBackend::Host>
DNDS_DEVICE void RecGradient_GGRec_Kernel_BndVal (EvaluatorDeviceView< B > &self_view, typename Evaluator_impl< B >::RecGradient_Arg::Portable &arg, index iBnd, index iBndEnd, int nVars, int nVarsScalar)
 Generates boundary ghost values for Green-Gauss gradient reconstruction.
 
template DNDS_DEVICE void RecGradient_GGRec_Kernel_BndVal< DeviceBackend::Host > (EvaluatorDeviceView< DeviceBackend::Host > &self_view, typename Evaluator_impl< DeviceBackend::Host >::RecGradient_Arg::Portable &arg, index iBnd, index iBndEnd, int nVars, int nVarsScalar)
 
template<DeviceBackend B = DeviceBackend::Host>
DNDS_DEVICE void RecGradient_GGRec_Kernel_GG (EvaluatorDeviceView< B > &self_view, typename Evaluator_impl< B >::RecGradient_Arg::Portable &arg, index iCell, index iCellEnd, int nVars, int nVarsScalar)
 Green-Gauss gradient reconstruction kernel (per-cell).
 
template DNDS_DEVICE void RecGradient_GGRec_Kernel_GG< DeviceBackend::Host > (EvaluatorDeviceView< DeviceBackend::Host > &self_view, typename Evaluator_impl< DeviceBackend::Host >::RecGradient_Arg::Portable &arg, index iCell, index iCellEnd, int nVars, int nVarsScalar)
 
template<DeviceBackend B = DeviceBackend::Host>
DNDS_DEVICE void RecGradient_BarthLimiter_Kernel_FlowPart (EvaluatorDeviceView< B > &self_view, typename Evaluator_impl< B >::RecGradient_Arg::Portable &arg, index iCell, index iCellEnd, int nVars, int nVarsScalar)
 Barth-Jespersen gradient limiter kernel for flow variables (per-cell).
 
template DNDS_DEVICE void RecGradient_BarthLimiter_Kernel_FlowPart< DeviceBackend::Host > (EvaluatorDeviceView< DeviceBackend::Host > &self_view, typename Evaluator_impl< DeviceBackend::Host >::RecGradient_Arg::Portable &arg, index iCell, index iCellEnd, int nVars, int nVarsScalar)
 
template<DeviceBackend B = DeviceBackend::Host>
DNDS_DEVICE void RecGradient_BarthLimiter_Kernel_ScalarPart (EvaluatorDeviceView< B > &self_view, typename Evaluator_impl< B >::RecGradient_Arg::Portable &arg, index iCell, index iCellEnd, int nVars, int nVarsScalar)
 Barth-Jespersen gradient limiter kernel for transported scalar variables (per-cell).
 
template DNDS_DEVICE void RecGradient_BarthLimiter_Kernel_ScalarPart< DeviceBackend::Host > (EvaluatorDeviceView< DeviceBackend::Host > &self_view, typename Evaluator_impl< DeviceBackend::Host >::RecGradient_Arg::Portable &arg, index iCell, index iCellEnd, int nVars, int nVarsScalar)
 
template<DeviceBackend B = DeviceBackend::Host>
DNDS_DEVICE void Cons2PrimMu_Kernel (EvaluatorDeviceView< B > &self_view, typename Evaluator_impl< B >::Cons2PrimMu_Arg::Portable &arg, index iPt, index iPtEnd, int nVars, int nVarsScalar)
 Conservative-to-primitive conversion with gradient transformation and viscosity (per-cell).
 
template DNDS_DEVICE void Cons2PrimMu_Kernel< DeviceBackend::Host > (EvaluatorDeviceView< DeviceBackend::Host > &self_view, typename Evaluator_impl< DeviceBackend::Host >::Cons2PrimMu_Arg::Portable &arg, index iPt, index iPtEnd, int nVars, int nVarsScalar)
 
template<DeviceBackend B = DeviceBackend::Host>
DNDS_DEVICE void Cons2Prim_Kernel (EvaluatorDeviceView< B > &self_view, typename Evaluator_impl< B >::Cons2Prim_Arg::Portable &arg, index iPt, index iPtEnd, int nVars, int nVarsScalar)
 Conservative-to-primitive conversion without gradient transformation or viscosity (per-cell).
 
template DNDS_DEVICE void Cons2Prim_Kernel< DeviceBackend::Host > (EvaluatorDeviceView< DeviceBackend::Host > &self_view, typename Evaluator_impl< DeviceBackend::Host >::Cons2Prim_Arg::Portable &arg, index iPt, index iPtEnd, int nVars, int nVarsScalar)
 
template<DeviceBackend B = DeviceBackend::Host>
DNDS_DEVICE void EstEigenDt_GetFaceLam_Kernel (EvaluatorDeviceView< B > &self_view, typename Evaluator_impl< B >::EstEigenDt_Arg::Portable &arg, index iFace, index iFaceEnd, int nVars, int nVarsScalar)
 Per-face eigenvalue estimation kernel for time-step computation.
 
template DNDS_DEVICE void EstEigenDt_GetFaceLam_Kernel< DeviceBackend::Host > (EvaluatorDeviceView< DeviceBackend::Host > &self_view, typename Evaluator_impl< DeviceBackend::Host >::EstEigenDt_Arg::Portable &arg, index iFace, index iFaceEnd, int nVars, int nVarsScalar)
 
template<DeviceBackend B = DeviceBackend::Host>
DNDS_DEVICE void EstEigenDt_FaceLam2CellDt_Kernel (EvaluatorDeviceView< B > &self_view, typename Evaluator_impl< B >::EstEigenDt_Arg::Portable &arg, index iCell, index iCellEnd, int nVars, int nVarsScalar)
 Per-cell kernel converting face eigenvalues to a local CFL time step.
 
template DNDS_DEVICE void EstEigenDt_FaceLam2CellDt_Kernel< DeviceBackend::Host > (EvaluatorDeviceView< DeviceBackend::Host > &self_view, typename Evaluator_impl< DeviceBackend::Host >::EstEigenDt_Arg::Portable &arg, index iCell, index iCellEnd, int nVars, int nVarsScalar)
 
template<DeviceBackend B = DeviceBackend::Host>
DNDS_DEVICE void RecFace2nd_Kernel (EvaluatorDeviceView< B > &self_view, typename Evaluator_impl< B >::RecFace2nd_Arg::Portable &arg, index iFace, index iFaceEnd, int nVars, int nVarsScalar)
 2nd-order face value reconstruction kernel (per-face).
 
template DNDS_DEVICE void RecFace2nd_Kernel< DeviceBackend::Host > (EvaluatorDeviceView< DeviceBackend::Host > &self_view, typename Evaluator_impl< DeviceBackend::Host >::RecFace2nd_Arg::Portable &arg, index iFace, index iFaceEnd, int nVars, int nVarsScalar)
 
template<DeviceBackend B = DeviceBackend::Host>
DNDS_DEVICE void Flux2nd_Kernel_FluxFace (EvaluatorDeviceView< B > &self_view, typename Evaluator_impl< B >::Flux2nd_Arg::Portable &arg, index iFace, index iFaceEnd, int nVars, int nVarsScalar)
 2nd-order Roe inviscid flux computation kernel (per-face).
 
template DNDS_DEVICE void Flux2nd_Kernel_FluxFace< DeviceBackend::Host > (EvaluatorDeviceView< DeviceBackend::Host > &self_view, typename Evaluator_impl< DeviceBackend::Host >::Flux2nd_Arg::Portable &arg, index iFace, index iFaceEnd, int nVars, int nVarsScalar)
 
template<DeviceBackend B = DeviceBackend::Host>
DNDS_DEVICE void Flux2nd_Kernel_Face2Cell (EvaluatorDeviceView< B > &self_view, typename Evaluator_impl< B >::Flux2nd_Arg::Portable &arg, index iCell, index iCellEnd, int nVars, int nVarsScalar)
 Scatters face fluxes to cell RHS residual (per-cell).
 
template DNDS_DEVICE void Flux2nd_Kernel_Face2Cell< DeviceBackend::Host > (EvaluatorDeviceView< DeviceBackend::Host > &self_view, typename Evaluator_impl< DeviceBackend::Host >::Flux2nd_Arg::Portable &arg, index iCell, index iCellEnd, int nVars, int nVarsScalar)
 
template<typename TU , typename TF , class TVecN >
DNDS_DEVICE_CALLABLE void GasInviscidFlux_XY (TU &&U, int nVars, real vn, real vgn, TVecN &&n, real p, TF &F)
 Computes the inviscid (Euler) flux projected onto a face normal direction.
 
void pybind11_Physics_bind (py::module_ &m)
 Registers pybind11 bindings for the Physics class.
 

Detailed Description

Namespace for the EulerP alternative evaluator module with GPU support.

Provides device-callable (Host + CUDA) compressible Navier-Stokes solvers using scalar loops instead of Eigen compile-time templates, enabling GPU offloading.

Typedef Documentation

◆ TDiffU

using DNDS::EulerP::TDiffU = typedef Eigen::Matrix<real, 3, nVarsFlow>

Fixed-size 3x5 spatial gradient of the conservative state (one row per spatial dimension).

Definition at line 34 of file EulerP.hpp.

◆ TDiffUFull

using DNDS::EulerP::TDiffUFull = typedef Eigen::Matrix<real, 3, Eigen::Dynamic>

Dynamic-size 3xN spatial gradient for extended variables.

Definition at line 37 of file EulerP.hpp.

◆ TDiffUFullMap

using DNDS::EulerP::TDiffUFullMap = typedef Eigen::Map<TDiffUFull>

Eigen::Map wrapper for non-owning access to a dynamic-size gradient matrix.

Definition at line 40 of file EulerP.hpp.

◆ TFiniteVolume

Definition at line 27 of file EulerP_Evaluator.hpp.

◆ TU

using DNDS::EulerP::TU = typedef Eigen::Vector<real, nVarsFlow>

Fixed-size 5-component conservative state vector (rho, rhoU, rhoV, rhoW, E).

Definition at line 33 of file EulerP.hpp.

◆ TUDof

using DNDS::EulerP::TUDof = typedef ArrayDof<nVarsFlow, 1>

Distributed array of 5-component flow state vectors (one per DOF).

Definition at line 42 of file EulerP.hpp.

◆ TUFull

using DNDS::EulerP::TUFull = typedef Eigen::Vector<real, Eigen::Dynamic>

Dynamic-size state vector for extended variables (flow + turbulence model variables).

Definition at line 36 of file EulerP.hpp.

◆ TUFullMap

using DNDS::EulerP::TUFullMap = typedef Eigen::Map<TUFull>

Eigen::Map wrapper for non-owning access to a dynamic-size state vector.

Definition at line 39 of file EulerP.hpp.

◆ TUGrad

using DNDS::EulerP::TUGrad = typedef ArrayDof<3, nVarsFlow>

Distributed array of 3x5 gradient matrices (one per DOF).

Definition at line 43 of file EulerP.hpp.

◆ TUScalar

using DNDS::EulerP::TUScalar = typedef ArrayDof<1, 1>

Distributed array of scalar values (e.g., pressure, temperature).

Definition at line 45 of file EulerP.hpp.

◆ TUScalar2

using DNDS::EulerP::TUScalar2 = typedef ArrayDof<2, 1>

Distributed array of 2-component scalar values (e.g., paired quantities).

Definition at line 48 of file EulerP.hpp.

◆ TUScalarGrad

using DNDS::EulerP::TUScalarGrad = typedef ArrayDof<3, 1>

Distributed array of 3-component scalar gradients.

Definition at line 46 of file EulerP.hpp.

◆ TUVec

using DNDS::EulerP::TUVec = typedef ArrayDof<3, 1>

Distributed array of 3-component vectors (e.g., velocity, coordinates).

Definition at line 50 of file EulerP.hpp.

◆ TUVecGrad

using DNDS::EulerP::TUVecGrad = typedef ArrayDof<3, 3>

Distributed array of 3x3 vector gradient tensors.

Definition at line 51 of file EulerP.hpp.

Enumeration Type Documentation

◆ BCType

enum class DNDS::EulerP::BCType : uint8_t
strong

Enumeration of boundary condition types for the EulerP module.

Mirrors the BC types in the Euler module. JSON-serializable via DNDS_DEFINE_ENUM_JSON.

Enumerator
Unknown 

Uninitialized or unrecognized BC type.

Far 

Farfield BC (characteristic-based).

Wall 

No-slip viscous wall (adiabatic).

WallInvis 

Inviscid (slip) wall BC.

WallIsothermal 

No-slip viscous wall with fixed temperature.

Out 

Supersonic/subsonic outflow BC.

OutP 

Outflow with specified back-pressure.

In 

Supersonic/subsonic inflow BC.

InPsTs 

Inflow with specified stagnation pressure and temperature.

Sym 

Symmetry plane BC.

Special 

Special-purpose BC for benchmark cases (e.g., DMR, Riemann).

Definition at line 31 of file EulerP_BC.hpp.

Function Documentation

◆ Cons2Prim_Kernel()

template<DeviceBackend B = DeviceBackend::Host>
DNDS_DEVICE void DNDS::EulerP::Cons2Prim_Kernel ( EvaluatorDeviceView< B > &  self_view,
typename Evaluator_impl< B >::Cons2Prim_Arg::Portable &  arg,
index  iPt,
index  iPtEnd,
int  nVars,
int  nVarsScalar 
)

Conservative-to-primitive conversion without gradient transformation or viscosity (per-cell).

Simplified version of Cons2PrimMu_Kernel that only converts conservative to primitive variables and computes thermodynamic scalars (T, p, a, gamma). No gradient transformation or viscosity computation is performed.

Template Parameters
BDeviceBackend (Host or CUDA), defaulting to Host.
Parameters
self_viewDevice view providing mesh, BC handler, and physics access.
argPortable struct with input/output device views.
iPtPoint (cell) index to process.
iPtEndUpper bound of point index range.
nVarsTotal number of variables (flow + scalar).
nVarsScalarNumber of additional transported scalar variables.

Definition at line 633 of file EulerP_Evaluator_impl_common.hxx.

Here is the caller graph for this function:

◆ Cons2Prim_Kernel< DeviceBackend::Host >()

template DNDS_DEVICE void DNDS::EulerP::Cons2Prim_Kernel< DeviceBackend::Host > ( EvaluatorDeviceView< DeviceBackend::Host > &  self_view,
typename Evaluator_impl< DeviceBackend::Host >::Cons2Prim_Arg::Portable &  arg,
index  iPt,
index  iPtEnd,
int  nVars,
int  nVarsScalar 
)

◆ Cons2PrimMu_Kernel()

template<DeviceBackend B = DeviceBackend::Host>
DNDS_DEVICE void DNDS::EulerP::Cons2PrimMu_Kernel ( EvaluatorDeviceView< B > &  self_view,
typename Evaluator_impl< B >::Cons2PrimMu_Arg::Portable &  arg,
index  iPt,
index  iPtEnd,
int  nVars,
int  nVarsScalar 
)

Conservative-to-primitive conversion with gradient transformation and viscosity (per-cell).

For each cell, converts conservative variables (rho, rho*u, rho*E, ...) to primitive variables (rho, u, v, w, p_or_E_internal, ...) and transforms the conservative gradients into primitive gradients using the Jacobian of the transformation. Also computes thermodynamic quantities (p, T, a, gamma) and total viscosity (mu).

Template Parameters
BDeviceBackend (Host or CUDA), defaulting to Host.
Parameters
self_viewDevice view providing mesh, BC handler, and physics access.
argPortable struct with all input/output device views.
iPtPoint (cell) index to process.
iPtEndUpper bound of point index range.
nVarsTotal number of variables (flow + scalar).
nVarsScalarNumber of additional transported scalar variables.

Definition at line 493 of file EulerP_Evaluator_impl_common.hxx.

Here is the caller graph for this function:

◆ Cons2PrimMu_Kernel< DeviceBackend::Host >()

template DNDS_DEVICE void DNDS::EulerP::Cons2PrimMu_Kernel< DeviceBackend::Host > ( EvaluatorDeviceView< DeviceBackend::Host > &  self_view,
typename Evaluator_impl< DeviceBackend::Host >::Cons2PrimMu_Arg::Portable &  arg,
index  iPt,
index  iPtEnd,
int  nVars,
int  nVarsScalar 
)

◆ DNDS_DEFINE_ENUM_JSON()

DNDS::EulerP::DNDS_DEFINE_ENUM_JSON ( BCType  ,
{ {BCType::Unknown, nullptr}, {BCType::Far, "Far"}, {BCType::Wall, "Wall"}, {BCType::WallInvis, "WallInvis"}, {BCType::WallIsothermal, "WallIsothermal"}, {BCType::Out, "Out"}, {BCType::OutP, "OutP"}, {BCType::In, "In"}, {BCType::InPsTs, "InPsTs"}, {BCType::Sym, "Sym"}, {BCType::Special, "Special"}, }   
)

Type alias for the device-resident BC parameter storage vector.

◆ EstEigenDt_FaceLam2CellDt_Kernel()

template<DeviceBackend B = DeviceBackend::Host>
DNDS_DEVICE void DNDS::EulerP::EstEigenDt_FaceLam2CellDt_Kernel ( EvaluatorDeviceView< B > &  self_view,
typename Evaluator_impl< B >::EstEigenDt_Arg::Portable &  arg,
index  iCell,
index  iCellEnd,
int  nVars,
int  nVarsScalar 
)

Per-cell kernel converting face eigenvalues to a local CFL time step.

For each cell, sums the maximum face eigenvalue (convective + viscous) weighted by face area across all cell faces. The local time step is computed as: dt = V * smoothScaleRatio / sum(lambda_face * A_face). Also computes the per-cell maximum deltaLamFace for H-correction.

Template Parameters
BDeviceBackend (Host or CUDA), defaulting to Host.
Parameters
self_viewDevice view providing mesh and FV geometry access.
argPortable struct with device views of eigenvalue and dt arrays.
iCellCell index to process.
iCellEndUpper bound of cell index range.
nVarsTotal number of variables (flow + scalar).
nVarsScalarNumber of additional transported scalar variables.

Definition at line 821 of file EulerP_Evaluator_impl_common.hxx.

Here is the caller graph for this function:

◆ EstEigenDt_FaceLam2CellDt_Kernel< DeviceBackend::Host >()

template DNDS_DEVICE void DNDS::EulerP::EstEigenDt_FaceLam2CellDt_Kernel< DeviceBackend::Host > ( EvaluatorDeviceView< DeviceBackend::Host > &  self_view,
typename Evaluator_impl< DeviceBackend::Host >::EstEigenDt_Arg::Portable &  arg,
index  iCell,
index  iCellEnd,
int  nVars,
int  nVarsScalar 
)

◆ EstEigenDt_GetFaceLam_Kernel()

template<DeviceBackend B = DeviceBackend::Host>
DNDS_DEVICE void DNDS::EulerP::EstEigenDt_GetFaceLam_Kernel ( EvaluatorDeviceView< B > &  self_view,
typename Evaluator_impl< B >::EstEigenDt_Arg::Portable &  arg,
index  iFace,
index  iFaceEnd,
int  nVars,
int  nVarsScalar 
)

Per-face eigenvalue estimation kernel for time-step computation.

For each face, estimates the convective eigenvalues (lam-, lam0, lam+) as the maximum of left/right normal velocity +/- speed of sound. Also estimates the viscous eigenvalue contribution based on face-averaged viscosity and the face-area-to-volume ratio. Computes deltaLamFace for H-correction.

Template Parameters
BDeviceBackend (Host or CUDA), defaulting to Host.
Parameters
self_viewDevice view providing mesh, BC handler, and physics access.
argPortable struct with device views of state, eigenvalue, and dt arrays.
iFaceFace index to process.
iFaceEndUpper bound of face index range.
nVarsTotal number of variables (flow + scalar).
nVarsScalarNumber of additional transported scalar variables.

need to verify this

Definition at line 732 of file EulerP_Evaluator_impl_common.hxx.

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

◆ EstEigenDt_GetFaceLam_Kernel< DeviceBackend::Host >()

template DNDS_DEVICE void DNDS::EulerP::EstEigenDt_GetFaceLam_Kernel< DeviceBackend::Host > ( EvaluatorDeviceView< DeviceBackend::Host > &  self_view,
typename Evaluator_impl< DeviceBackend::Host >::EstEigenDt_Arg::Portable &  arg,
index  iFace,
index  iFaceEnd,
int  nVars,
int  nVarsScalar 
)

◆ Flux2nd_Kernel_Face2Cell()

template<DeviceBackend B = DeviceBackend::Host>
DNDS_DEVICE void DNDS::EulerP::Flux2nd_Kernel_Face2Cell ( EvaluatorDeviceView< B > &  self_view,
typename Evaluator_impl< B >::Flux2nd_Arg::Portable &  arg,
index  iCell,
index  iCellEnd,
int  nVars,
int  nVarsScalar 
)

Scatters face fluxes to cell RHS residual (per-cell).

For each cell, accumulates the face-level numerical fluxes into the cell's RHS vector. Each face flux is weighted by -Area / Volume * sign, where sign accounts for the face normal orientation relative to the cell. Also accumulates scalar fluxes into rhsScalar.

Template Parameters
BDeviceBackend (Host or CUDA), defaulting to Host.
Parameters
self_viewDevice view providing mesh and FV geometry access.
argPortable struct with device views of flux and RHS arrays.
iCellCell index to process.
iCellEndUpper bound of cell index range.
nVarsTotal number of variables (flow + scalar).
nVarsScalarNumber of additional transported scalar variables.

Definition at line 1215 of file EulerP_Evaluator_impl_common.hxx.

Here is the caller graph for this function:

◆ Flux2nd_Kernel_Face2Cell< DeviceBackend::Host >()

template DNDS_DEVICE void DNDS::EulerP::Flux2nd_Kernel_Face2Cell< DeviceBackend::Host > ( EvaluatorDeviceView< DeviceBackend::Host > &  self_view,
typename Evaluator_impl< DeviceBackend::Host >::Flux2nd_Arg::Portable &  arg,
index  iCell,
index  iCellEnd,
int  nVars,
int  nVarsScalar 
)

◆ Flux2nd_Kernel_FluxFace()

template<DeviceBackend B = DeviceBackend::Host>
DNDS_DEVICE void DNDS::EulerP::Flux2nd_Kernel_FluxFace ( EvaluatorDeviceView< B > &  self_view,
typename Evaluator_impl< B >::Flux2nd_Arg::Portable &  arg,
index  iFace,
index  iFaceEnd,
int  nVars,
int  nVarsScalar 
)

2nd-order Roe inviscid flux computation kernel (per-face).

For each face, computes the Roe-averaged state (velocity, enthalpy, speed of sound), evaluates the three characteristic eigenvalues with H-correction fixing, and assembles the numerical flux via RoeFluxFlow. The result is stored in fluxFF.

Note
Scalar flux and viscous flux are not yet implemented (marked TODO).
Template Parameters
BDeviceBackend (Host or CUDA), defaulting to Host.
Parameters
self_viewDevice view providing mesh, BC handler, and physics access.
argPortable struct with device views of all state, face, thermodynamic, and flux arrays.
iFaceFace index to process.
iFaceEndUpper bound of face index range.
nVarsTotal number of variables (flow + scalar).
nVarsScalarNumber of additional transported scalar variables.

Definition at line 1035 of file EulerP_Evaluator_impl_common.hxx.

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

◆ Flux2nd_Kernel_FluxFace< DeviceBackend::Host >()

template DNDS_DEVICE void DNDS::EulerP::Flux2nd_Kernel_FluxFace< DeviceBackend::Host > ( EvaluatorDeviceView< DeviceBackend::Host > &  self_view,
typename Evaluator_impl< DeviceBackend::Host >::Flux2nd_Arg::Portable &  arg,
index  iFace,
index  iFaceEnd,
int  nVars,
int  nVarsScalar 
)

◆ GasInviscidFlux_XY()

template<typename TU , typename TF , class TVecN >
DNDS_DEVICE_CALLABLE void DNDS::EulerP::GasInviscidFlux_XY ( TU &&  U,
int  nVars,
real  vn,
real  vgn,
TVecN &&  n,
real  p,
TF &  F 
)
inline

Computes the inviscid (Euler) flux projected onto a face normal direction.

Accumulates the inviscid flux contribution into F:

  • F_i += U_i * (vn - vgn) for all variables (advection with grid velocity correction)
  • F_{d+1} += p * n_d for d = 0,1,2 (pressure contribution to momentum)
  • F_{I4} += vn * p (pressure work on energy)
Note
This function adds to F (does not zero it first). Caller must initialize F.
Template Parameters
TUConservative state vector type (deduced).
TFFlux vector type (deduced).
TVecNNormal vector type (deduced).
Parameters
UConservative state vector (rho, rhoU, rhoV, rhoW, E, ...).
nVarsTotal number of variables.
vnNormal velocity v dot n.
vgnGrid normal velocity (for ALE/moving mesh; 0 for static mesh).
nOutward unit face normal vector (3 components).
pPressure at the face.
[in,out]FFlux vector to accumulate into.

Definition at line 427 of file EulerP_Physics.hpp.

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

◆ pybind11_BC_bind()

void DNDS::EulerP::pybind11_BC_bind ( py::module_ &  m)
inline

Top-level binding function for all EulerP boundary condition Python API.

Calls pybind11_BCType_define, pybind11_BCInput_define, pybind11_BC_define, and pybind11_BCHandler_define to register the complete BC Python interface.

Parameters
mPybind11 module to register bindings into.

Definition at line 140 of file EulerP_BC_bind.hpp.

Here is the call graph for this function:

◆ pybind11_BC_define()

void DNDS::EulerP::pybind11_BC_define ( py::module_ &  m)
inline

Registers pybind11 bindings for the BC class.

Exposes BC with:

  • Default constructor
  • id (read/write property via getId/setId)
  • type (read/write property via getType/setType)
  • values (read/write property via getValues/setValues)
Parameters
mPybind11 module to register bindings into.

Definition at line 96 of file EulerP_BC_bind.hpp.

Here is the caller graph for this function:

◆ pybind11_BCHandler_define()

void DNDS::EulerP::pybind11_BCHandler_define ( py::module_ &  m)
inline

Registers pybind11 bindings for the BCHandler class.

Exposes BCHandler with:

  • Constructor taking (bc_inputs: list[BCInput], name2id: AutoAppendName2ID)
  • id2bc(id) method for BC lookup by zone ID
Parameters
mPybind11 module to register bindings into.

Definition at line 119 of file EulerP_BC_bind.hpp.

Here is the caller graph for this function:

◆ pybind11_BCInput_define()

void DNDS::EulerP::pybind11_BCInput_define ( py::module_ &  m)
inline

Registers pybind11 bindings for the BCInput struct.

Exposes BCInput with:

  • Default constructor
  • to_dict() / from_dict() for JSON round-trip
  • name, type, value as read/write properties
Parameters
mPybind11 module to register bindings into.

Definition at line 59 of file EulerP_BC_bind.hpp.

Here is the caller graph for this function:

◆ pybind11_BCType_define()

void DNDS::EulerP::pybind11_BCType_define ( py::module_ &  m)
inline

Registers the BCType enum as a Python enum class.

Exposes all BCType values: Wall, WallInvis, WallIsothermal, Far, Sym, In, InPsTs, Out, OutP, Special, Unknown.

Parameters
mPybind11 module to register the enum into.

Definition at line 30 of file EulerP_BC_bind.hpp.

Here is the caller graph for this function:

◆ pybind11_Evaluator_bind()

void DNDS::EulerP::pybind11_Evaluator_bind ( py::module_ &  m)
inline

Top-level binding function for the EulerP Evaluator Python API.

Calls pybind11_Evaluator_define to register all Evaluator bindings.

Parameters
mPybind11 module to register bindings into.

Definition at line 103 of file EulerP_Evaluator_bind.hpp.

Here is the call graph for this function:

◆ pybind11_Evaluator_define()

void DNDS::EulerP::pybind11_Evaluator_define ( py::module_ &  m)
inline

Defines pybind11 bindings for the Evaluator class and its kernel argument structs.

Registers the Evaluator class with:

  • Constructor taking (fv, bcHandler, physics)
  • setConfig / getConfig for JSON configuration round-trip
  • fv, bcHandler, physics as read/write properties
  • to_host, to_device, device for device management
  • PrintDataVTKHDF for VTK output
  • Kernel methods and their packed argument structs via DNDS_EULERP_EVALUATOR_BIND_STANDARD_PACKED_API
Parameters
mPybind11 module to register bindings into.

Definition at line 35 of file EulerP_Evaluator_bind.hpp.

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

◆ pybind11_Physics_bind()

void DNDS::EulerP::pybind11_Physics_bind ( py::module_ &  m)
inline

Registers pybind11 bindings for the Physics class.

Exposes Physics as a shared_ptr-held class with:

  • Physics() default constructor
  • to_dict() → Python dict (via nlohmann_json serialization)
  • Physics.from_dict(dict)Physics (static factory via JSON deserialization)
Parameters
mPybind11 module to register bindings into.

Definition at line 28 of file EulerP_Physics_bind.hpp.

◆ RecFace2nd_Kernel()

template<DeviceBackend B = DeviceBackend::Host>
DNDS_DEVICE void DNDS::EulerP::RecFace2nd_Kernel ( EvaluatorDeviceView< B > &  self_view,
typename Evaluator_impl< B >::RecFace2nd_Arg::Portable &  arg,
index  iFace,
index  iFaceEnd,
int  nVars,
int  nVarsScalar 
)

2nd-order face value reconstruction kernel (per-face).

For each face, extrapolates cell-centered values using the cell gradient: u_face = u_cell + grad(u) . (x_face - x_cell). For internal faces, computes both left and right states and a face-averaged gradient with a correction term proportional to the jump (uR - uL) / distance. For boundary faces, applies the BC handler to generate the right (ghost) state. Also reconstructs transported scalar face values and gradients.

Template Parameters
BDeviceBackend (Host or CUDA), defaulting to Host.
Parameters
self_viewDevice view providing mesh, BC handler, and physics access.
argPortable struct with device views of state, gradient, and face output arrays.
iFaceFace index to process.
iFaceEndUpper bound of face index range.
nVarsTotal number of variables (flow + scalar).
nVarsScalarNumber of additional transported scalar variables.

Definition at line 886 of file EulerP_Evaluator_impl_common.hxx.

Here is the caller graph for this function:

◆ RecFace2nd_Kernel< DeviceBackend::Host >()

template DNDS_DEVICE void DNDS::EulerP::RecFace2nd_Kernel< DeviceBackend::Host > ( EvaluatorDeviceView< DeviceBackend::Host > &  self_view,
typename Evaluator_impl< DeviceBackend::Host >::RecFace2nd_Arg::Portable &  arg,
index  iFace,
index  iFaceEnd,
int  nVars,
int  nVarsScalar 
)

◆ RecGradient_BarthLimiter_Kernel_FlowPart()

template<DeviceBackend B = DeviceBackend::Host>
DNDS_DEVICE void DNDS::EulerP::RecGradient_BarthLimiter_Kernel_FlowPart ( EvaluatorDeviceView< B > &  self_view,
typename Evaluator_impl< B >::RecGradient_Arg::Portable &  arg,
index  iCell,
index  iCellEnd,
int  nVars,
int  nVarsScalar 
)

Barth-Jespersen gradient limiter kernel for flow variables (per-cell).

Limits cell gradients to prevent the reconstructed face values from creating new extrema. Operates on density, total energy, and velocity magnitude independently. Also enforces positivity of internal energy by computing a secondary limiting factor based on the minimum reconstructed internal energy across cell faces.

Template Parameters
BDeviceBackend (Host or CUDA), defaulting to Host.
Parameters
self_viewDevice view providing mesh, BC handler, and physics access.
argPortable struct with device views of state and gradient arrays.
iCellCell index to process.
iCellEndUpper bound of cell index range.
nVarsTotal number of variables (flow + scalar).
nVarsScalarNumber of additional transported scalar variables.

Definition at line 252 of file EulerP_Evaluator_impl_common.hxx.

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

◆ RecGradient_BarthLimiter_Kernel_FlowPart< DeviceBackend::Host >()

template DNDS_DEVICE void DNDS::EulerP::RecGradient_BarthLimiter_Kernel_FlowPart< DeviceBackend::Host > ( EvaluatorDeviceView< DeviceBackend::Host > &  self_view,
typename Evaluator_impl< DeviceBackend::Host >::RecGradient_Arg::Portable &  arg,
index  iCell,
index  iCellEnd,
int  nVars,
int  nVarsScalar 
)

◆ RecGradient_BarthLimiter_Kernel_ScalarPart()

template<DeviceBackend B = DeviceBackend::Host>
DNDS_DEVICE void DNDS::EulerP::RecGradient_BarthLimiter_Kernel_ScalarPart ( EvaluatorDeviceView< B > &  self_view,
typename Evaluator_impl< B >::RecGradient_Arg::Portable &  arg,
index  iCell,
index  iCellEnd,
int  nVars,
int  nVarsScalar 
)

Barth-Jespersen gradient limiter kernel for transported scalar variables (per-cell).

Applies the same min/max bounding logic as the flow limiter but for additional transported scalar fields. Processes scalars in batches of bufSize (nVarsFlow) to reuse the TU-sized buffers.

Template Parameters
BDeviceBackend (Host or CUDA), defaulting to Host.
Parameters
self_viewDevice view providing mesh, BC handler, and physics access.
argPortable struct with device views of scalar and scalar gradient arrays.
iCellCell index to process.
iCellEndUpper bound of cell index range.
nVarsTotal number of variables (flow + scalar).
nVarsScalarNumber of additional transported scalar variables.

Definition at line 400 of file EulerP_Evaluator_impl_common.hxx.

Here is the caller graph for this function:

◆ RecGradient_BarthLimiter_Kernel_ScalarPart< DeviceBackend::Host >()

template DNDS_DEVICE void DNDS::EulerP::RecGradient_BarthLimiter_Kernel_ScalarPart< DeviceBackend::Host > ( EvaluatorDeviceView< DeviceBackend::Host > &  self_view,
typename Evaluator_impl< DeviceBackend::Host >::RecGradient_Arg::Portable &  arg,
index  iCell,
index  iCellEnd,
int  nVars,
int  nVarsScalar 
)

◆ RecGradient_GGRec_Kernel_BndVal()

template<DeviceBackend B = DeviceBackend::Host>
DNDS_DEVICE void DNDS::EulerP::RecGradient_GGRec_Kernel_BndVal ( EvaluatorDeviceView< B > &  self_view,
typename Evaluator_impl< B >::RecGradient_Arg::Portable &  arg,
index  iBnd,
index  iBndEnd,
int  nVars,
int  nVarsScalar 
)

Generates boundary ghost values for Green-Gauss gradient reconstruction.

For a given boundary face, applies the boundary condition to produce a ghost-cell state stored in the face BC buffer. Internal faces are skipped. The ghost values are later used by the Green-Gauss kernel in place of a neighbor cell's state.

Template Parameters
BDeviceBackend (Host or CUDA), defaulting to Host.
Parameters
self_viewDevice view providing mesh, BC handler, and physics access.
argPortable struct with device views of state arrays and face BC buffers.
iBndBoundary index to process.
iBndEndUpper bound of the boundary index range (used for bounds checking).
nVarsTotal number of variables (flow + scalar).
nVarsScalarNumber of additional transported scalar variables.

Definition at line 52 of file EulerP_Evaluator_impl_common.hxx.

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

◆ RecGradient_GGRec_Kernel_BndVal< DeviceBackend::Host >()

template DNDS_DEVICE void DNDS::EulerP::RecGradient_GGRec_Kernel_BndVal< DeviceBackend::Host > ( EvaluatorDeviceView< DeviceBackend::Host > &  self_view,
typename Evaluator_impl< DeviceBackend::Host >::RecGradient_Arg::Portable &  arg,
index  iBnd,
index  iBndEnd,
int  nVars,
int  nVarsScalar 
)

◆ RecGradient_GGRec_Kernel_GG()

template<DeviceBackend B = DeviceBackend::Host>
DNDS_DEVICE void DNDS::EulerP::RecGradient_GGRec_Kernel_GG ( EvaluatorDeviceView< B > &  self_view,
typename Evaluator_impl< B >::RecGradient_Arg::Portable &  arg,
index  iCell,
index  iCellEnd,
int  nVars,
int  nVarsScalar 
)

Green-Gauss gradient reconstruction kernel (per-cell).

Computes the cell gradient of the conservative state using the Green-Gauss theorem: grad(u)_cell = (1/V) * sum_faces [ 0.5 * (u_neighbor - u_cell) (x) n * A ]. Boundary faces use ghost values from the face BC buffer. Also computes gradients for transported scalar variables.

On CUDA, uses shared memory for coalesced global writes via CUDA_Local2GlobalAssign.

Template Parameters
BDeviceBackend (Host or CUDA), defaulting to Host.
Parameters
self_viewDevice view providing mesh, BC handler, and physics access.
argPortable struct with device views of state, gradient, and BC buffer arrays.
iCellCell index to process.
iCellEndUpper bound of cell index range (used for bounds checking).
nVarsTotal number of variables (flow + scalar).
nVarsScalarNumber of additional transported scalar variables.

Definition at line 127 of file EulerP_Evaluator_impl_common.hxx.

Here is the caller graph for this function:

◆ RecGradient_GGRec_Kernel_GG< DeviceBackend::Host >()

template DNDS_DEVICE void DNDS::EulerP::RecGradient_GGRec_Kernel_GG< DeviceBackend::Host > ( EvaluatorDeviceView< DeviceBackend::Host > &  self_view,
typename Evaluator_impl< DeviceBackend::Host >::RecGradient_Arg::Portable &  arg,
index  iCell,
index  iCellEnd,
int  nVars,
int  nVarsScalar 
)

◆ RoeAverageNS()

template<DeviceBackend B, class TUL , class TUR , class TULPrim , class TURPrim >
DNDS_DEVICE_CALLABLE void DNDS::EulerP::RoeAverageNS ( TUL &&  UL,
TUR &&  UR,
TULPrim &&  ULPrim,
TURPrim &&  URPrim,
int  nVars,
real  pL,
real  pR,
PhysicsDeviceView< B > &  phy,
Geom::tPoint veloRoe,
real vsqrRoe,
real HRoe,
real rhoRoe,
real aSqrRoe 
)

Computes Roe-averaged quantities from left and right states.

Calculates the Roe-averaged velocity, specific total enthalpy, density, and speed of sound squared using the density-weighted averaging formula:

  • phi_Roe = (sqrt(rho_L) * phi_L + sqrt(rho_R) * phi_R) / (sqrt(rho_L) + sqrt(rho_R))
  • rho_Roe = sqrt(rho_L * rho_R)
Note
Speed of sound squared uses IdealGas::RoeSpeedOfSoundSqr, which assumes a perfect gas.
Template Parameters
BDevice backend (Host or CUDA).
TULLeft conservative state type (deduced).
TURRight conservative state type (deduced).
TULPrimLeft primitive state type (deduced).
TURPrimRight primitive state type (deduced).
Parameters
ULLeft conservative state vector.
URRight conservative state vector.
ULPrimLeft primitive state vector.
URPrimRight primitive state vector.
nVarsTotal number of variables.
pLLeft pressure.
pRRight pressure.
phyPhysics device view for thermodynamic computations.
[out]veloRoeRoe-averaged velocity vector (3 components).
[out]vsqrRoeRoe-averaged velocity magnitude squared.
[out]HRoeRoe-averaged specific total enthalpy.
[out]rhoRoeRoe-averaged density.
[out]aSqrRoeRoe-averaged speed of sound squared.

Definition at line 73 of file EulerP_ARS.hpp.

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

◆ RoeEigenValueFixer()

DNDS_DEVICE_CALLABLE void DNDS::EulerP::RoeEigenValueFixer ( real  aL,
real  aR,
real  vnL,
real  vnR,
real  dLambda,
real  fixScale,
real lam0,
real lam123,
real lam4 
)
inline

Applies entropy fix to Roe eigenvalues to prevent expansion shocks.

Thin wrapper delegating to IdealGas::EntropyFix_HCorrHY, which is shared with the Euler module. Modifies the three characteristic eigenvalues (lam0, lam123, lam4) in-place.

Parameters
aLSpeed of sound on the left state.
aRSpeed of sound on the right state.
vnLNormal velocity on the left state.
vnRNormal velocity on the right state.
dLambdaEntropy fix threshold parameter.
fixScaleScaling factor for the entropy fix.
[in,out]lam0Eigenvalue for the u-a characteristic wave.
[in,out]lam123Eigenvalue for the u (contact/shear) waves.
[in,out]lam4Eigenvalue for the u+a characteristic wave.

Definition at line 34 of file EulerP_ARS.hpp.

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

◆ RoeFluxFlow()

template<DeviceBackend B>
DNDS_DEVICE_CALLABLE void DNDS::EulerP::RoeFluxFlow ( const TU UL,
const TU UR,
real  pL,
real  pR,
const Geom::tPoint veloRoe,
real  vsqrRoe,
real  vgn,
const Geom::tPoint n,
real  asqrRoe,
real  aRoe,
real  HRoe,
PhysicsDeviceView< B > &  phy,
real  lam0,
real  lam123,
real  lam4,
TU F 
)

Computes the complete Roe numerical flux for the 5-equation flow system.

Implements the Roe approximate Riemann solver: F = 0.5 * (F_L + F_R - |A_Roe| * dU)

The procedure:

  1. Decomposes the state jump dU = U_R - U_L into Roe characteristic variables (alpha decomposition) via IdealGas::RoeAlphaDecomposition.
  2. Scales each alpha by the corresponding entropy-fixed eigenvalue (lam0, lam123, lam4).
  3. Accumulates left and right inviscid fluxes via GasInviscidFlux_XY.
  4. Subtracts the upwind dissipation from the averaged flux.
  5. Multiplies the result by 0.5 for the final Roe flux.
Template Parameters
BDevice backend (Host or CUDA).
Parameters
ULLeft conservative state.
URRight conservative state.
pLLeft pressure.
pRRight pressure.
veloRoeRoe-averaged velocity vector (from RoeAverageNS).
vsqrRoeRoe-averaged velocity magnitude squared.
vgnGrid normal velocity (for ALE; 0 for static mesh).
nOutward unit face normal vector.
asqrRoeRoe-averaged speed of sound squared.
aRoeRoe-averaged speed of sound.
HRoeRoe-averaged specific total enthalpy.
phyPhysics device view.
lam0Entropy-fixed eigenvalue for the u-a wave.
lam123Entropy-fixed eigenvalue for the contact/shear waves.
lam4Entropy-fixed eigenvalue for the u+a wave.
[out]FOutput Roe flux vector (must be zero-initialized by caller).

Definition at line 131 of file EulerP_ARS.hpp.

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

◆ U012()

template<class TU >
DNDS_DEVICE_CALLABLE DNDS_FORCEINLINE auto DNDS::EulerP::U012 ( TU &&  v)

Extracts the first 3 components (indices 0,1,2) from a vector as a 3x1 block.

Template Parameters
TUEigen vector type (deduced).
Parameters
vVector of at least 3 components.
Returns
A 3x1 Eigen block expression referencing components (0, 1, 2).

Definition at line 72 of file EulerP.hpp.

◆ U123()

template<class TU >
DNDS_DEVICE_CALLABLE DNDS_FORCEINLINE auto DNDS::EulerP::U123 ( TU &&  v)

Extracts the momentum components (indices 1,2,3) from a state vector as a 3x1 block.

Template Parameters
TUEigen vector type (deduced).
Parameters
vState vector of at least 4 components (rho, rhoU, rhoV, rhoW, ...).
Returns
A 3x1 Eigen block expression referencing (rhoU, rhoV, rhoW).

Definition at line 60 of file EulerP.hpp.

Here is the caller graph for this function: