Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
ShallowWaterSolver/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",
44  "LaxFriedrichs Riemann solver");
45 
46  LaxFriedrichsSolver::LaxFriedrichsSolver() : NonlinearSWESolver()
47  {
48 
49  }
50 
51  /**
52  * @brief Lax-Friedrichs Riemann solver
53  *
54  * @param hL Water depth left state.
55  * @param hR Water depth right state.
56  * @param huL x-momentum component left state.
57  * @param huR x-momentum component right state.
58  * @param hvL y-momentum component left state.
59  * @param hvR y-momentum component right state.
60  * @param hf Computed Riemann flux for density.
61  * @param huf Computed Riemann flux for x-momentum component
62  * @param hvf Computed Riemann flux for y-momentum component
63  */
65  double hL, double huL, double hvL,
66  double hR, double huR, double hvR,
67  double &hf, double &huf, double &hvf)
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  }
112 }