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

#include <LEELaxFriedrichsSolver.h>

Inheritance diagram for Nektar::LEELaxFriedrichsSolver:
[legend]

Static Public Member Functions

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

Static Public Attributes

static std::string solverName
 

Protected Member Functions

 LEELaxFriedrichsSolver (const LibUtilities::SessionReaderSharedPtr &pSession)
 
void v_PointSolve (NekDouble pL, NekDouble rhoL, NekDouble rhouL, NekDouble rhovL, NekDouble rhowL, NekDouble pR, NekDouble rhoR, NekDouble rhouR, NekDouble rhovR, NekDouble rhowR, NekDouble c0sqL, NekDouble rho0L, NekDouble u0L, NekDouble v0L, NekDouble w0L, NekDouble c0sqR, NekDouble rho0R, NekDouble u0R, NekDouble v0R, NekDouble w0R, NekDouble &pF, NekDouble &rhoF, NekDouble &rhouF, NekDouble &rhovF, NekDouble &rhowF) override
 Lax-Friedrichs Riemann solver. More...
 
- Protected Member Functions inherited from Nektar::LEESolver
 LEESolver (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
 
- Protected Member Functions inherited from Nektar::AcousticSolver
 AcousticSolver (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_PointSolve (NekDouble pL, NekDouble rhoL, NekDouble uL, NekDouble vL, NekDouble wL, NekDouble pR, NekDouble rhoR, NekDouble uR, NekDouble vR, NekDouble wR, NekDouble c0sqL, NekDouble rho0L, NekDouble u0L, NekDouble v0L, NekDouble w0L, NekDouble c0sqR, NekDouble rho0R, NekDouble u0R, NekDouble v0R, NekDouble w0R, NekDouble &pF, NekDouble &rhoF, NekDouble &uF, NekDouble &vF, NekDouble &wF)=0
 
void GetRotBasefield (Array< OneD, Array< OneD, NekDouble > > &bfFwd, Array< OneD, Array< OneD, NekDouble > > &bfBwd)
 
- 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)
 

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 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 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...
 

Detailed Description

Definition at line 47 of file LEELaxFriedrichsSolver.h.

Constructor & Destructor Documentation

◆ LEELaxFriedrichsSolver()

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

Definition at line 51 of file LEELaxFriedrichsSolver.cpp.

53 : LEESolver(pSession)
54{
55}
LEESolver(const LibUtilities::SessionReaderSharedPtr &pSession)
Definition: LEESolver.cpp:46

Referenced by create().

Member Function Documentation

◆ create()

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

Definition at line 50 of file LEELaxFriedrichsSolver.h.

52 {
54 }
LEELaxFriedrichsSolver(const LibUtilities::SessionReaderSharedPtr &pSession)
std::shared_ptr< RiemannSolver > RiemannSolverSharedPtr
A shared pointer to an EquationSystem object.

References LEELaxFriedrichsSolver().

◆ v_PointSolve()

void Nektar::LEELaxFriedrichsSolver::v_PointSolve ( NekDouble  pL,
NekDouble  rhoL,
NekDouble  rhouL,
NekDouble  rhovL,
NekDouble  rhowL,
NekDouble  pR,
NekDouble  rhoR,
NekDouble  rhouR,
NekDouble  rhovR,
NekDouble  rhowR,
NekDouble  c0sqL,
NekDouble  rho0L,
NekDouble  u0L,
NekDouble  v0L,
NekDouble  w0L,
NekDouble  c0sqR,
NekDouble  rho0R,
NekDouble  u0R,
NekDouble  v0R,
NekDouble  w0R,
NekDouble pF,
NekDouble rhoF,
NekDouble rhouF,
NekDouble rhovF,
NekDouble rhowF 
)
overrideprotectedvirtual

Lax-Friedrichs Riemann solver.

Parameters
pLPerturbation pressure left state
rhoLPerturbation density left state
pRPerturbation pressure right state
rhoRPerturbation density right state
rhouLx perturbation velocity component left state
rhouRx perturbation velocity component right state
rhovLy perturbation velocity component left state
rhovRy perturbation velocity component right state
rhowLz perturbation velocity component left state
rhowRz perturbation velocity component right state
c0sqLBase pressure left state
c0sqRBase pressure right state
rho0LBase density left state
rho0RBase density right state
u0LBase x velocity component left state
u0RBase x velocity component right state
v0LBase y velocity component left state
v0RBase y velocity component right state
w0LBase z velocity component left state
w0RBase z velocity component right state
pFComputed Riemann flux for perturbation pressure
rhoFComputed Riemann flux for perturbation density
rhouFComputed Riemann flux for x perturbation velocity component
rhovFComputed Riemann flux for y perturbation velocity component
rhowFComputed Riemann flux for z perturbation velocity component

Implements Nektar::AcousticSolver.

Definition at line 86 of file LEELaxFriedrichsSolver.cpp.

96{
97 // Speed of sound
98 NekDouble cL = sqrt(c0sqL);
99 NekDouble cR = sqrt(c0sqR);
100
101 // max absolute eigenvalue of the jacobian of F_n1
102 NekDouble a_1_max = 0;
103 a_1_max = std::max(a_1_max, std::abs(u0L - cL));
104 a_1_max = std::max(a_1_max, std::abs(u0R - cR));
105 a_1_max = std::max(a_1_max, std::abs(u0L + cL));
106 a_1_max = std::max(a_1_max, std::abs(u0R + cR));
107
108 NekDouble pFL = rhouL * c0sqL + u0L * pL;
109 NekDouble rhoFL = rhouL + rhoL * u0L;
110 NekDouble rhouFL = pL + rhouL * u0L;
111 NekDouble rhovFL = 0;
112 NekDouble rhowFL = 0;
113
114 NekDouble pFR = rhouR * c0sqR + u0R * pR;
115 NekDouble rhoFR = rhouR + rhoR * u0R;
116 NekDouble rhouFR = pR + rhouR * u0R;
117 NekDouble rhovFR = 0;
118 NekDouble rhowFR = 0;
119
120 // assemble the face-normal fluxes
121 pF = 0.5 * (pFL + pFR - a_1_max * (pR - pL));
122 rhoF = 0.5 * (rhoFL + rhoFR - a_1_max * (rhoR - rhoL));
123 rhouF = 0.5 * (rhouFL + rhouFR - a_1_max * (rhouR - rhouL));
124 rhovF = 0.5 * (rhovFL + rhovFR - a_1_max * (rhovR - rhovL));
125 rhowF = 0.5 * (rhowFL + rhowFR - a_1_max * (rhowR - rhowL));
126}
double NekDouble
scalarT< T > abs(scalarT< T > in)
Definition: scalar.hpp:298
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294

References tinysimd::abs(), and tinysimd::sqrt().

Member Data Documentation

◆ solverName

std::string Nektar::LEELaxFriedrichsSolver::solverName
static
Initial value:
=
"LEELaxFriedrichs", LEELaxFriedrichsSolver::create,
"Lax-Friedrichs Solver for LEE")
static RiemannSolverSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession)
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:197
RiemannSolverFactory & GetRiemannSolverFactory()

Definition at line 56 of file LEELaxFriedrichsSolver.h.