PDF version (NAG web site
, 64-bit version, 64-bit version)
NAG Toolbox: nag_pde_1d_parab_euler_roe (d03pu)
Purpose
nag_pde_1d_parab_euler_roe (d03pu) 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
nag_pde_1d_parab_convdiff (d03pf),
nag_pde_1d_parab_convdiff_dae (d03pl) or
nag_pde_1d_parab_convdiff_remesh (d03ps), but may also be applicable to other conservative upwind schemes requiring numerical flux functions.
Syntax
Description
nag_pde_1d_parab_euler_roe (d03pu) 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
nag_pde_1d_parab_convdiff (d03pf),
nag_pde_1d_parab_convdiff_dae (d03pl) and
nag_pde_1d_parab_convdiff_remesh (d03ps), 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
nag_pde_1d_parab_euler_roe (d03pu).
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).
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
Parameters
Compulsory Input Parameters
- 1:
– double array
-
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:
– double array
-
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 scalar
-
The ratio of specific heats, .
Constraint:
.
Optional Input Parameters
None.
Output Parameters
- 1:
– double array
-
contains the numerical flux component , for .
- 2:
– int64int32nag_int scalar
unless the function detects an error (see
Error Indicators and Warnings).
Error Indicators and Warnings
Errors or warnings detected by the function:
-
-
-
-
On entry, | the left and/or right density or pressure value is less than . |
-
An unexpected error has been triggered by this routine. Please
contact
NAG.
-
Your licence key may have expired or may not have been installed correctly.
-
Dynamic memory allocation failed.
Accuracy
nag_pde_1d_parab_euler_roe (d03pu) performs an exact calculation of the Roe numerical flux function, and so the result will be accurate to
machine precision.
Further Comments
nag_pde_1d_parab_euler_roe (d03pu) 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.
Example
See
Example in
nag_pde_1d_parab_convdiff_dae (d03pl).
Open in the MATLAB editor:
d03pu_example
function d03pu_example
fprintf('d03pu example results\n\n');
uleft = [1 0 2.5];
uright = [2 1 0.5];
gamma = 1.9;
[flux, ifail] = d03pu(uleft, uright, gamma);
fprintf(' Given uleft = [%7.1f %7.1f %7.1f]\n',uleft);
fprintf(' and uright = [%7.1f %7.1f %7.1f]\n',uright);
fprintf(' with gamma = %7.4f\n',gamma);
fprintf('\n flux = [%7.4f %7.4f %7.4f]\n',flux);
d03pu example results
Given uleft = [ 1.0 0.0 2.5]
and uright = [ 2.0 1.0 0.5]
with gamma = 1.9000
flux = [ 0.8548 1.3149 1.5161]
PDF version (NAG web site
, 64-bit version, 64-bit version)
© The Numerical Algorithms Group Ltd, Oxford, UK. 2009–2015