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:
Inheritance graph
[legend]
Collaboration diagram for Nektar::LaxFriedrichsSolver:
Collaboration graph
[legend]

Static Public Member Functions

static RiemannSolverSharedPtr create ()
 
static RiemannSolverSharedPtr create ()
 
static RiemannSolverSharedPtr create ()
 

Static Public Attributes

static std::string solverName
 

Protected Member Functions

 LaxFriedrichsSolver ()
 
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 p0, NekDouble rho0, NekDouble u0, NekDouble v0, NekDouble w0, NekDouble &pF, NekDouble &rhoF, NekDouble &uF, NekDouble &vF, NekDouble &wF)
 Lax-Friedrichs Riemann solver. More...
 
 LaxFriedrichsSolver ()
 
virtual void v_PointSolve (NekDouble rhoL, NekDouble rhouL, NekDouble rhovL, NekDouble rhowL, NekDouble EL, NekDouble rhoR, NekDouble rhouR, NekDouble rhovR, NekDouble rhowR, NekDouble ER, NekDouble &rhof, NekDouble &rhouf, NekDouble &rhovf, NekDouble &rhowf, NekDouble &Ef)
 Lax-Friedrichs Riemann solver. More...
 
virtual void v_PointSolveVisc (NekDouble rhoL, NekDouble rhouL, NekDouble rhovL, NekDouble rhowL, NekDouble EL, NekDouble EpsL, NekDouble rhoR, NekDouble rhouR, NekDouble rhovR, NekDouble rhowR, NekDouble ER, NekDouble EpsR, NekDouble &rhof, NekDouble &rhouf, NekDouble &rhovf, NekDouble &rhowf, NekDouble &Ef, NekDouble &Epsf)
 
 LaxFriedrichsSolver ()
 
virtual void v_PointSolve (double hL, double huL, double hvL, double hR, double huR, double hvR, double &hf, double &huf, double &hvf)
 Lax-Friedrichs Riemann solver. More...
 
- Protected Member Functions inherited from Nektar::NonlinearSWESolver
 NonlinearSWESolver ()
 
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)
 
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)
 
- Protected Member Functions inherited from Nektar::SolverUtils::RiemannSolver
SOLVER_UTILS_EXPORT RiemannSolver ()
 
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...
 
- Protected Member Functions inherited from Nektar::CompressibleSolver
 CompressibleSolver ()
 
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)
 
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)
 
- Protected Member Functions inherited from Nektar::APESolver
 APESolver ()
 
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)
 
Array< OneD, Array< OneD, NekDouble > > GetRotBasefield ()
 

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)
 
std::map< std::string, RSScalarFuncType > & GetScalars ()
 
std::map< std::string, RSVecFuncType > & GetVectors ()
 
std::map< std::string, RSParamFuncType > & GetParams ()
 
- Public Attributes inherited from Nektar::SolverUtils::RiemannSolver
int m_spacedim
 
- Protected Attributes inherited from Nektar::NonlinearSWESolver
bool m_pointSolve
 
- 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...
 
- Protected Attributes inherited from Nektar::CompressibleSolver
bool m_pointSolve
 

Detailed Description

Definition at line 48 of file APESolver/RiemannSolvers/LaxFriedrichsSolver.h.

Constructor & Destructor Documentation

Nektar::LaxFriedrichsSolver::LaxFriedrichsSolver ( )
protected

Definition at line 53 of file APESolver/RiemannSolvers/LaxFriedrichsSolver.cpp.

Referenced by create().

53  :
54  APESolver()
55 {
56 
57 }
Nektar::LaxFriedrichsSolver::LaxFriedrichsSolver ( )
protected
Nektar::LaxFriedrichsSolver::LaxFriedrichsSolver ( )
protected

Member Function Documentation

static RiemannSolverSharedPtr Nektar::LaxFriedrichsSolver::create ( )
inlinestatic

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

References LaxFriedrichsSolver().

47  {
49  }
boost::shared_ptr< RiemannSolver > RiemannSolverSharedPtr
A shared pointer to an EquationSystem object.
static RiemannSolverSharedPtr Nektar::LaxFriedrichsSolver::create ( )
inlinestatic

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

References LaxFriedrichsSolver().

