D03PXF (PDF version)
D03 Chapter Contents
D03 Chapter Introduction
NAG Library Manual

NAG Library Routine Document

D03PXF

Note:  before using this routine, please read the Users' Note for your implementation to check the interpretation of bold italicised terms and other implementation-dependent details.

+ Contents

    1  Purpose
    7  Accuracy

1  Purpose

D03PXF calculates a numerical flux function using an Exact 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

SUBROUTINE D03PXF ( ULEFT, URIGHT, GAMMA, TOL, NITER, FLUX, IFAIL)
INTEGER  NITER, IFAIL
REAL (KIND=nag_wp)  ULEFT(3), URIGHT(3), GAMMA, TOL, FLUX(3)

3  Description

D03PXF calculates a numerical flux function at a single spatial point using an Exact Riemann Solver (see Toro (1996) and Toro (1989)) 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 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 D03PXF.
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-m22 ρ meρ+mργ-1 e-m22ρ , (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 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.
The algorithm is termed an Exact Riemann Solver although it does in fact calculate an approximate solution to a true Riemann problem, as opposed to an Approximate Riemann Solver which involves some form of alternative modelling of the Riemann problem. The approximation part of the Exact Riemann Solver is a Newton–Raphson iterative procedure to calculate the pressure, and you must supply a tolerance TOL and a maximum number of iterations NITER. Default values for these parameters can be chosen.
A solution cannot be found by this routine if there is a vacuum state in the Riemann problem (loosely characterised by zero density), or if such a state is generated by the interaction of two non-vacuum data states. In this case a Riemann solver which can handle vacuum states has to be used (see Toro (1996)).

4  References

Toro E F (1989) A weighted average flux method for hyperbolic conservation laws Proc. Roy. Soc. Lond. A423 401–418
Toro E F (1996) Riemann Solvers and Upwind Methods for Fluid Dynamics Springer–Verlag

5  Parameters

1:     ULEFT(3) – 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.
2:     URIGHT(3) – 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.
3:     GAMMA – REAL (KIND=nag_wp)Input
On entry: the ratio of specific heats, γ.
Constraint: GAMMA>0.0.
4:     TOL – REAL (KIND=nag_wp)Input
On entry: the tolerance to be used in the Newton–Raphson procedure to calculate the pressure. If TOL is set to zero then the default value of 1.0×10-6 is used.
Constraint: TOL0.0.
5:     NITER – INTEGERInput
On entry: the maximum number of Newton–Raphson iterations allowed. If NITER is set to zero then the default value of 20 is used.
Constraint: NITER0.
6:     FLUX(3) – REAL (KIND=nag_wp) arrayOutput
On exit: FLUXi contains the numerical flux component F^i, for i=1,2,3.
7:     IFAIL – INTEGERInput/Output
On entry: IFAIL must be set to 0, -1​ or ​1. If you are unfamiliar with this parameter you should refer to Section 3.3 in the Essential Introduction 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 parameter, 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 parameter 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,
orTOL<0.0,
orNITER<0.
IFAIL=2
On entry,the left and/or right density or derived pressure value is less than 0.0.
IFAIL=3
A vacuum condition has been detected therefore a solution cannot be found using this routine. You are advised to check your problem formulation.
IFAIL=4
The internal Newton–Raphson iterative procedure used to solve for the pressure has failed to converge. The value of TOL or NITER may be too small, but if the problem persists try an Approximate Riemann Solver (D03PUF, D03PVF or D03PWF).

7  Accuracy

The algorithm is exact apart from the calculation of the pressure which uses a Newton–Raphson iterative procedure, the accuracy of which is controlled by the parameter TOL. In some cases the initial guess for the Newton–Raphson procedure is exact and no further iterations are required.

8  Further Comments

D03PXF 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.
For some problems the routine may fail or be highly inefficient in comparison with an Approximate Riemann Solver (e.g., D03PUF, D03PVF or D03PWF). Hence it is advisable to try more than one Riemann solver and to compare the performance and the results.
The time taken by the routine is independent of all input parameters other than TOL.

9  Example

This example uses D03PLF and D03PXF to solve the Euler equations in the domain 0x1 for 0<t0.035 with initial conditions for the primitive variables ρx,t, ux,t and px,t given by
ρx,0=5.99924, ux,0=-19.5975, px,0=460.894,   for ​x<0.5, ρx,0=5.99242, ux,0=-6.19633, px,0=046.095,   for ​x>0.5.
This test problem is taken from Toro (1996) and its solution represents the collision of two strong shocks travelling in opposite directions, consisting of a left facing shock (travelling slowly to the right), a right travelling contact discontinuity and a right travelling shock wave. There is an exact solution to this problem (see Toro (1996)) but the calculation is lengthy and has therefore been omitted.

9.1  Program Text

Program Text (d03pxfe.f90)

9.2  Program Data

Program Data (d03pxfe.d)

9.3  Program Results

Program Results (d03pxfe.r)

Produced by GNUPLOT 4.4 patchlevel 0 0 0.005 0.01 0.015 0.02 0.025 0.03 0.035 0 0.2 0.4 0.6 0.8 1 5 10 15 20 25 30 35 Example Program 1 Euler Equation Solution Showing Collision of Two Strong Shocks DENSITY Time x
Produced by GNUPLOT 4.4 patchlevel 0 0 0.005 0.01 0.015 0.02 0.025 0.03 0.035 0 0.2 0.4 0.6 0.8 1 -10 -5 0 5 10 15 20 Euler Equation Solution Showing Collision of Two Strong Shocks VELOCITY Time x
Produced by GNUPLOT 4.4 patchlevel 0 0 0.005 0.01 0.015 0.02 0.025 0.03 0.035 0 0.2 0.4 0.6 0.8 1 0 200 400 600 800 1000 1200 1400 1600 1800 Euler Equation Solution Showing Collision of Two Strong Shocks PRESSURE Time x

D03PXF (PDF version)
D03 Chapter Contents
D03 Chapter Introduction
NAG Library Manual

© The Numerical Algorithms Group Ltd, Oxford, UK. 2012