NAG Library Routine Document

d03puf  (dim1_parab_euler_roe)

 Contents

    1  Purpose
    7  Accuracy
    10  Example

1
Purpose

d03puf 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 d03pff, d03plf or d03psf, but may also be applicable to other conservative upwind schemes requiring numerical flux functions.

2
Specification

Fortran Interface
Subroutine d03puf ( uleft, uright, gamma, flux, ifail)
Integer, Intent (Inout):: ifail
Real (Kind=nag_wp), Intent (In):: uleft(3), uright(3), gamma
Real (Kind=nag_wp), Intent (Out):: flux(3)
C Header Interface
#include nagmk26.h
void  d03puf_ ( const double uleft[], const double uright[], const double *gamma, double flux[], Integer *ifail)

3
Description

d03puf 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 routines d03pff, d03plf and d03psf, the left and right solution values are derived automatically from the solution values at adjacent spatial points and supplied to the subroutine argument numflx from which you may call d03puf.
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 routine 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:     uleft3 – Real (Kind=nag_wp) arrayInput
On entry: ulefti must contain the left value of the component Ui, for i=1,2,3. That is, uleft1 must contain the left value of ρ, uleft2 must contain the left value of m and uleft3 must contain the left value of e.
Constraints:
  • uleft10.0;
  • Left pressure, pl0.0, where pl is calculated using (3).
2:     uright3 – Real (Kind=nag_wp) arrayInput
On entry: urighti must contain the right value of the component Ui, for i=1,2,3. That is, uright1 must contain the right value of ρ, uright2 must contain the right value of m and uright3 must contain the right value of e.
Constraints:
  • uright10.0;
  • Right pressure, pr0.0, where pr is calculated using (3).
3:     gamma – Real (Kind=nag_wp)Input
On entry: the ratio of specific heats, γ.
Constraint: gamma>0.0.
4:     flux3 – Real (Kind=nag_wp) arrayOutput
On exit: fluxi contains the numerical flux component F^i, for i=1,2,3.
5:     ifail – IntegerInput/Output
On entry: ifail must be set to 0, -1​ or ​1. If you are unfamiliar with this argument you should refer to Section 3.4 in How to Use the NAG Library and its Documentation for details.
For environments where it might be inappropriate to halt program execution when an error is detected, the value -1​ or ​1 is recommended. If the output of error messages is undesirable, then the value 1 is recommended. Otherwise, if you are not familiar with this argument, the recommended value is 0. When the value -1​ or ​1 is used it is essential to test the value of ifail on exit.
On exit: ifail=0 unless the routine detects an error or a warning has been flagged (see Section 6).
Note: if the left and/or right values of ρ or p (from (3)) are found to be negative, then the routine will terminate with an error exit (ifail=2). If the routine is being called from the numflx etc., then a soft fail option (ifail=1 or -1) is recommended so that a recalculation of the current time step can be forced using the numflx argument ires (see d03pff or d03plf).

6
Error Indicators and Warnings

If on entry ifail=0 or -1, explanatory error messages are output on the current error message unit (as defined by x04aaf).
Errors or warnings detected by the routine:
ifail=1
On entry,gamma0.0.
ifail=2
On entry,the left and/or right density or pressure value is less than 0.0.
ifail=-99
An unexpected error has been triggered by this routine. Please contact NAG.
See Section 3.9 in How to Use the NAG Library and its Documentation for further information.
ifail=-399
Your licence key may have expired or may not have been installed correctly.
See Section 3.8 in How to Use the NAG Library and its Documentation for further information.
ifail=-999
Dynamic memory allocation failed.
See Section 3.7 in How to Use the NAG Library and its Documentation for further information.

7
Accuracy

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

8
Parallelism and Performance

d03puf is not thread safe and should not be called from a multithreaded user program. Please see Section 3.12.1 in How to Use the NAG Library and its Documentation for more information on thread safety.
d03puf is not threaded in any implementation.

9
Further Comments

d03puf must only be used to calculate the numerical flux for the Euler equations in exactly the form given by (2), with ulefti and urighti 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 routine does not detect or correct any entropy violation. The time taken is independent of the input arguments.

10
Example

See Section 10 in d03plf.
© The Numerical Algorithms Group Ltd, Oxford, UK. 2017