Nektar++
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 ()
 
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 47 of file CompressibleSolver.h.

Member Typedef Documentation

◆ ND

Definition at line 61 of file CompressibleSolver.h.

Constructor & Destructor Documentation

◆ CompressibleSolver() [1/2]

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

Session ctor.

Definition at line 41 of file CompressibleSolver.cpp.

43  : RiemannSolver(pSession), m_pointSolve(true)
44  {
45  m_requiresRotation = true;
46 
47  // Create equation of state object
48  std::string eosType;
49  pSession->LoadSolverInfo("EquationOfState",
50  eosType, "IdealGas");
52  .CreateInstance(eosType, pSession);
53  // Check if using ideal gas
54  m_idealGas = boost::iequals(eosType, "IdealGas");
55  }
EquationOfStateSharedPtr m_eos
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:145
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 57 of file CompressibleSolver.cpp.

58  : RiemannSolver(), m_idealGas(true), m_pointSolve(true)
59  {
60  m_requiresRotation = true;
61  }

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 113 of file CompressibleSolver.cpp.

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

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::AverageSolver, and Nektar::RoeSolver.

Definition at line 69 of file CompressibleSolver.h.

73  {
74  boost::ignore_unused(Fwd, Bwd, flux);
76  "This function should be defined by subclasses.");
77  }
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:209

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

Reimplemented in Nektar::HLLCSolver, Nektar::ExactSolverToro, Nektar::RoeSolver, Nektar::LaxFriedrichsSolver, Nektar::HLLSolver, Nektar::AUSM3Solver, Nektar::AUSM2Solver, Nektar::AUSM1Solver, and Nektar::AUSM0Solver.

Definition at line 79 of file CompressibleSolver.h.

83  {
84  boost::ignore_unused(rhoL, rhouL, rhovL, rhowL, EL,
85  rhoR, rhouR, rhovR, rhowR, ER,
86  rhof, rhouf, rhovf, rhowf, Ef);
88  "This function should be defined by subclasses.");
89  }

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

Referenced by v_Solve().

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

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

Referenced by CompressibleSolver(), and GetRoeSoundSpeed().

◆ m_pointSolve

bool Nektar::CompressibleSolver::m_pointSolve
protected

Definition at line 50 of file CompressibleSolver.h.

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