Nektar++
Static Public Member Functions | Static Public Attributes | Protected Member Functions | List of all members
Nektar::LaxFriedrichsSolver Class Reference

#include <LaxFriedrichsSolver.h>

Inheritance diagram for Nektar::LaxFriedrichsSolver:
[legend]

Static Public Member Functions

static RiemannSolverSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession)
 
static RiemannSolverSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession)
 

Static Public Attributes

static std::string solverName
 

Protected Member Functions

 LaxFriedrichsSolver (const LibUtilities::SessionReaderSharedPtr &pSession)
 
void v_PointSolve (double rhoL, double rhouL, double rhovL, double rhowL, double EL, double rhoR, double rhouR, double rhovR, double rhowR, double ER, double &rhof, double &rhouf, double &rhovf, double &rhowf, double &Ef) override
 Lax-Friedrichs Riemann solver. More...
 
 LaxFriedrichsSolver (const LibUtilities::SessionReaderSharedPtr &pSession)
 
void v_PointSolve (NekDouble hL, NekDouble huL, NekDouble hvL, NekDouble hR, NekDouble huR, NekDouble hvR, NekDouble &hf, NekDouble &huf, NekDouble &hvf) override
 Lax-Friedrichs Riemann solver. More...
 
- Protected Member Functions inherited from Nektar::CompressibleSolver
 CompressibleSolver (const LibUtilities::SessionReaderSharedPtr &pSession)
 Session ctor. More...
 
 CompressibleSolver ()
 Programmatic ctor. More...
 
void v_Solve (const int nDim, const Array< OneD, const Array< OneD, ND > > &Fwd, const Array< OneD, const Array< OneD, ND > > &Bwd, Array< OneD, Array< OneD, ND > > &flux) override
 
virtual void v_ArraySolve (const Array< OneD, const Array< OneD, ND > > &Fwd, const Array< OneD, const Array< OneD, ND > > &Bwd, Array< OneD, Array< OneD, ND > > &flux)
 
virtual void v_PointSolve (ND rhoL, ND rhouL, ND rhovL, ND rhowL, ND EL, ND rhoR, ND rhouR, ND rhovR, ND rhowR, ND ER, ND &rhof, ND &rhouf, ND &rhovf, ND &rhowf, ND &Ef)
 
ND GetRoeSoundSpeed (ND rhoL, ND pL, ND eL, ND HL, ND srL, ND rhoR, ND pR, ND eR, ND HR, ND srR, ND HRoe, ND URoe2, ND srLR)
 
- Protected Member Functions inherited from Nektar::SolverUtils::RiemannSolver
SOLVER_UTILS_EXPORT RiemannSolver ()
 
SOLVER_UTILS_EXPORT RiemannSolver (const LibUtilities::SessionReaderSharedPtr &pSession)
 
virtual SOLVER_UTILS_EXPORT ~RiemannSolver ()
 
virtual void v_Solve (const int nDim, const Array< OneD, const Array< OneD, NekDouble > > &Fwd, const Array< OneD, const Array< OneD, NekDouble > > &Bwd, Array< OneD, Array< OneD, NekDouble > > &flux)=0
 
SOLVER_UTILS_EXPORT void GenerateRotationMatrices (const Array< OneD, const Array< OneD, NekDouble > > &normals)
 Generate rotation matrices for 3D expansions. More...
 
void FromToRotation (Array< OneD, const NekDouble > &from, Array< OneD, const NekDouble > &to, NekDouble *mat)
 A function for creating a rotation matrix that rotates a vector from into another vector to. More...
 
SOLVER_UTILS_EXPORT void rotateToNormal (const Array< OneD, const Array< OneD, NekDouble > > &inarray, const Array< OneD, const Array< OneD, NekDouble > > &normals, const Array< OneD, const Array< OneD, NekDouble > > &vecLocs, Array< OneD, Array< OneD, NekDouble > > &outarray)
 Rotate a vector field to trace normal. More...
 