47  {
49  new LaxFriedrichsSolver());
50  }
boost::shared_ptr< RiemannSolver > RiemannSolverSharedPtr
A shared pointer to an EquationSystem object.
static RiemannSolverSharedPtr Nektar::LaxFriedrichsSolver::create ( )
inlinestatic

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

52  {
54  }
boost::shared_ptr< RiemannSolver > RiemannSolverSharedPtr
A shared pointer to an EquationSystem object.
void Nektar::LaxFriedrichsSolver::v_PointSolve ( NekDouble  rhoL,
NekDouble  rhouL,
NekDouble  rhovL,
NekDouble  rhowL,
NekDouble  EL,
NekDouble  rhoR,
NekDouble  rhouR,
NekDouble  rhovR,
NekDouble  rhowR,
NekDouble  ER,
NekDouble rhof,
NekDouble rhouf,
NekDouble rhovf,
NekDouble rhowf,
NekDouble Ef 
)
protectedvirtual

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.

References Nektar::SolverUtils::RiemannSolver::m_params, and sign.

73  {
74  static NekDouble gamma = m_params["gamma"]();
75 
76  // Left and right velocities
77  NekDouble uL = rhouL / rhoL;
78  NekDouble vL = rhovL / rhoL;
79  NekDouble wL = rhowL / rhoL;
80  NekDouble uR = rhouR / rhoR;
81  NekDouble vR = rhovR / rhoR;
82  NekDouble wR = rhowR / rhoR;
83 
84  // Left and right pressures
85  NekDouble pL = (gamma - 1.0) *
86  (EL - 0.5 * (rhouL*uL + rhovL*vL + rhowL*wL));
87 
88  NekDouble pR = (gamma - 1.0) *
89  (ER - 0.5 * (rhouR*uR + rhovR*vR + rhowR*wR));
90 
91  // Left and right speeds of sound
92  NekDouble cL = sqrt(gamma * pL / rhoL);
93  NekDouble cR = sqrt(gamma * pR / rhoR);
94 
95  // Left and right entalpies
96  NekDouble hL = (EL + pL) / rhoL;
97  NekDouble hR = (ER + pR) / rhoR;
98 
99  // Square root of rhoL and rhoR.
100  NekDouble srL = sqrt(rhoL);
101  NekDouble srR = sqrt(rhoR);
102  NekDouble srLR = srL + srR;
103 
104  // Velocity Roe averages
105  NekDouble uRoe = (srL * uL + srR * uR) / srLR;
106  NekDouble vRoe = (srL * vL + srR * vR) / srLR;
107  NekDouble wRoe = (srL * wL + srR * wR) / srLR;
108  NekDouble hRoe = (srL * hL + srR * hR) / srLR;
109  NekDouble cRoe = sqrt((gamma - 1.0)*(hRoe - 0.5 *
110  (uRoe * uRoe + vRoe * vRoe + wRoe * wRoe)));
111 
112  // Minimum and maximum wave speeds
113  NekDouble S = std::max(uRoe+cRoe, std::max(uR+cR, -uL+cL));
114  NekDouble sign = 1.0;
115 
116  if(S == -uL+cL)
117  {
118  sign = -1.0;
119  }
120 
121  // Lax-Friedrichs Riemann rho flux
122  rhof = 0.5 * ((rhouL + rhouR) - sign * S * (rhoR -rhoL));
123 
124  // Lax-Friedrichs Riemann rhou flux
125  rhouf = 0.5 * ((rhoL * uL * uL + pL + rhoR * uR * uR + pR) -
126  sign * S * (rhouR - rhouL));
127 
128  // Lax-Friedrichs Riemann rhov flux
129  rhovf = 0.5 * ((rhoL * uL * vL + rhoR * uR * vR) -
130  sign * S * (rhovR - rhovL));
131 
132  // Lax-Friedrichs Riemann rhow flux
133  rhowf = 0.5 * ((rhoL * uL * wL + rhoR * uR * wR) -
134  sign * S * (rhowR - rhowL));
135 
136  // Lax-Friedrichs Riemann E flux
137  Ef = 0.5 * ((uL * (EL + pL) + uR * (ER + pR)) -
138  sign * S * (ER - EL));
139  }
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:22
double NekDouble
std::map< std::string, RSParamFuncType > m_params
Map of parameter function types.
void Nektar::LaxFriedrichsSolver::v_PointSolve ( double  hL,
double  huL,
double  hvL,
double  hR,
double  huR,
double  hvR,
double &  hf,
double &  huf,
double &  hvf 
)
protectedvirtual

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 64 of file ShallowWaterSolver/RiemannSolvers/LaxFriedrichsSolver.cpp.

