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