SOLVER_UTILS_EXPORT void rotateFromNormal (const Array< OneD, const Array< OneD, NekDouble > > &inarray, const Array< OneD, const Array< OneD, NekDouble > > &normals, const Array< OneD, const Array< OneD, NekDouble > > &vecLocs, Array< OneD, Array< OneD, NekDouble > > &outarray)
 Rotate a vector field from trace normal. More...
 
SOLVER_UTILS_EXPORT bool CheckScalars (std::string name)
 Determine whether a scalar has been defined in m_scalars. More...
 
SOLVER_UTILS_EXPORT bool CheckVectors (std::string name)
 Determine whether a vector has been defined in m_vectors. More...
 
SOLVER_UTILS_EXPORT bool CheckParams (std::string name)
 Determine whether a parameter has been defined in m_params. More...
 
SOLVER_UTILS_EXPORT bool CheckAuxScal (std::string name)
 Determine whether a scalar has been defined in m_auxScal. More...
 
SOLVER_UTILS_EXPORT bool CheckAuxVec (std::string name)
 Determine whether a vector has been defined in m_auxVec. More...
 
virtual SOLVER_UTILS_EXPORT void v_CalcFluxJacobian (const int nDim, const Array< OneD, const Array< OneD, NekDouble > > &Fwd, const Array< OneD, const Array< OneD, NekDouble > > &Bwd, const Array< OneD, const Array< OneD, NekDouble > > &normals, DNekBlkMatSharedPtr &FJac, DNekBlkMatSharedPtr &BJac)
 
- Protected Member Functions inherited from Nektar::NonlinearSWESolver
 NonlinearSWESolver (const LibUtilities::SessionReaderSharedPtr &pSession)
 
void v_Solve (const int nDim, const Array< OneD, const Array< OneD, NekDouble > > &Fwd, const Array< OneD, const Array< OneD, NekDouble > > &Bwd, Array< OneD, Array< OneD, NekDouble > > &flux) override
 
virtual void v_ArraySolve (const Array< OneD, const Array< OneD, NekDouble > > &Fwd, const Array< OneD, const Array< OneD, NekDouble > > &Bwd, Array< OneD, Array< OneD, NekDouble > > &flux)
 
virtual void v_PointSolve (NekDouble hL, NekDouble huL, NekDouble hvL, NekDouble hR, NekDouble huR, NekDouble hvR, NekDouble &hf, NekDouble &huf, NekDouble &hvf)
 

Additional Inherited Members

- Public Member Functions inherited from Nektar::SolverUtils::RiemannSolver
SOLVER_UTILS_EXPORT void Solve (const int nDim, const Array< OneD, const Array< OneD, NekDouble > > &Fwd, const Array< OneD, const Array< OneD, NekDouble > > &Bwd, Array< OneD, Array< OneD, NekDouble > > &flux)
 Perform the Riemann solve given the forwards and backwards spaces. More...
 
template<typename FuncPointerT , typename ObjectPointerT >
void SetScalar (std::string name, FuncPointerT func, ObjectPointerT obj)
 
void SetScalar (std::string name, RSScalarFuncType fp)
 
template<typename FuncPointerT , typename ObjectPointerT >
void SetVector (std::string name, FuncPointerT func, ObjectPointerT obj)
 
void SetVector (std::string name, RSVecFuncType fp)
 
template<typename FuncPointerT , typename ObjectPointerT >
void SetParam (std::string name, FuncPointerT func, ObjectPointerT obj)
 
void SetALEFlag (bool &ALE)
 
void SetParam (std::string name, RSParamFuncType fp)
 
template<typename FuncPointerT , typename ObjectPointerT >
void SetAuxScal (std::string name, FuncPointerT func, ObjectPointerT obj)
 
template<typename FuncPointerT , typename ObjectPointerT >
void SetAuxVec (std::string name, FuncPointerT func, ObjectPointerT obj)
 
