NAG CL Interface
d03pxc (dim1_​parab_​euler_​exact)

Settings help

CL Name Style:


1 Purpose

d03pxc 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 d03pfc, d03plc or d03psc, but may also be applicable to other conservative upwind schemes requiring numerical flux functions.

2 Specification

#include <nag.h>
void  d03pxc (const double uleft[], const double uright[], double gamma, double tol, Integer niter, double flux[], Nag_D03_Save *saved, NagError *fail)
The function may be called by the names: d03pxc, nag_pde_dim1_parab_euler_exact or nag_pde_parab_1d_euler_exact.

3 Description

d03pxc 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 d03pfc, d03plc and 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 d03pxc.
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 function calculates the numerical flux function F(UL,UR)=F(U*(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 U(y,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 function 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: uleft[3] const double Input
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.
2: uright[3] const double Input
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.
3: gamma double Input
On entry: the ratio of specific heats, γ.
Constraint: gamma>0.0.
4: tol double 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: flux[3] double Output
On exit: flux[i-1] contains the numerical flux component F^i, for i=1,2,3.
7: saved Nag_D03_Save * Communication Structure
saved may contain data concerning the computation required by d03pxc as passed through to numflx from one of the integrator functions d03pfc, d03plc or d03psc. You should not change the components of saved.
8: fail NagError * Input/Output
The NAG error argument (see Section 7 in the Introduction to the NAG Library CL Interface).

6 Error Indicators and Warnings

NE_ALLOC_FAIL
Dynamic memory allocation failed.
See Section 3.1.2 in the Introduction to the NAG Library CL Interface for further information.
NE_BAD_PARAM
On entry, argument value had an illegal value.
NE_INT
On entry, niter=value.
Constraint: niter0.
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 7.5 in the Introduction to the NAG Library CL Interface for further information.
NE_ITER_FAIL_CONV
Newton–Raphson iteration failed to converge.
NE_NO_LICENCE
Your licence key may have expired or may not have been installed correctly.
See Section 8 in the Introduction to the NAG Library CL Interface for further information.
NE_REAL
Left pressure value pl<0.0: pl=value.
On entry, gamma=value.
Constraint: gamma>0.0.
On entry, tol=value.
Constraint: tol0.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.
NE_VACUUM
A vacuum condition has been detected.

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

Background information to multithreading can be found in the Multithreading documentation.
d03pxc is not threaded in any implementation.

9 Further Comments

d03pxc 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.
For some problems the function may fail or be highly inefficient in comparison with an Approximate Riemann Solver (e.g., d03puc, d03pvc or d03pwc). Hence it is advisable to try more than one Riemann solver and to compare the performance and the results.
The time taken by the function is independent of all input arguments other than tol.

10 Example

This example uses d03plc and d03pxc to solve the Euler equations in the domain 0x1 for 0<t0.035 with initial conditions for the primitive variables ρ(x,t), u(x,t) and p(x,t) given by
ρ(x,0)=5.99924, u(x,0)=-19.5975, p(x,0)=460.894,   for ​x<0.5, ρ(x,0)=5.99242, u(x,0)=-6.19633, p(x,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 (d03pxce.c)

10.2 Program Data

Program Data (d03pxce.d)

10.3 Program Results

Program Results (d03pxce.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