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,d03plford03psf, but may also be applicable to other conservative upwind schemes requiring numerical flux functions.
The routine may be called by the names d03puf or nagf_pde_dim1_parab_euler_roe.
3Description
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,d03plfandd03psf, 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:
(1)
with
(2)
where is the density, is the momentum, is the specific total energy, and is the (constant) ratio of specific heats. The pressure is given by
(3)
where is the velocity.
The routine calculates the Roe approximation to the numerical flux function , where and are the left and right solution values, and is the intermediate state arising from the similarity solution of the Riemann problem defined by
(4)
with and as in (2), and initial piecewise constant values for and for . The spatial domain is , where 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).
4References
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 Fluids18 555–574
Roe P L (1981) Approximate Riemann solvers, parameter vectors, and difference schemes J. Comput. Phys.43 357–372
5Arguments
1: – Real (Kind=nag_wp) arrayInput
On entry: must contain the left value of the component , for . That is, must contain the left value of , must contain the left value of and must contain the left value of .
On entry: must contain the right value of the component , for . That is, must contain the right value of , must contain the right value of and must contain the right value of .
On exit: contains the numerical flux component , for .
5: – IntegerInput/Output
On entry: ifail must be set to , or to set behaviour on detection of an error; these values have no effect when no error is detected.
A value of causes the printing of an error message and program execution will be halted; otherwise program execution continues. A value of means that an error message is printed while a value of means that it is not.
If halting is not appropriate, the value or is recommended. If message printing is undesirable, then the value is recommended. Otherwise, the value is recommended. When the value or is used it is essential to test the value of ifail on exit.
On exit: 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 (from (3)) are found to be negative, then the routine will terminate with an error exit (). If the routine is being called from the numflx etc., then a soft fail option ( or ) is recommended so that a recalculation of the current time step can be forced using the numflx argument ires (see d03pfford03plf).
6Error Indicators and Warnings
If on entry or , explanatory error messages are output on the current error message unit (as defined by x04aaf).
Errors or warnings detected by the routine:
On entry, .
Constraint: .
Left pressure value : .
On entry, .
Constraint: .
On entry, .
Constraint: .
Right pressure value : .
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.
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.
Dynamic memory allocation failed.
See Section 9 in the Introduction to the NAG Library FL Interface for further information.
7Accuracy
d03puf performs an exact calculation of the Roe numerical flux function, and so the result will be accurate to machine precision.
8Parallelism and Performance
Background information to multithreading can be found in the Multithreading documentation.
d03puf 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.
d03puf is not threaded in any implementation.
9Further Comments
d03puf must only be used to calculate the numerical flux for the Euler equations in exactly the form given by (2), with and containing the left and right values of and , for , 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.