Nektar++
Protected Types | Protected Member Functions | Protected Attributes | List of all members
Nektar::CompressibleSolver Class Reference

#include <CompressibleSolver.h>

Inheritance diagram for Nektar::CompressibleSolver:
[legend]

Protected Types

using ND = NekDouble
 

Protected Member Functions

 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 Attributes

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

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
 

Detailed Description

Definition at line 45 of file CompressibleSolver.h.

Member Typedef Documentation

◆ ND

Definition at line 58 of file CompressibleSolver.h.

Constructor & Destructor Documentation

◆ CompressibleSolver() [1/2]

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

Session ctor.

Definition at line 40 of file CompressibleSolver.cpp.

42 : RiemannSolver(pSession), m_pointSolve(true)
43{
44 m_requiresRotation = true;
45
46 // Create equation of state object
47 std::string eosType;
48 pSession->LoadSolverInfo("EquationOfState", eosType, "IdealGas");
49 m_eos = GetEquationOfStateFactory().CreateInstance(eosType, pSession);
50 // Check if using ideal gas
51 m_idealGas = boost::iequals(eosType, "IdealGas");
52}
EquationOfStateSharedPtr m_eos
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:143
bool m_requiresRotation
Indicates whether the Riemann solver requires a rotation to be applied to the velocity fields.
SOLVER_UTILS_EXPORT RiemannSolver()
EquationOfStateFactory & GetEquationOfStateFactory()
Declaration of the equation of state factory singleton.

References Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), Nektar::GetEquationOfStateFactory(), m_eos, m_idealGas, and Nektar::SolverUtils::RiemannSolver::m_requiresRotation.

◆ CompressibleSolver() [2/2]

Nektar::CompressibleSolver::CompressibleSolver ( )
protected

Programmatic ctor.

Definition at line 54 of file CompressibleSolver.cpp.

55 : RiemannSolver(), m_pointSolve(true), m_idealGas(true)
56{
57 m_requiresRotation = true;
58}

References Nektar::SolverUtils::RiemannSolver::m_requiresRotation.

Member Function Documentation

◆ GetRoeSoundSpeed()

NekDouble Nektar::CompressibleSolver::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

Definition at line 112 of file CompressibleSolver.cpp.

117{
118 static NekDouble gamma = m_params["gamma"]();
119 NekDouble cRoe;
120 if (m_idealGas)
121 {
122 cRoe = sqrt((gamma - 1.0) * (HRoe - 0.5 * URoe2));
123 }
124 else
125 {
126 // Calculate static enthalpy of left and right states
127 NekDouble hL = eL + pL / rhoL;
128 NekDouble hR = eR + pR / rhoR;
129
130 // Get partial derivatives of P(rho,e)
131 NekDouble dpdeL = m_eos->GetDPDe_rho(rhoL, eL);
132 NekDouble dpdeR = m_eos->GetDPDe_rho(rhoR, eR);
133 NekDouble dpdrhoL = m_eos->GetDPDrho_e(rhoL, eL);
134 NekDouble dpdrhoR = m_eos->GetDPDrho_e(rhoR, eR);
135
136 // Evaluate chi and kappa parameters
137 NekDouble chiL = dpdrhoL - eL / rhoL * dpdeL;
138 NekDouble kappaL = dpdeL / rhoL;
139 NekDouble chiR = dpdrhoR - eR / rhoR * dpdeR;
140 NekDouble kappaR = dpdeR / rhoR;
141
142 //
143 // Calculate interface speed of sound using procedure from
144 // Vinokur, M.; Montagné, J.-L. "Generalized Flux-Vector
145 // Splitting and Roe Average for an Equilibrium Real Gas",
146 // JCP (1990).
147 //
148
149 // Calculate averages
150 NekDouble avgChi = 0.5 * (chiL + chiR);
151 NekDouble avgKappa = 0.5 * (kappaL + kappaR);
152 NekDouble avgKappaH = 0.5 * (kappaL * hL + kappaR * hR);
153
154 // Calculate jumps
155 NekDouble deltaP = pR - pL;
156 NekDouble deltaRho = rhoR - rhoL;
157 NekDouble deltaRhoe = rhoR * eR - rhoL * eL;
158
159 // Evaluate dP: equation (64) from Vinokur-Montagné
160 NekDouble dP = deltaP - avgChi * deltaRho - avgKappa * deltaRhoe;
161 // s (eq 66)
162 NekDouble s = avgChi + avgKappaH;
163 // D (eq 65)
164 NekDouble D = (s * deltaRho) * (s * deltaRho) + deltaP * deltaP;
165 // chiRoe and kappaRoe (eq 66)
166 NekDouble chiRoe, kappaRoe;
167 NekDouble fac = D - deltaP * deltaRho;
169 {
170 chiRoe = (D * avgChi + s * s * deltaRho * dP) / fac;
171 kappaRoe = D * avgKappa / fac;
172 }
173 else
174 {
175 chiRoe = avgChi;
176 kappaRoe = avgKappa;
177 }
178 // Speed of sound (eq 53)
179 cRoe = sqrt(chiRoe + kappaRoe * (HRoe - 0.5 * URoe2));
180 }
181 return cRoe;
182}
std::map< std::string, RSParamFuncType > m_params
Map of parameter function types.
static const NekDouble kNekZeroTol
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(), Nektar::NekConstants::kNekZeroTol, m_eos, m_idealGas, Nektar::SolverUtils::RiemannSolver::m_params, and tinysimd::sqrt().