void SetAuxVec (std::string name, RSVecFuncType fp)
 
std::map< std::string, RSScalarFuncType > & GetScalars ()
 
std::map< std::string, RSVecFuncType > & GetVectors ()
 
std::map< std::string, RSParamFuncType > & GetParams ()
 
SOLVER_UTILS_EXPORT void CalcFluxJacobian (const int nDim, const Array< OneD, const Array< OneD, NekDouble > > &Fwd, const Array< OneD, const Array< OneD, NekDouble > > &Bwd, DNekBlkMatSharedPtr &FJac, DNekBlkMatSharedPtr &BJac)
 Calculate the flux jacobian of Fwd and Bwd. More...
 
- Public Attributes inherited from Nektar::SolverUtils::RiemannSolver
int m_spacedim
 
- Protected Types inherited from Nektar::CompressibleSolver
using ND = NekDouble
 
- Protected Attributes inherited from Nektar::CompressibleSolver
bool m_pointSolve
 
EquationOfStateSharedPtr m_eos
 
bool m_idealGas
 
- Protected Attributes inherited from Nektar::SolverUtils::RiemannSolver
bool m_requiresRotation
 Indicates whether the Riemann solver requires a rotation to be applied to the velocity fields. More...
 
std::map< std::string, RSScalarFuncTypem_scalars
 Map of scalar function types. More...
 
std::map< std::string, RSVecFuncTypem_vectors
 Map of vector function types. More...
 
std::map< std::string, RSParamFuncTypem_params
 Map of parameter function types. More...
 
std::map< std::string, RSScalarFuncTypem_auxScal
 Map of auxiliary scalar function types. More...
 
std::map< std::string, RSVecFuncTypem_auxVec
 Map of auxiliary vector function types. More...
 
Array< OneD, Array< OneD, NekDouble > > m_rotMat
 Rotation matrices for each trace quadrature point. More...
 
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_rotStorage
 Rotation storage. More...
 
bool m_ALESolver = false
 Flag if using the ALE formulation. More...
 
- Protected Attributes inherited from Nektar::NonlinearSWESolver
bool m_pointSolve
 

Detailed Description

Definition at line 42 of file CompressibleFlowSolver/RiemannSolvers/LaxFriedrichsSolver.h.

Constructor & Destructor Documentation

◆ LaxFriedrichsSolver() [1/2]

Nektar::LaxFriedrichsSolver::LaxFriedrichsSolver ( const LibUtilities::SessionReaderSharedPtr pSession)
protected

Definition at line 44 of file CompressibleFlowSolver/RiemannSolvers/LaxFriedrichsSolver.cpp.

46 : CompressibleSolver(pSession)
47{
48}
CompressibleSolver()
Programmatic ctor.

Referenced by create().

◆ LaxFriedrichsSolver() [2/2]

Nektar::LaxFriedrichsSolver::LaxFriedrichsSolver ( const LibUtilities::SessionReaderSharedPtr pSession)
protected

Member Function Documentation

◆ create() [1/2]

static RiemannSolverSharedPtr Nektar::LaxFriedrichsSolver::create ( const LibUtilities::SessionReaderSharedPtr pSession)
inlinestatic

Definition at line 45 of file CompressibleFlowSolver/RiemannSolvers/LaxFriedrichsSolver.h.

47 {
48 return RiemannSolverSharedPtr(new LaxFriedrichsSolver(pSession));
49 }
LaxFriedrichsSolver(const LibUtilities::SessionReaderSharedPtr &pSession)
std::shared_ptr< RiemannSolver > RiemannSolverSharedPtr
A shared pointer to an EquationSystem object.

References LaxFriedrichsSolver().

◆ create() [2/2]

static RiemannSolverSharedPtr Nektar::LaxFriedrichsSolver::create ( const LibUtilities::SessionReaderSharedPtr pSession)
inlinestatic

