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
38namespace Nektar
39{
42 "LinearHLL", LinearHLLSolver::create, "Linear HLL Riemann solver");
43
46 : LinearSWESolver(pSession)
47{
48}
49
50/**
51 * @brief HLL Riemann solver for the linear Shallow Water Equations
52 *
53 * @param etaL Free surface elevation left state.
54 * @param etaR Free surface elevation right state.
55 * @param uL x-velocity left state.
56 * @param uR x-velocity right state.
57 * @param vL y-velocity left state.
58 * @param vR y-velocity right state.
59 * @param dL still water depth component left state.
60 * @param dR still water depth component right state.
61 * @param etaf Computed Riemann flux for continuity.
62 * @param uf Computed Riemann flux for x-momentum component
63 * @param vf Computed Riemann flux for y-momentum component
64 */
66 NekDouble dL, NekDouble etaR, NekDouble uR,
67 NekDouble vR, NekDouble dR, NekDouble &etaf,
68 NekDouble &uf, NekDouble &vf)
69{
70 static NekDouble g = m_params["gravity"]();
71
72 ASSERTL0(dL == dR, "not constant depth in LinearHLL");
73
74 // Left and right wave speeds
75 NekDouble cL = sqrt(g * dL);
76 NekDouble cR = sqrt(g * dR);
77
78 NekDouble SL = uL - cL;
79 NekDouble SR = uR + cR;
80
81 if (SL >= 0)
82 {
83 etaf = dL * uL;
84 uf = g * etaL;
85 vf = 0.0;
86 }
87 else if (SR <= 0)
88 {
89 etaf = dR * uR;
90 uf = g * etaR;
91 vf = 0.0;
92 }
93 else
94 {
95 etaf =
96 (SR * dL * uL - SL * dR * uR + SL * SR * (etaR - etaL)) / (SR - SL);
97 uf = (SR * g * etaL - SL * g * etaR + SL * SR * (uR - uL)) / (SR - SL);
98 vf = (SL * SR * (vR - vL)) / (SR - SL);
99 }
100}
101} // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
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) override
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()
double NekDouble
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294