Nektar++
CompressibleFlowSolver/RiemannSolvers/LaxFriedrichsSolver.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File: LaxFriedrichsSolver.cpp
4 //
5 // For more information, please see: http://www.nektar.info
6 //
7 // The MIT License
8 //
9 // Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10 // Department of Aeronautics, Imperial College London (UK), and Scientific
11 // Computing and Imaging Institute, University of Utah (USA).
12 //
13 // License for the specific language governing rights and limitations under
14 // Permission is hereby granted, free of charge, to any person obtaining a
15 // copy of this software and associated documentation files (the "Software"),
16 // to deal in the Software without restriction, including without limitation
17 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
18 // and/or sell copies of the Software, and to permit persons to whom the
19 // Software is furnished to do so, subject to the following conditions:
20 //
21 // The above copyright notice and this permission notice shall be included
22 // in all copies or substantial portions of the Software.
23 //
24 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30 // DEALINGS IN THE SOFTWARE.
31 //
32 // Description: Lax-Friedrichs Riemann solver.
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
37 
38 namespace Nektar
39 {
42  "LaxFriedrichs", LaxFriedrichsSolver::create,
43  "Lax-Friedrichs Riemann solver");
44 
45  LaxFriedrichsSolver::LaxFriedrichsSolver() : CompressibleSolver()
46  {
47 
48  }
49 
50  /**
51  * @brief Lax-Friedrichs Riemann solver
52  *
53  * @param rhoL Density left state.
54  * @param rhoR Density right state.
55  * @param rhouL x-momentum component left state.
56  * @param rhouR x-momentum component right state.
57  * @param rhovL y-momentum component left state.
58  * @param rhovR y-momentum component right state.
59  * @param rhowL z-momentum component left state.
60  * @param rhowR z-momentum component right state.
61  * @param EL Energy left state.
62  * @param ER Energy right state.
63  * @param rhof Computed Riemann flux for density.
64  * @param rhouf Computed Riemann flux for x-momentum component
65  * @param rhovf Computed Riemann flux for y-momentum component
66  * @param rhowf Computed Riemann flux for z-momentum component
67  * @param Ef Computed Riemann flux for energy.
68  */
70  NekDouble rhoL, NekDouble rhouL, NekDouble rhovL, NekDouble rhowL, NekDouble EL,
71  NekDouble rhoR, NekDouble rhouR, NekDouble rhovR, NekDouble rhowR, NekDouble ER,
72  NekDouble &rhof, NekDouble &rhouf, NekDouble &rhovf, NekDouble &rhowf, NekDouble &Ef)
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  }
140 
142  NekDouble rhoL, NekDouble rhouL, NekDouble rhovL, NekDouble rhowL, NekDouble EL, NekDouble EpsL,
143  NekDouble rhoR, NekDouble rhouR, NekDouble rhovR, NekDouble rhowR, NekDouble ER, NekDouble EpsR,
144  NekDouble &rhof, NekDouble &rhouf, NekDouble &rhovf, NekDouble &rhowf, NekDouble &Ef, NekDouble &Epsf)
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  }
215 
216 }
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:22
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)
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.
RiemannSolverFactory & GetRiemannSolverFactory()
double NekDouble
std::map< std::string, RSParamFuncType > m_params
Map of parameter function types.
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, tDescription pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:215