References Nektar::SolverUtils::RiemannSolver::m_params, and sign.

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 = 0.5 * ((hL * uL * uL + 0.5 * g * hL * hL +
104  hR * uR * uR + 0.5 * g * hR * hR) -
105  sign * S * (huR - huL));
106 
107  // Lax-Friedrichs Riemann hv flux
108  hvf = 0.5 * ((hL * uL * vL + hR * uR * vR) -
109  sign * S * (hvR - hvL));
110 
111  }
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:22
double NekDouble
std::map< std::string, RSParamFuncType > m_params
Map of parameter function types.
void Nektar::LaxFriedrichsSolver::v_PointSolve ( NekDouble  pL,
NekDouble  rhoL,
NekDouble  uL,
NekDouble  vL,
NekDouble  wL,
NekDouble  pR,
NekDouble  rhoR,
NekDouble  uR,
NekDouble  vR,
NekDouble  wR,
NekDouble  p0,
NekDouble  rho0,
NekDouble  u0,
NekDouble  v0,
NekDouble  w0,
NekDouble pF,
NekDouble rhoF,
NekDouble uF,
NekDouble vF,
NekDouble wF 
)
protectedvirtual

Lax-Friedrichs Riemann solver.

Parameters
pLPerturbation pressure left state
rhoLPerturbation density left state
pRPerturbation pressure right state
rhoRPerturbation density right state
uLx perturbation velocity component left state
uRx perturbation velocity component right state
vLy perturbation velocity component left state
vRy perturbation velocity component right state
wLz perturbation velocity component left state
wRz perturbation velocity component right state
p0Base pressure
rho0Base density
u0Base x velocity component
v0Base y velocity component
w0Base z velocity component
pFComputed Riemann flux for perturbation pressure
rhoFComputed Riemann flux for perturbation density
uFComputed Riemann flux for x perturbation velocity component
vFComputed Riemann flux for y perturbation velocity component
wFComputed Riemann flux for z perturbation velocity component

Reimplemented from Nektar::APESolver.

Definition at line 83 of file APESolver/RiemannSolvers/LaxFriedrichsSolver.cpp.

References ASSERTL1, Nektar::SolverUtils::RiemannSolver::CheckParams(), and Nektar::SolverUtils::RiemannSolver::m_params.

88 {
89  ASSERTL1(CheckParams("Gamma"), "Gamma not defined.");
90  const NekDouble &gamma = m_params["Gamma"]();
91 
92  // Speed of sound
93  NekDouble c = sqrt(gamma * p0 / rho0);
94 
95  // max absolute eigenvalue of the jacobian of F_n1
96  NekDouble a_1_max = std::max(std::abs(u0 - c), std::abs(u0 + c));
97 
98  NekDouble pFL = rho0*uL + u0*pL/(c*c);
99  NekDouble uFL = pL/rho0 + u0*uL + v0*vL + w0*wL;
100  NekDouble vFL = 0;
101  NekDouble wFL = 0;
102 
103  NekDouble pFR = rho0*uR + u0*pR/(c*c);
104  NekDouble uFR = pR/rho0 + u0*uR + v0*vR + w0*wR;
105  NekDouble vFR = 0;
106  NekDouble wFR = 0;
107 
108  // assemble the face-normal fluxes
109  pF = 0.5*(pFL + pFR - a_1_max*(pR/(c*c) - pL/(c*c)));
110  uF = 0.5*(uFL + uFR - a_1_max*(uR - uL));
111  vF = 0.5*(vFL + vFR - a_1_max*(vR - vL));
112  wF = 0.5*(wFL + wFR - a_1_max*(wR - wL));
113 }
double NekDouble
std::map< std::string, RSParamFuncType > m_params
Map of parameter function types.
SOLVER_UTILS_EXPORT bool CheckParams(std::string name)
Determine whether a parameter has been defined in m_params.
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191
void Nektar::LaxFriedrichsSolver::v_PointSolveVisc ( NekDouble  rhoL,
NekDouble  rhouL,
NekDouble  rhovL,
NekDouble  rhowL,
NekDouble  EL,
NekDouble  EpsL,
NekDouble  rhoR,
NekDouble  rhouR,
NekDouble  rhovR,
NekDouble  rhowR,
NekDouble  ER,
NekDouble  EpsR,
NekDouble rhof,
NekDouble rhouf,
NekDouble rhovf,
NekDouble rhowf,
NekDouble Ef,
NekDouble Epsf 
)
protectedvirtual