Referenced by Nektar::HLLSolver::v_PointSolve(), Nektar::LaxFriedrichsSolver::v_PointSolve(), and Nektar::HLLCSolver::v_PointSolve().

◆ v_ArraySolve()

virtual void Nektar::CompressibleSolver::v_ArraySolve ( const Array< OneD, const Array< OneD, ND > > &  Fwd,
const Array< OneD, const Array< OneD, ND > > &  Bwd,
Array< OneD, Array< OneD, ND > > &  flux 
)
inlineprotectedvirtual

Reimplemented in Nektar::RoeSolver, and Nektar::AverageSolver.

Definition at line 64 of file CompressibleSolver.h.

68 {
70 "This function should be defined by subclasses.");
71 }
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:202

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by v_Solve().

◆ v_PointSolve()

virtual void Nektar::CompressibleSolver::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 
)
inlineprotectedvirtual

◆ v_Solve()

void Nektar::CompressibleSolver::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 
)
overrideprotectedvirtual

Implements Nektar::SolverUtils::RiemannSolver.

Reimplemented in Nektar::RoeSolverSIMD.

Definition at line 63 of file CompressibleSolver.cpp.

67{
68 if (m_pointSolve)
69 {
70 size_t expDim = nDim;
71 NekDouble rhouf{}, rhovf{};
72
73 if (expDim == 1)
74 {
75 for (size_t i = 0; i < Fwd[0].size(); ++i)
76 {
77 v_PointSolve(Fwd[0][i], Fwd[1][i], 0.0, 0.0, Fwd[2][i],
78 Bwd[0][i], Bwd[1][i], 0.0, 0.0, Bwd[2][i],
79 flux[0][i], flux[1][i], rhouf, rhovf, flux[2][i]);
80 }
81 }
82 else if (expDim == 2)
83 {
84 for (size_t i = 0; i < Fwd[0].size(); ++i)
85 {
86 v_PointSolve(Fwd[0][i], Fwd[1][i], Fwd[2][i], 0.0, Fwd[3][i],
87 Bwd[0][i], Bwd[1][i], Bwd[2][i], 0.0, Bwd[3][i],
88 flux[0][i], flux[1][i], flux[2][i], rhovf,
89 flux[3][i]);
90 }
91 }
92 else if (expDim == 3)
93 {
94 for (size_t i = 0; i < Fwd[0].size(); ++i)
95 {
96 v_PointSolve(Fwd[0][i], Fwd[1][i], Fwd[2][i], Fwd[3][i],
97 Fwd[4][i], Bwd[0][i], Bwd[1][i], Bwd[2][i],
98 Bwd[3][i], Bwd[4][i], flux[0][i], flux[1][i],
99 flux[2][i], flux[3][i], flux[4][i]);
100 }
101 }
102 }
103 else
104 {
105 v_ArraySolve(Fwd, Bwd, flux);
106 }
107}
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)
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)

References m_pointSolve, v_ArraySolve(), and v_PointSolve().

Member Data Documentation

◆ m_eos

EquationOfStateSharedPtr Nektar::CompressibleSolver::m_eos
protected

◆ m_idealGas

bool Nektar::CompressibleSolver::m_idealGas
protected

Definition at line 50 of file CompressibleSolver.h.

Referenced by CompressibleSolver(), and GetRoeSoundSpeed().

◆ m_pointSolve

bool Nektar::CompressibleSolver::m_pointSolve
protected

Definition at line 48 of file CompressibleSolver.h.

Referenced by Nektar::AverageSolver::AverageSolver(), and v_Solve().