Nektar++
ExactSolverToro.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File: ExactSolver.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: Exact Riemann solver (Toro 2009).
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
36 
37 #define TOL 1e-6
38 #define NRITER 20
39 
40 namespace Nektar
41 {
42  std::string ExactSolverToro::solverName =
44  "ExactToro",
46  "Exact Riemann solver");
47 
50  : CompressibleSolver(pSession)
51  {
52 
53  }
54 
55  /**
56  * @brief Use either PVRS, two-rarefaction or two-shock Riemann solvers to
57  * calculate an initial pressure for the Newton-Raphson scheme.
58  *
59  * @param g Array of calculated gamma values.
60  * @param rhoL Density left state.
61  * @param rhoR Density right state.
62  * @param uL x-velocity component left state.
63  * @param uR x-velocity component right state.
64  * @param pL Pressure component left state.
65  * @param pR Pressure component right state.
66  * @param cL Sound speed component left state.
67  * @param cR Sound speed component right state.
68  * @return Computed initial guess for the Newton-Raphson scheme.
69  */
70  inline NekDouble guessp(
71  NekDouble g[],
72  NekDouble rhoL, NekDouble uL, NekDouble pL, NekDouble cL,
73  NekDouble rhoR, NekDouble uR, NekDouble pR, NekDouble cR)
74  {
75  const NekDouble quser = 2.0;
76  NekDouble cup, ppv, pmin, pmax, qmax;
77 
78  cup = 0.25*(rhoL + rhoR)*(cL + cR);
79  ppv = 0.5 *(pL + pR) + 0.5*(uL - uR)*cup;
80  ppv = std::max(0.0, ppv);
81  pmin = std::min(pL, pR);
82  pmax = std::max(pL, pR);
83  qmax = pmax/pmin;
84 
85  if (qmax <= quser && pmin <= ppv && ppv <= pmax)
86  {
87  // Select PVRS Riemann solver.
88  return ppv;
89  }
90  else if (ppv < pmin)
91  {
92  // Select two-rarefaction Riemann solver.
93  NekDouble pq = pow(pL/pR, g[1]);
94  NekDouble um = (pq*uL/cL + uR/cR + g[4]*(pq - 1.0))/(pq/cL + 1.0/cR);
95  NekDouble ptL = 1.0 + g[7]*(uL - um)/cL;
96  NekDouble ptR = 1.0 + g[7]*(um - uR)/cR;
97  return 0.5*(pL*pow(ptL, g[3]) + pR*pow(ptR, g[3]));
98  }
99  else
100  {
101  // Select two-shock Riemann solver with PVRS as estimate.
102  NekDouble geL = sqrt((g[5]/rhoL)/(g[6]*pL + ppv));
103  NekDouble geR = sqrt((g[5]/rhoR)/(g[6]*pR + ppv));
104  return (geL*pL + geR*pR - (uR - uL))/(geL + geR);
105  }
106  }
107 
108  /**
109  * @brief Evaluate pressure functions fL and fR in Newton iteration of
110  * Riemann solver (see equation 4.85 and 4.86 of Toro 2009).
111  *
112  * @param g Array of gamma parameters.
113  * @param p Pressure at current iteration.
114  * @param dk Density (left or right state)
115  * @param pk Pressure (left or right state)
116  * @param ck Sound speed (left or right state)
117  * @param f Computed pressure (shock).
118  * @param fd Computed pressure (rarefaction).
119  */
120  inline void prefun(NekDouble *g, NekDouble p, NekDouble dk, NekDouble pk,
121  NekDouble ck, NekDouble &f, NekDouble &fd)
122  {
123  if (p <= pk)
124  {
125  // Rarefaction wave
126  NekDouble prat = p/pk;
127  f = g[4]*ck*(pow(prat, g[1])-1.0);
128  fd = pow(prat, -g[2])/(dk*ck);
129  }
130  else
131  {
132  // Shock wave
133  NekDouble ak = g[5]/dk;
134  NekDouble bk = g[6]*pk;
135  NekDouble qrt = sqrt(ak/(bk+p));
136  f = (p-pk)*qrt;
137  fd = (1.0-0.5*(p-pk)/(bk+p))*qrt;
138  }
139  }
140 
141  /**
142  * @brief Exact Riemann solver for the Euler equations.
143  *
144  * This algorithm is transcribed from:
145  *
146  * "Riemann Solvers and Numerical Methods for Fluid Dynamics: A Practical
147  * Introduction", E. F. Toro (3rd edition, 2009).
148  *
149  * The full Fortran 77 routine can be found at the end of chapter 4 (section
150  * 4.9). This transcription is essentially the functions STARPU and SAMPLE
151  * glued together, and variable names are kept mostly the same. See the
152  * preceding chapter which explains the derivation of the solver. The
153  * routines PREFUN and GUESSP are kept separate and are reproduced above.
154  *
155  * @param rhoL Density left state.
156  * @param rhoR Density right state.
157  * @param rhouL x-momentum component left state.
158  * @param rhouR x-momentum component right state.
159  * @param rhovL y-momentum component left state.
160  * @param rhovR y-momentum component right state.
161  * @param rhowL z-momentum component left state.
162  * @param rhowR z-momentum component right state.
163  * @param EL Energy left state.
164  * @param ER Energy right state.
165  * @param rhof Computed Riemann flux for density.
166  * @param rhouf Computed Riemann flux for x-momentum component
167  * @param rhovf Computed Riemann flux for y-momentum component
168  * @param rhowf Computed Riemann flux for z-momentum component
169  * @param Ef Computed Riemann flux for energy.
170  */
172  NekDouble rhoL, NekDouble rhouL, NekDouble rhovL, NekDouble rhowL, NekDouble EL,
173  NekDouble rhoR, NekDouble rhouR, NekDouble rhovR, NekDouble rhowR, NekDouble ER,
174  NekDouble &rhof, NekDouble &rhouf, NekDouble &rhovf, NekDouble &rhowf, NekDouble &Ef)
175  {
176  static NekDouble gamma = m_params["gamma"]();
177 
178  // Left and right variables.
179  NekDouble uL = rhouL / rhoL;
180  NekDouble vL = rhovL / rhoL;
181  NekDouble wL = rhowL / rhoL;
182  NekDouble uR = rhouR / rhoR;
183  NekDouble vR = rhovR / rhoR;
184  NekDouble wR = rhowR / rhoR;
185 
186  // Left and right pressure.
187  NekDouble pL = (gamma - 1.0) *
188  (EL - 0.5 * (rhouL*uL + rhovL*vL + rhowL*wL));
189  NekDouble pR = (gamma - 1.0) *
190  (ER - 0.5 * (rhouR*uR + rhovR*vR + rhowR*wR));
191 
192  // Compute gammas.
193  NekDouble g[] = { gamma,
194  (gamma-1.0)/(2.0*gamma),
195  (gamma+1.0)/(2.0*gamma),
196  2.0*gamma/(gamma-1.0),
197  2.0/(gamma-1.0),
198  2.0/(gamma+1.0),
199  (gamma-1.0)/(gamma+1.0),
200  0.5*(gamma-1.0),
201  gamma-1.0 };
202 
203  // Compute sound speeds.
204  NekDouble cL = sqrt(gamma * pL / rhoL);
205  NekDouble cR = sqrt(gamma * pR / rhoR);
206 
207  ASSERTL0(g[4]*(cL+cR) > (uR-uL), "Vacuum is generated by given data.");
208 
209  // Guess initial pressure.
210  NekDouble pOld = guessp(g, rhoL, uL, pL, cL, rhoR, uR, pR, cR);
211  NekDouble uDiff = uR - uL;
212  NekDouble p, fL, fR, fLd, fRd, change;
213  int k;
214 
215  // Newton-Raphson iteration for pressure in star region.
216  for (k = 0; k < NRITER; ++k)
217  {
218  prefun(g, pOld, rhoL, pL, cL, fL, fLd);
219  prefun(g, pOld, rhoR, pR, cR, fR, fRd);
220  p = pOld - (fL+fR+uDiff) / (fLd+fRd);
221  change = 2 * fabs((p-pOld)/(p+pOld));
222 
223  if (change <= TOL)
224  {
225  break;
226  }
227 
228  if (p < 0.0)
229  {
230  p = TOL;
231  }
232 
233  pOld = p;
234  }
235 
236  ASSERTL0(k < NRITER, "Divergence in Newton-Raphson scheme");
237 
238  // Compute velocity in star region.
239  NekDouble u = 0.5*(uL+uR+fR-fL);
240 
241  // -- SAMPLE ROUTINE --
242  // The variable S aligns with the Fortran parameter S of the SAMPLE
243  // routine, but is hard-coded as 0 (and should be optimised out by the
244  // compiler). Since we are using a Godunov scheme we pick this as 0 (see
245  // chapter 6 of Toro 2009).
246  const NekDouble S = 0.0;
247 
248  // Computed primitive variables.
249  NekDouble outRho, outU, outV, outW, outP;
250 
251  if (S <= u)
252  {
253  if (p <= pL)
254  {
255  // Left rarefaction
256  NekDouble shL = uL - cL;
257  if (S <= shL)
258  {
259  // Sampled point is left data state.
260  outRho = rhoL;
261  outU = uL;
262  outV = vL;
263  outW = wL;
264  outP = pL;
265  }
266  else
267  {
268  NekDouble cmL = cL*pow(p/pL, g[1]);
269  NekDouble stL = u - cmL;
270 
271  if (S > stL)
272  {
273  // Sampled point is star left state
274  outRho = rhoL*pow(p/pL, 1.0/gamma);
275  outU = u;
276  outV = vL;
277  outW = wL;
278  outP = p;
279  }
280  else
281  {
282  // Sampled point is inside left fan
283  NekDouble c = g[5]*(cL + g[7]*(uL - S));
284  outRho = rhoL*pow(c/cL, g[4]);
285  outU = g[5]*(cL + g[7]*uL + S);
286  outV = vL;
287  outW = wL;
288  outP = pL*pow(c/cL, g[3]);
289  }
290  }
291  }
292  else
293  {
294  // Left shock
295  NekDouble pmL = p / pL;
296  NekDouble SL = uL - cL*sqrt(g[2]*pmL + g[1]);
297  if (S <= SL)
298  {
299  // Sampled point is left data state
300  outRho = rhoL;
301  outU = uL;
302  outV = vL;
303  outW = wL;
304  outP = pL;
305  }
306  else
307  {
308  // Sampled point is star left state
309  outRho = rhoL*(pmL + g[6])/(pmL*g[6] + 1.0);
310  outU = u;
311  outV = vL;
312  outW = wL;
313  outP = p;
314  }
315  }
316  }
317  else
318  {
319  if (p > pR)
320  {
321  // Right shock
322  NekDouble pmR = p/pR;
323  NekDouble SR = uR + cR*sqrt(g[2]*pmR + g[1]);
324  if (S >= SR)
325  {
326  // Sampled point is right data state
327  outRho = rhoR;
328  outU = uR;
329  outV = vR;
330  outW = wR;
331  outP = pR;
332  }
333  else
334  {
335  // Sampled point is star right state
336  outRho = rhoR*(pmR + g[6])/(pmR*g[6] + 1.0);
337  outU = u;
338  outV = vR;
339  outW = wR;
340  outP = p;
341  }
342  }
343  else
344  {
345  // Right rarefaction
346  NekDouble shR = uR + cR;
347 
348  if (S >= shR)
349  {
350  // Sampled point is right data state
351  outRho = rhoR;
352  outU = uR;
353  outV = vR;
354  outW = wR;
355  outP = pR;
356  }
357  else
358  {
359  NekDouble cmR = cR*pow(p/pR, g[1]);
360  NekDouble stR = u + cmR;
361 
362  if (S <= stR)
363  {
364  // Sampled point is star right state
365  outRho = rhoR*pow(p/pR, 1.0/gamma);
366  outU = u;
367  outV = vR;
368  outW = wR;
369  outP = p;
370  }
371  else
372  {
373  // Sampled point is inside left fan
374  NekDouble c = g[5]*(cR - g[7]*(uR-S));
375  outRho = rhoR*pow(c/cR, g[4]);
376  outU = g[5]*(-cR + g[7]*uR + S);
377  outV = vR;
378  outW = wR;
379  outP = pR*pow(c/cR, g[3]);
380  }
381  }
382  }
383  }
384 
385  // Transform computed primitive variables to fluxes.
386  rhof = outRho * outU;
387  rhouf = outP + outRho*outU*outU;
388  rhovf = outRho * outU * outV;
389  rhowf = outRho * outU * outW;
390  Ef = outU*(outP/(gamma-1.0) + 0.5*outRho*
391  (outU*outU + outV*outV + outW*outW) + outP);
392  }
393 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
virtual void v_PointSolve(NekDouble rhoL, NekDouble rhouL, NekDouble rhovL, NekDouble rhowL, NekDouble EL, NekDouble rhoR, NekDouble rhouR, NekDouble rhovR, NekDouble rhowR, NekDouble ER, NekDouble &rhof, NekDouble &rhouf, NekDouble &rhovf, NekDouble &rhowf, NekDouble &Ef)
Exact Riemann solver for the Euler equations.
#define TOL
RiemannSolverFactory & GetRiemannSolverFactory()
double NekDouble
NekDouble guessp(NekDouble g[], NekDouble rhoL, NekDouble uL, NekDouble pL, NekDouble cL, NekDouble rhoR, NekDouble uR, NekDouble pR, NekDouble cR)
Use either PVRS, two-rarefaction or two-shock Riemann solvers to calculate an initial pressure for th...
ExactSolverToro(const LibUtilities::SessionReaderSharedPtr &pSession)
#define NRITER
static RiemannSolverSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession)
std::map< std::string, RSParamFuncType > m_params
Map of parameter function types.
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:199
std::shared_ptr< SessionReader > SessionReaderSharedPtr
void prefun(NekDouble *g, NekDouble p, NekDouble dk, NekDouble pk, NekDouble ck, NekDouble &f, NekDouble &fd)
Evaluate pressure functions fL and fR in Newton iteration of Riemann solver (see equation 4...
static std::string solverName