Definition at line 45 of file ShallowWaterSolver/RiemannSolvers/LaxFriedrichsSolver.h.

47 {
48 return RiemannSolverSharedPtr(new LaxFriedrichsSolver(pSession));
49 }

References LaxFriedrichsSolver().

◆ v_PointSolve() [1/2]

void Nektar::LaxFriedrichsSolver::v_PointSolve ( double  rhoL,
double  rhouL,
double  rhovL,
double  rhowL,
double  EL,
double  rhoR,
double  rhouR,
double  rhovR,
double  rhowR,
double  ER,
double &  rhof,
double &  rhouf,
double &  rhovf,
double &  rhowf,
double &  Ef 
)
overrideprotectedvirtual

Lax-Friedrichs Riemann solver.

Parameters
rhoLDensity left state.
rhoRDensity right state.
rhouLx-momentum component left state.
rhouRx-momentum component right state.
rhovLy-momentum component left state.
rhovRy-momentum component right state.
rhowLz-momentum component left state.
rhowRz-momentum component right state.
ELEnergy left state.
EREnergy right state.
rhofComputed Riemann flux for density.
rhoufComputed Riemann flux for x-momentum component
rhovfComputed Riemann flux for y-momentum component
rhowfComputed Riemann flux for z-momentum component
EfComputed Riemann flux for energy.

Reimplemented from Nektar::CompressibleSolver.

Definition at line 69 of file CompressibleFlowSolver/RiemannSolvers/LaxFriedrichsSolver.cpp.

74{
75 // Left and right velocities
76 NekDouble uL = rhouL / rhoL;
77 NekDouble vL = rhovL / rhoL;
78 NekDouble wL = rhowL / rhoL;
79 NekDouble uR = rhouR / rhoR;
80 NekDouble vR = rhovR / rhoR;
81 NekDouble wR = rhowR / rhoR;
82
83 // Internal energy (per unit mass)
84 NekDouble eL = (EL - 0.5 * (rhouL * uL + rhovL * vL + rhowL * wL)) / rhoL;
85 NekDouble eR = (ER - 0.5 * (rhouR * uR + rhovR * vR + rhowR * wR)) / rhoR;
86
87 // Pressure
88 NekDouble pL = m_eos->GetPressure(rhoL, eL);
89 NekDouble pR = m_eos->GetPressure(rhoR, eR);
90
91 // Left and right total enthalpy
92 NekDouble HL = (EL + pL) / rhoL;
93 NekDouble HR = (ER + pR) / rhoR;
94
95 // Square root of rhoL and rhoR.
96 NekDouble srL = sqrt(rhoL);
97 NekDouble srR = sqrt(rhoR);
98 NekDouble srLR = srL + srR;
99
100 // Roe average state
101 NekDouble uRoe = (srL * uL + srR * uR) / srLR;
102 NekDouble vRoe = (srL * vL + srR * vR) / srLR;
103 NekDouble wRoe = (srL * wL + srR * wR) / srLR;
104 NekDouble URoe2 = uRoe * uRoe + vRoe * vRoe + wRoe * wRoe;
105 NekDouble HRoe = (srL * HL + srR * HR) / srLR;
106 NekDouble cRoe = GetRoeSoundSpeed(rhoL, pL, eL, HL, srL, rhoR, pR, eR, HR,
107 srR, HRoe, URoe2, srLR);
108
109 // Maximum eigenvalue
110 NekDouble URoe = fabs(uRoe) + cRoe;
111
112 // Lax-Friedrichs flux formula
113 rhof = 0.5 * (rhouL + rhouR - URoe * (rhoR - rhoL));
114 rhouf = 0.5 * (pL + rhouL * uL + pR + rhouR * uR - URoe * (rhouR - rhouL));
115 rhovf = 0.5 * (rhouL * vL + rhouR * vR - URoe * (rhovR - rhovL));
116 rhowf = 0.5 * (rhouL * wL + rhouR * wR - URoe * (rhowR - rhowL));
117 Ef = 0.5 * (uL * (EL + pL) + uR * (ER + pR) - URoe * (ER - EL));
118}
ND GetRoeSoundSpeed(ND rhoL, ND pL, ND eL, ND HL, ND srL, ND rhoR, ND pR, ND eR, ND HR, ND srR, ND HRoe, ND URoe2, ND srLR)
EquationOfStateSharedPtr m_eos
double NekDouble
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294

