Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
LinearHLLSolver.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File: LinearHLLSolver.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: HLL Riemann solver for the Linear Shallow Water Equations.
33 // Only valid for constant depth
34 //
35 ///////////////////////////////////////////////////////////////////////////////
36 
38 
39 namespace Nektar
40 {
41  std::string LinearHLLSolver::solverName =
43  "LinearHLL",
45  "Linear HLL Riemann solver");
46 
48  {
49 
50  }
51 
52  /**
53  * @brief HLL Riemann solver for the linear Shallow Water Equations
54  *
55  * @param etaL Free surface elevation left state.
56  * @param etaR Free surface elevation right state.
57  * @param uL x-velocity left state.
58  * @param uR x-velocity right state.
59  * @param vL y-velocity left state.
60  * @param vR y-velocity right state.
61  * @param dL still water depth component left state.
62  * @param dR still water depth component right state.
63  * @param etaf Computed Riemann flux for continuity.
64  * @param uf Computed Riemann flux for x-momentum component
65  * @param vf Computed Riemann flux for y-momentum component
66  */
68  double etaL, double uL, double vL, double dL,
69  double etaR, double uR, double vR, double dR,
70  double &etaf, double &uf, double &vf)
71  {
72  static NekDouble g = m_params["gravity"]();
73 
74  if (dL != dR)
75  {
76  ASSERTL0(false,"not constant depth in LinearHLL");
77  }
78 
79  // Left and right wave speeds
80  NekDouble cL = sqrt(g * dL);
81  NekDouble cR = sqrt(g * dR);
82 
83  NekDouble SL = uL - cL;
84  NekDouble SR = uR + cR;
85 
86  if (SL >= 0)
87  {
88  etaf = dL * uL;
89  uf = g * etaL;
90  vf = 0.0;
91  }
92  else if (SR <= 0)
93  {
94  etaf = dR * uR;
95  uf = g * etaR;
96  vf = 0.0;
97  }
98  else
99  {
100  etaf = (SR * dL * uL - SL * dR * uR +
101  SL * SR * (etaR - etaL)) / (SR - SL);
102  uf = (SR * g * etaL - SL * g * etaR +
103  SL * SR * (uR - uL)) / (SR - SL);
104  vf = (SL * SR * (vR - vL)) / (SR - SL);
105  }
106  }
107 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
static std::string solverName
static RiemannSolverSharedPtr create()
RiemannSolverFactory & GetRiemannSolverFactory()
double NekDouble
virtual void v_PointSolve(double etaL, double uL, double vL, double dL, double etaR, double uR, double vR, double dR, double &etaf, double &uf, double &vf)
HLL Riemann solver for the linear Shallow Water Equations.
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