NAG CL Interface
d03puc (dim1_parab_euler_roe)
1
Purpose
d03puc 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
d03pfc,
d03plc or
d03psc, but may also be applicable to other conservative upwind schemes requiring numerical flux functions.
2
Specification
The function may be called by the names: d03puc, nag_pde_dim1_parab_euler_roe or nag_pde_parab_1d_euler_roe.
3
Description
d03puc 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 functions
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
d03puc.
The Euler equations for a perfect gas in conservative form are:
with
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
where
is the velocity.
The function 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
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).
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:
– const double
Input
-
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 .
Constraints:
- ;
- Left pressure, , where is calculated using (3).
-
2:
– const double
Input
-
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 .
Constraints:
- ;
- Right pressure, , where is calculated using (3).
-
3:
– double
Input
-
On entry: the ratio of specific heats, .
Constraint:
.
-
4:
– double
Output
-
On exit: contains the numerical flux component , for .
-
5:
– Nag_D03_Save *
Communication Structure
-
saved may contain data concerning the computation required by
d03puc as passed through to
numflx from one of the integrator functions
d03pfc,
d03plc or
d03psc. You should not change the components of
saved.
-
6:
– 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 had an illegal value.
- 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_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 : .
On entry, .
Constraint: .
On entry, .
Constraint: .
On entry, .
Constraint: .
Right pressure value : .
7
Accuracy
d03puc performs an exact calculation of the Roe numerical flux function, and so the result will be accurate to machine precision.
8
Parallelism and Performance
Background information to multithreading can be found in the
Multithreading documentation.
d03puc is not threaded in any implementation.
d03puc 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 function does not detect or correct any entropy violation. The time taken is independent of the input arguments.
10
Example