NAG FL Interface
d03pxf (dim1_​parab_​euler_​exact)

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

Fortran Interface
Subroutine d03pxf ( uleft, uright, gamma, tol, niter, flux, ifail)
Integer, Intent (In) :: niter
Integer, Intent (Inout) :: ifail
Real (Kind=nag_wp), Intent (In) :: uleft(3), uright(3), gamma, tol
Real (Kind=nag_wp), Intent (Out) :: flux(3)
C Header Interface
#include <nag.h>
void  d03pxf_ (const double uleft[], const double uright[], const double *gamma, const double *tol, const Integer *niter, double flux[], Integer *ifail)
The routine may be called by the names d03pxf or nagf_pde_dim1_parab_euler_exact.

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 arguments 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 Arguments

1: uleft3 Real (Kind=nag_wp) array Input
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: uright3 Real (Kind=nag_wp) array Input
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 Integer Input
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: flux3 Real (Kind=nag_wp) array Output
On exit: fluxi contains the numerical flux component F^i, for i=1,2,3.
7: ifail Integer Input/Output
On entry: ifail must be set to 0, -1 or 1. If you are unfamiliar with this argument you should refer to Section 4 in the Introduction to the NAG Library FL Interface 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, gamma=value.
Constraint: gamma>0.0.
On entry, niter=value.
Constraint: niter0.
On entry, tol=value.
Constraint: tol0.0.
ifail=2
Left pressure value pl<0.0: pl=value.
On entry, uleft1=value.
Constraint: uleft10.0.
On entry, uright1=value.
Constraint: uright10.0.
Right pressure value pr<0.0: pr=value.
ifail=3
A vacuum condition has been detected.
ifail=4
Newton–Raphson iteration failed to converge.
ifail=-99
An unexpected error has been triggered by this routine. Please contact NAG.
See Section 7 in the Introduction to the NAG Library FL Interface for further information.
ifail=-399
Your licence key may have expired or may not have been installed correctly.
See Section 8 in the Introduction to the NAG Library FL Interface for further information.
ifail=-999
Dynamic memory allocation failed.
See Section 9 in the Introduction to the NAG Library FL Interface for further information.

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 argument tol. In some cases the initial guess for the Newton–Raphson procedure is exact and no further iterations are required.

8 Parallelism and Performance

d03pxf is not thread safe and should not be called from a multithreaded user program. Please see Section 1 in FL Interface Multithreading for more information on thread safety.
d03pxf is not threaded in any implementation.

9 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 arguments other than tol.

10 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.

10.1 Program Text

Program Text (d03pxfe.f90)

10.2 Program Data

Program Data (d03pxfe.d)

10.3 Program Results

Program Results (d03pxfe.r)
GnuplotProduced by GNUPLOT 5.0 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 gnuplot_plot_1 Time x
GnuplotProduced by GNUPLOT 5.0 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 gnuplot_plot_1 Time x
GnuplotProduced by GNUPLOT 5.0 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 gnuplot_plot_1 Time x