Nektar++
ShallowWaterSolver/RiemannSolvers/LaxFriedrichsSolver.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: LaxFriedrichsSolver.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: Lax-Friedrichs Riemann solver.
32//
33///////////////////////////////////////////////////////////////////////////////
34
36
37namespace Nektar
38{
39
42 "LaxFriedrichs", LaxFriedrichsSolver::create,
43 "LaxFriedrichs Riemann solver");
44
47 : NonlinearSWESolver(pSession)
48{
49}
50
51/**
52 * @brief Lax-Friedrichs Riemann solver
53 *
54 * @param hL Water depth left state.
55 * @param hR Water depth right state.
56 * @param huL x-momentum component left state.
57 * @param huR x-momentum component right state.
58 * @param hvL y-momentum component left state.
59 * @param hvR y-momentum component right state.
60 * @param hf Computed Riemann flux for density.
61 * @param huf Computed Riemann flux for x-momentum component
62 * @param hvf Computed Riemann flux for y-momentum component
63 */
64void LaxFriedrichsSolver::v_PointSolve(NekDouble hL, NekDouble huL,
65 NekDouble hvL, NekDouble hR,
66 NekDouble huR, NekDouble hvR,
67 NekDouble &hf, NekDouble &huf,
68 NekDouble &hvf)
69{
70 static NekDouble g = m_params["gravity"]();
71
72 // Left and right velocities
73 NekDouble uL = huL / hL;
74 NekDouble vL = hvL / hL;
75 NekDouble uR = huR / hR;
76 NekDouble vR = hvR / hR;
77
78 // Left and right wave speeds
79 NekDouble cL = sqrt(g * hL);
80 NekDouble cR = sqrt(g * hR);
81
82 // Square root of hL and hR.
83 NekDouble srL = sqrt(hL);
84 NekDouble srR = sqrt(hR);
85 NekDouble srLR = srL + srR;
86
87 // Velocity Roe averages
88 NekDouble uRoe = (srL * uL + srR * uR) / srLR;
89 NekDouble cRoe = sqrt(0.5 * (cL * cL + cR * cR));
90
91 // Minimum and maximum wave speeds
92 NekDouble S = std::max(uRoe + cRoe, std::max(uR + cR, -(uL - cL)));
93 NekDouble sign = 1.0;
94
95 if (S == -(uL - cL))
96 {
97 sign = -1.0;
98 }
99
100 // Lax-Friedrichs Riemann h flux
101 hf = 0.5 * ((huL + huR) - sign * S * (hR - hL));
102
103 // Lax-Friedrichs Riemann hu flux
104 huf =
105 0.5 *
106 ((hL * uL * uL + 0.5 * g * hL * hL + hR * uR * uR + 0.5 * g * hR * hR) -
107 sign * S * (huR - huL));
108
109 // Lax-Friedrichs Riemann hv flux
110 hvf = 0.5 * ((hL * uL * vL + hR * uR * vR) - sign * S * (hvR - hvL));
111}
112
113} // namespace Nektar
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:47
static RiemannSolverSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession)
LaxFriedrichsSolver(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