Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
solvers/APESolver/RiemannSolvers/UpwindSolver.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File: UpwindSolver.cpp
4 //
5 // For more information, please see: http://www.nektar.info
6 //
7 // The MIT License
8 //
9 // Copyright (c) 2014 Kilian Lackhove
10 // Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
11 // Department of Aeronautics, Imperial College London (UK), and Scientific
12 // Computing and Imaging Institute, University of Utah (USA).
13 //
14 // License for the specific language governing rights and limitations under
15 // Permission is hereby granted, free of charge, to any person obtaining a
16 // copy of this software and associated documentation files (the "Software"),
17 // to deal in the Software without restriction, including without limitation
18 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
19 // and/or sell copies of the Software, and to permit persons to whom the
20 // Software is furnished to do so, subject to the following conditions:
21 //
22 // The above copyright notice and this permission notice shall be included
23 // in all copies or substantial portions of the Software.
24 //
25 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
26 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
27 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
28 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
29 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
30 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
31 // DEALINGS IN THE SOFTWARE.
32 //
33 // Description: Upwind Riemann solver for the APE equations.
34 //
35 ///////////////////////////////////////////////////////////////////////////////
36 
38 
39 using namespace std;
40 
41 namespace Nektar
42 {
43 
44 std::string UpwindSolver::solverName = SolverUtils::GetRiemannSolverFactory().
45  RegisterCreatorFunction("APEUpwind", UpwindSolver::create, "Upwind solver for the APE equation");
46 
47 UpwindSolver::UpwindSolver() :
48  APESolver()
49 {
50 
51 }
52 
53 
54 /**
55  * @brief Upwind Riemann solver
56  *
57  * @param pL Perturbation pressure left state.
58  * @param pR Perturbation pressure right state.
59  * @param uL x perturbation verlocity component left state.
60  * @param uR x perturbation verlocity component right state.
61  * @param vL y perturbation verlocity component left state.
62  * @param vR y perturbation verlocity component right state.
63  * @param wL z perturbation verlocity component left state.
64  * @param wR z perturbation verlocity component right state.
65  * @param p0 Base pressure.
66  * @param u0 Base x verlocity component
67  * @param v0 Base y verlocity component
68  * @param w0 Base z verlocity component
69  * @param pF Computed Riemann flux for perturbation pressure.
70  * @param uF Computed Riemann flux for x perturbation verlocity component
71  * @param vF Computed Riemann flux for y perturbation verlocity component
72  * @param wF Computed Riemann flux for z perturbation verlocity component
73  */
75  NekDouble pL, NekDouble uL, NekDouble vL, NekDouble wL,
76  NekDouble pR, NekDouble uR, NekDouble vR, NekDouble wR,
77  NekDouble p0, NekDouble u0, NekDouble v0, NekDouble w0,
78  NekDouble &pF, NekDouble &uF, NekDouble &vF, NekDouble &wF)
79 {
80  // fetch params
81  ASSERTL1(CheckParams("Gamma"), "Gamma not defined.");
82  ASSERTL1(CheckParams("Rho"), "Rho not defined.");
83  const NekDouble &gamma= m_params["Gamma"]();
84  const NekDouble &rho = m_params["Rho"]();
85 
86  Array<OneD, NekDouble> characteristic(4);
87  Array<OneD, NekDouble> W(2);
88  Array<OneD, NekDouble> lambda(2);
89 
90  // compute the wave speeds
91  lambda[0]=u0 + sqrt(p0*gamma*rho)/rho;
92  lambda[1]=u0 - sqrt(p0*gamma*rho)/rho;
93 
94  // calculate the caracteristic variables
95  //left characteristics
96  characteristic[0] = pL/2 + uL*sqrt(p0*gamma*rho)/2;
97  characteristic[1] = pL/2 - uL*sqrt(p0*gamma*rho)/2;
98  //right characteristics
99  characteristic[2] = pR/2 + uR*sqrt(p0*gamma*rho)/2;
100  characteristic[3] = pR/2 - uR*sqrt(p0*gamma*rho)/2;
101 
102  //take left or right value of characteristic variable
103  for (int j = 0; j < 2; j++)
104  {
105  if (lambda[j]>=0)
106  {
107  W[j]=characteristic[j];
108  }
109  if(lambda[j]<0)
110  {
111  W[j]=characteristic[j+2];
112  }
113  }
114 
115  //calculate conservative variables from characteristics
116  NekDouble p = W[0]+W[1];
117  NekDouble u = (W[0]-W[1])/sqrt(p0*gamma*rho);
118  // do not divide by zero
119  if (p0*gamma*rho == 0)
120  {
121  u = 0.0;
122  }
123 
124  // assemble the fluxes
125  pF = gamma*p0 * u + u0*p;
126  uF = p/rho + u0*u + v0*vL + w0*wL;
127  vF = 0.0;
128  wF = 0.0;
129 }
130 
131 }
132