Nektar++
ShallowWaterSolver/RiemannSolvers/HLLSolver.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: HLLSolver.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 Nonlinear Shallow Water Equations.
32//
33///////////////////////////////////////////////////////////////////////////////
34
36
37namespace Nektar
38{
39
40std::string HLLSolver::solverName =
42 "HLL", HLLSolver::create, "HLL Riemann solver");
43
45 : NonlinearSWESolver(pSession)
46{
47}
48
49/**
50 * @brief HLL Riemann solver for the Nonlinear Shallow Water Equations
51 *
52 * @param hL Water depth left state.
53 * @param hR Water depth right state.
54 * @param huL x-momentum component left state.
55 * @param huR x-momentum component right state.
56 * @param hvL y-momentum component left state.
57 * @param hvR y-momentum component right state.
58 * @param hf Computed Riemann flux for density.
59 * @param huf Computed Riemann flux for x-momentum component
60 * @param hvf Computed Riemann flux for y-momentum component
61 */
62void HLLSolver::v_PointSolve(NekDouble hL, NekDouble huL, NekDouble hvL,
63 NekDouble hR, NekDouble huR, NekDouble hvR,
64 NekDouble &hf, NekDouble &huf, NekDouble &hvf)
65{
66 static NekDouble g = m_params["gravity"]();
67
68 // Left and Right velocities
69 NekDouble uL = huL / hL;
70 NekDouble vL = hvL / hL;
71 NekDouble uR = huR / hR;
72 NekDouble vR = hvR / hR;
73
74 // Left and right wave speeds
75 NekDouble cL = sqrt(g * hL);
76 NekDouble cR = sqrt(g * hR);
77
78 // the two-rarefaction wave assumption
79 NekDouble hstar, fL, fR;
80 hstar = 0.5 * (cL + cR) + 0.25 * (uL - uR);
81 hstar *= hstar;
82 hstar *= (1.0 / g);
83
84 // Compute SL
85 NekDouble SL;
86 if (hstar > hL)
87 {
88 SL = uL - cL * sqrt(0.5 * ((hstar * hstar + hstar * hL) / (hL * hL)));
89 }
90 else
91 {
92 SL = uL - cL;
93 }
94
95 // Compute SR
96 NekDouble SR;
97 if (hstar > hR)
98 {
99 SR = uR + cR * sqrt(0.5 * ((hstar * hstar + hstar * hR) / (hR * hR)));
100 }
101 else
102 {
103 SR = uR + cR;
104 }
105
106 if (SL >= 0)
107 {
108 hf = hL * uL;
109 huf = uL * uL * hL + 0.5 * g * hL * hL;
110 hvf = hL * uL * vL;
111 }
112 else if (SR <= 0)
113 {
114 hf = hR * uR;
115 huf = uR * uR * hR + 0.5 * g * hR * hR;
116 hvf = hR * uR * vR;
117 }
118 else
119 {
120 hf = (SR * hL * uL - SL * hR * uR + SL * SR * (hR - hL)) / (SR - SL);
121 fL = uL * uL * hL + 0.5 * g * hL * hL;
122 fR = uR * uR * hR + 0.5 * g * hR * hR;
123 huf = (SR * fL - SL * fR + SL * SR * (hR * uR - hL * uL)) / (SR - SL);
124 fL = uL * vL * hL;
125 fR = uR * vR * hR;
126 hvf = (SR * fL - SL * fR + SL * SR * (hR * vR - hL * vL)) / (SR - SL);
127 }
128}
129
130} // namespace Nektar
HLLSolver(const LibUtilities::SessionReaderSharedPtr &pSession)
static RiemannSolverSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession)
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
RiemannSolverFactory & GetRiemannSolverFactory()
double NekDouble
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:285