Reimplemented from Nektar::CompressibleSolver.

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

References Nektar::SolverUtils::RiemannSolver::m_params, and sign.

145  {
146  static NekDouble gamma = m_params["gamma"]();
147 
148  // Left and right velocities
149  NekDouble uL = rhouL / rhoL;
150  NekDouble vL = rhovL / rhoL;
151  NekDouble wL = rhowL / rhoL;
152  NekDouble uR = rhouR / rhoR;
153  NekDouble vR = rhovR / rhoR;
154  NekDouble wR = rhowR / rhoR;
155 
156  // Left and right pressures
157  NekDouble pL = (gamma - 1.0) *
158  (EL - 0.5 * (rhouL*uL + rhovL*vL + rhowL*wL));
159 
160  NekDouble pR = (gamma - 1.0) *
161  (ER - 0.5 * (rhouR*uR + rhovR*vR + rhowR*wR));
162 
163  // Left and right speeds of sound
164  NekDouble cL = sqrt(gamma * pL / rhoL);
165  NekDouble cR = sqrt(gamma * pR / rhoR);
166 
167  // Left and right entalpies
168  NekDouble hL = (EL + pL) / rhoL;
169  NekDouble hR = (ER + pR) / rhoR;
170 
171  // Square root of rhoL and rhoR.
172  NekDouble srL = sqrt(rhoL);
173  NekDouble srR = sqrt(rhoR);
174  NekDouble srLR = srL + srR;
175 
176  // Velocity Roe averages
177  NekDouble uRoe = (srL * uL + srR * uR) / srLR;
178  NekDouble vRoe = (srL * vL + srR * vR) / srLR;
179  NekDouble wRoe = (srL * wL + srR * wR) / srLR;
180  NekDouble hRoe = (srL * hL + srR * hR) / srLR;
181  NekDouble cRoe = sqrt((gamma - 1.0)*(hRoe - 0.5 *
182  (uRoe * uRoe + vRoe * vRoe + wRoe * wRoe)));
183 
184  // Minimum and maximum wave speeds
185  NekDouble S = std::max(uRoe+cRoe, std::max(uR+cR, -uL+cL));
186  NekDouble sign = 1.0;
187 
188  if(S == -uL+cL)
189  {
190  sign = -1.0;
191  }
192 
193  // Lax-Friedrichs Riemann rho flux
194  rhof = 0.5 * ((rhouL + rhouR) - sign * S * (rhoR -rhoL));
195 
196  // Lax-Friedrichs Riemann rhou flux
197  rhouf = 0.5 * ((rhoL * uL * uL + pL + rhoR * uR * uR + pR) -
198  sign * S * (rhouR - rhouL));
199 
200  // Lax-Friedrichs Riemann rhov flux
201  rhovf = 0.5 * ((rhoL * uL * vL + rhoR * uR * vR) -
202  sign * S * (rhovR - rhovL));
203 
204  // Lax-Friedrichs Riemann rhow flux
205  rhowf = 0.5 * ((rhoL * uL * wL + rhoR * uR * wR) -
206  sign * S * (rhowR - rhowL));
207 
208  // Lax-Friedrichs Riemann E flux
209  Ef = 0.5 * ((uL * (EL + pL) + uR * (ER + pR)) -
210  sign * S * (ER - EL));
211 
212  Epsf = 0.5 * ((EpsL + EpsR) -
213  sign * S * (EpsR - EpsL));
214  }
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:22
double NekDouble
std::map< std::string, RSParamFuncType > m_params
Map of parameter function types.

Member Data Documentation

static std::string Nektar::LaxFriedrichsSolver::solverName
static
Initial value:
=
RegisterCreatorFunction("LaxFriedrichs", LaxFriedrichsSolver::create,
"Lax-Friedrichs Solver")

Definition at line 56 of file APESolver/RiemannSolvers/LaxFriedrichsSolver.h.