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{
41 "LaxFriedrichs", LaxFriedrichsSolver::create,
42 "LaxFriedrichs Riemann solver");
43
46 : NonlinearSWESolver(pSession)
47{
48}
49
50/**
51 * @brief Lax-Friedrichs Riemann solver
52 *
53 * @param hL Water depth left state.
54 * @param hR Water depth right state.
55 * @param huL x-momentum component left state.
56 * @param huR x-momentum component right state.
57 * @param hvL y-momentum component left state.
58 * @param hvR y-momentum component right state.
59 * @param hf Computed Riemann flux for density.
60 * @param huf Computed Riemann flux for x-momentum component
61 * @param hvf Computed Riemann flux for y-momentum component
62 */
63void LaxFriedrichsSolver::v_PointSolve(NekDouble hL, NekDouble huL,
64 NekDouble hvL, NekDouble hR,
65 NekDouble huR, NekDouble hvR,
66 NekDouble &hf, NekDouble &huf,
67 NekDouble &hvf)
68{
69 static NekDouble g = m_params["gravity"]();
70
71 // Left and right velocities
72 NekDouble uL = huL / hL;
73 NekDouble vL = hvL / hL;
74 NekDouble uR = huR / hR;
75 NekDouble vR = hvR / hR;
76
77 // Left and right wave speeds
78 NekDouble cL = sqrt(g * hL);
79 NekDouble cR = sqrt(g * hR);
80
81 // Square root of hL and hR.
82 NekDouble srL = sqrt(hL);
83 NekDouble srR = sqrt(hR);
84 NekDouble srLR = srL + srR;
85
86 // Velocity Roe averages
87 NekDouble uRoe = (srL * uL + srR * uR) / srLR;
88 NekDouble cRoe = sqrt(0.5 * (cL * cL + cR * cR));
89
90 // Minimum and maximum wave speeds
91 NekDouble S = std::max(uRoe + cRoe, std::max(uR + cR, -(uL - cL)));
92 NekDouble sign = 1.0;
93
94 if (S == -(uL - cL))
95 {
96 sign = -1.0;
97 }
98
99 // Lax-Friedrichs Riemann h flux
100 hf = 0.5 * ((huL + huR) - sign * S * (hR - hL));
101
102 // Lax-Friedrichs Riemann hu flux
103 huf =
104 0.5 *
105 ((hL * uL * uL + 0.5 * g * hL * hL + hR * uR * uR + 0.5 * g * hR * hR) -
106 sign * S * (huR - huL));
107
108 // Lax-Friedrichs Riemann hv flux
109 hvf = 0.5 * ((hL * uL * vL + hR * uR * vR) - sign * S * (hvR - hvL));
110}
111} // 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.
Definition: NekFactory.hpp:197
std::shared_ptr< SessionReader > SessionReaderSharedPtr
RiemannSolverFactory & GetRiemannSolverFactory()
double NekDouble
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294