NAG C Library Function Document

nag_pde_parab_1d_euler_roe (d03puc)

1
Purpose

nag_pde_parab_1d_euler_roe (d03puc) calculates a numerical flux function using Roe's Approximate Riemann Solver for the Euler equations in conservative form. It is designed primarily for use with the upwind discretization schemes nag_pde_parab_1d_cd (d03pfc), nag_pde_parab_1d_cd_ode (d03plc) or nag_pde_parab_1d_cd_ode_remesh (d03psc), but may also be applicable to other conservative upwind schemes requiring numerical flux functions.

2
Specification

#include <nag.h>
#include <nagd03.h>
void  nag_pde_parab_1d_euler_roe (const double uleft[], const double uright[], double gamma, double flux[], Nag_D03_Save *saved, NagError *fail)

3
Description

nag_pde_parab_1d_euler_roe (d03puc) calculates a numerical flux function at a single spatial point using Roe's Approximate Riemann Solver (see Roe (1981)) for the Euler equations (for a perfect gas) in conservative form. You must supply the left and right solution values at the point where the numerical flux is required, i.e., the initial left and right states of the Riemann problem defined below.
In the functions nag_pde_parab_1d_cd (d03pfc), nag_pde_parab_1d_cd_ode (d03plc) and nag_pde_parab_1d_cd_ode_remesh (d03psc), the left and right solution values are derived automatically from the solution values at adjacent spatial points and supplied to the function argument numflx from which you may call nag_pde_parab_1d_euler_roe (d03puc).
The Euler equations for a perfect gas in conservative form are:
U t + F x = 0 , (1)
with
U = ρ m e   and   F = m m2ρ + γ-1 e - m2 2ρ meρ + mρ γ-1 e - m2 2ρ , (2)
where ρ is the density, m is the momentum, e is the specific total energy, and γ is the (constant) ratio of specific heats. The pressure p is given by
p=γ-1 e-ρu22 , (3)
where u=m/ρ is the velocity.
The function calculates the Roe approximation to the numerical flux function FUL,UR=FU*UL,UR, where U=UL and U=UR are the left and right solution values, and U*UL,UR is the intermediate state ω0 arising from the similarity solution Uy,t=ωy/t of the Riemann problem defined by
U t + F y =0, (4)
with U and F as in (2), and initial piecewise constant values U=UL for y<0 and U=UR for y>0. The spatial domain is -<y<, where y=0 is the point at which the numerical flux is required. This implementation of Roe's scheme for the Euler equations uses the so-called argument-vector method described in Roe (1981).

4
References

LeVeque R J (1990) Numerical Methods for Conservation Laws Birkhäuser Verlag
Quirk J J (1994) A contribution to the great Riemann solver debate Internat. J. Numer. Methods Fluids 18 555–574
Roe P L (1981) Approximate Riemann solvers, parameter vectors, and difference schemes J. Comput. Phys. 43 357–372

5
Arguments

1:     uleft[3] const doubleInput
On entry: uleft[i-1] must contain the left value of the component Ui, for i=1,2,3. That is, uleft[0] must contain the left value of ρ, uleft[1] must contain the left value of m and uleft[2] must contain the left value of e.
Constraints:
  • uleft[0]0.0;
  • Left pressure, pl0.0, where pl is calculated using (3).
2:     uright[3] const doubleInput
On entry: uright[i-1] must contain the right value of the component Ui, for i=1,2,3. That is, uright[0] must contain the right value of ρ, uright[1] must contain the right value of m and uright[2] must contain the right value of e.
Constraints:
  • uright[0]0.0;
  • Right pressure, pr0.0, where pr is calculated using (3).
3:     gamma doubleInput
On entry: the ratio of specific heats, γ.
Constraint: gamma>0.0.
4:     flux[3] doubleOutput
On exit: flux[i-1] contains the numerical flux component F^i, for i=1,2,3.
5:     saved Nag_D03_Save *Communication Structure
saved may contain data concerning the computation required by nag_pde_parab_1d_euler_roe (d03puc) as passed through to numflx from one of the integrator functions nag_pde_parab_1d_cd (d03pfc), nag_pde_parab_1d_cd_ode (d03plc) or nag_pde_parab_1d_cd_ode_remesh (d03psc). You should not change the components of saved.
6:     fail NagError *Input/Output
The NAG error argument (see Section 3.7 in How to Use the NAG Library and its Documentation).

6
Error Indicators and Warnings

NE_ALLOC_FAIL
Dynamic memory allocation failed.
See Section 2.3.1.2 in How to Use the NAG Library and its Documentation for further information.
NE_BAD_PARAM
On entry, argument value had an illegal value.
NE_INTERNAL_ERROR
An internal error has occurred in this function. Check the function call and any array sizes. If the call is correct then please contact NAG for assistance.
See Section 2.7.6 in How to Use the NAG Library and its Documentation for further information.
NE_NO_LICENCE
Your licence key may have expired or may not have been installed correctly.
See Section 2.7.5 in How to Use the NAG Library and its Documentation for further information.
NE_REAL
Left pressure value pl<0.0: pl=value.
On entry, gamma=value.
Constraint: gamma>0.0.
On entry, uleft[0]=value.
Constraint: uleft[0]0.0.
On entry, uright[0]=value.
Constraint: uright[0]0.0.
Right pressure value pr<0.0: pr=value.

7
Accuracy

nag_pde_parab_1d_euler_roe (d03puc) performs an exact calculation of the Roe numerical flux function, and so the result will be accurate to machine precision.

8
Parallelism and Performance

nag_pde_parab_1d_euler_roe (d03puc) is not threaded in any implementation.

9
Further Comments

nag_pde_parab_1d_euler_roe (d03puc) must only be used to calculate the numerical flux for the Euler equations in exactly the form given by (2), with uleft[i-1] and uright[i-1] containing the left and right values of ρ,m and e, for i=1,2,3, respectively. It should be noted that Roe's scheme, in common with all Riemann solvers, may be unsuitable for some problems (see Quirk (1994) for examples). In particular Roe's scheme does not satisfy an ‘entropy condition’ which guarantees that the approximate solution of the PDE converges to the correct physical solution, and hence it may admit non-physical solutions such as expansion shocks. The algorithm used in this function does not detect or correct any entropy violation. The time taken is independent of the input arguments.

10
Example

See Section 10 in nag_pde_parab_1d_cd_ode (d03plc).