References Nektar::CompressibleSolver::GetRoeSoundSpeed(), Nektar::CompressibleSolver::m_eos, and tinysimd::sqrt().

◆ v_PointSolve() [2/2]

void Nektar::LaxFriedrichsSolver::v_PointSolve ( NekDouble  hL,
NekDouble  huL,
NekDouble  hvL,
NekDouble  hR,
NekDouble  huR,
NekDouble  hvR,
NekDouble hf,
NekDouble huf,
NekDouble hvf 
)
overrideprotectedvirtual

Lax-Friedrichs Riemann solver.

Parameters
hLWater depth left state.
hRWater depth right state.
huLx-momentum component left state.
huRx-momentum component right state.
hvLy-momentum component left state.
hvRy-momentum component right state.
hfComputed Riemann flux for density.
hufComputed Riemann flux for x-momentum component
hvfComputed Riemann flux for y-momentum component

Reimplemented from Nektar::NonlinearSWESolver.

Definition at line 63 of file ShallowWaterSolver/RiemannSolvers/LaxFriedrichsSolver.cpp.

68{
69 static NekDouble g = m_params["gravity"]();
70
71 // Left and right velocities
72 NekDouble uL = huL / hL;
73 NekDouble vL = hvL / hL;
74 NekDouble uR = huR / hR;
75 NekDouble vR = hvR / hR;
76
77 // Left and right wave speeds
78 NekDouble cL = sqrt(g * hL);
79 NekDouble cR = sqrt(g * hR);
80
81 // Square root of hL and hR.
82 NekDouble srL = sqrt(hL);
83 NekDouble srR = sqrt(hR);
84 NekDouble srLR = srL + srR;
85
86 // Velocity Roe averages
87 NekDouble uRoe = (srL * uL + srR * uR) / srLR;
88 NekDouble cRoe = sqrt(0.5 * (cL * cL + cR * cR));
89
90 // Minimum and maximum wave speeds
91 NekDouble S = std::max(uRoe + cRoe, std::max(uR + cR, -(uL - cL)));
92 NekDouble sign = 1.0;
93
94 if (S == -(uL - cL))
95 {
96 sign = -1.0;
97 }
98
99 // Lax-Friedrichs Riemann h flux
100 hf = 0.5 * ((huL + huR) - sign * S * (hR - hL));
101
102 // Lax-Friedrichs Riemann hu flux
103 huf =
104 0.5 *
105 ((hL * uL * uL + 0.5 * g * hL * hL + hR * uR * uR + 0.5 * g * hR * hR) -
106 sign * S * (huR - huL));
107
108 // Lax-Friedrichs Riemann hv flux
109 hvf = 0.5 * ((hL * uL * vL + hR * uR * vR) - sign * S * (hvR - hvL));
110}
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:47
std::map< std::string, RSParamFuncType > m_params
Map of parameter function types.

References sign, and tinysimd::sqrt().

Member Data Documentation

◆ solverName

static std::string Nektar::LaxFriedrichsSolver::solverName
static
Initial value:
=
"LaxFriedrichs", LaxFriedrichsSolver::create,
"Lax-Friedrichs Riemann solver")
static RiemannSolverSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession)
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
RiemannSolverFactory & GetRiemannSolverFactory()

Definition at line 51 of file CompressibleFlowSolver/RiemannSolvers/LaxFriedrichsSolver.h.