Nektar++
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 // Permission is hereby granted, free of charge, to any person obtaining a
14 // copy of this software and associated documentation files (the "Software"),
15 // to deal in the Software without restriction, including without limitation
16 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
17 // and/or sell copies of the Software, and to permit persons to whom the
18 // Software is furnished to do so, subject to the following conditions:
19 //
20 // The above copyright notice and this permission notice shall be included
21 // in all copies or substantial portions of the Software.
22 //
23 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29 // DEALINGS IN THE SOFTWARE.
30 //
31 // Description: HLL Riemann solver for the Linear Shallow Water Equations.
32 // Only valid for constant depth
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
37 
38 namespace Nektar
39 {
40  std::string LinearHLLSolver::solverName =
42  "LinearHLL",
44  "Linear HLL Riemann solver");
45 
48  : LinearSWESolver(pSession)
49  {
50 
51  }
52 
53  /**
54  * @brief HLL Riemann solver for the linear Shallow Water Equations
55  *
56  * @param etaL Free surface elevation left state.
57  * @param etaR Free surface elevation right state.
58  * @param uL x-velocity left state.
59  * @param uR x-velocity right state.
60  * @param vL y-velocity left state.
61  * @param vR y-velocity right state.
62  * @param dL still water depth component left state.
63  * @param dR still water depth component right state.
64  * @param etaf Computed Riemann flux for continuity.
65  * @param uf Computed Riemann flux for x-momentum component
66  * @param vf Computed Riemann flux for y-momentum component
67  */
69  NekDouble etaL, NekDouble uL, NekDouble vL, NekDouble dL,
70  NekDouble etaR, NekDouble uR, NekDouble vR, NekDouble dR,
71  NekDouble &etaf, NekDouble &uf, NekDouble &vf)
72  {
73  static NekDouble g = m_params["gravity"]();
74 
75  ASSERTL0(dL == dR, "not constant depth in LinearHLL");
76 
77  // Left and right wave speeds
78  NekDouble cL = sqrt(g * dL);
79  NekDouble cR = sqrt(g * dR);
80 
81  NekDouble SL = uL - cL;
82  NekDouble SR = uR + cR;
83 
84  if (SL >= 0)
85  {
86  etaf = dL * uL;
87  uf = g * etaL;
88  vf = 0.0;
89  }
90  else if (SR <= 0)
91  {
92  etaf = dR * uR;
93  uf = g * etaR;
94  vf = 0.0;
95  }
96  else
97  {
98  etaf = (SR * dL * uL - SL * dR * uR +
99  SL * SR * (etaR - etaL)) / (SR - SL);
100  uf = (SR * g * etaL - SL * g * etaR +
101  SL * SR * (uR - uL)) / (SR - SL);
102  vf = (SL * SR * (vR - vL)) / (SR - SL);
103  }
104  }
105 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:200
virtual void v_PointSolve(NekDouble etaL, NekDouble uL, NekDouble vL, NekDouble dL, NekDouble etaR, NekDouble uR, NekDouble vR, NekDouble dR, NekDouble &etaf, NekDouble &uf, NekDouble &vf)
HLL Riemann solver for the linear Shallow Water Equations.
static std::string solverName
LinearHLLSolver(const LibUtilities::SessionReaderSharedPtr &pSession)
static RiemannSolverSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession)
std::map< std::string, RSParamFuncType > m_params
Map of parameter function types.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
RiemannSolverFactory & GetRiemannSolverFactory()
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
double NekDouble
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:267