PDF version (NAG web site
, 64-bit version, 64-bit version)
NAG Toolbox: nag_pde_1d_parab_euler_osher (d03pv)
Purpose
nag_pde_1d_parab_euler_osher (d03pv) calculates a numerical flux function using Osher'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_osher (d03pv) calculates a numerical flux function at a single spatial point using Osher's Approximate Riemann Solver (see
Hemker and Spekreijse (1986) and
Pennington and Berzins (1994)) 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_osher (d03pv).
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 Osher 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. Osher's solver carries out an integration along a path in the phase space of
consisting of subpaths which are piecewise parallel to the eigenvectors of the Jacobian of the PDE system. There are two variants of the Osher solver termed O (original) and P (physical), which differ in the order in which the subpaths are taken. The P-variant is generally more efficient, but in some rare cases may fail (see
Hemker and Spekreijse (1986) for details). The argument
path specifies which variant is to be used. The algorithm for Osher's solver for the Euler equations is given in detail in the Appendix of
Pennington and Berzins (1994).
References
Hemker P W and Spekreijse S P (1986) Multiple grid and Osher's scheme for the efficient solution of the steady Euler equations Applied Numerical Mathematics 2 475–493
Pennington S V and Berzins M (1994) New NAG Library software for first-order partial differential equations ACM Trans. Math. Softw. 20 63–99
Quirk J J (1994) A contribution to the great Riemann solver debate Internat. J. Numer. Methods Fluids 18 555–574
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:
.
- 4:
– string (length ≥ 1)
-
The variant of the Osher scheme.
- Original.
- Physical.
Constraint:
or .
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, | , |
or | or . |
-
-
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_osher (d03pv) performs an exact calculation of the Osher numerical flux function, and so the result will be accurate to
machine precision.
Further Comments
nag_pde_1d_parab_euler_osher (d03pv) 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 Osher's scheme, in common with all Riemann solvers, may be unsuitable for some problems (see
Quirk (1994) for examples). The time taken depends on the input argument
path and on the left and right solution values, since inclusion of each subpath depends on the signs of the eigenvalues. In general this cannot be determined in advance.
Example
See
Example in
nag_pde_1d_parab_convdiff_dae (d03pl).
Open in the MATLAB editor:
d03pv_example
function d03pv_example
fprintf('d03pv example results\n\n');
uleft = [1 0 2.5];
uright = [2 1 0.5];
gamma = 1.9;
path = 'P';
[flux, ifail] = d03pv(uleft, uright, gamma, path);
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 and path = ''%s''\n',gamma,path);
fprintf('\n flux = [%7.4f %7.4f %7.4f]\n',flux);
d03pv example results
Given uleft = [ 1.0 0.0 2.5]
and uright = [ 2.0 1.0 0.5]
with gamma = 1.9000 and path = 'P'
flux = [ 0.5540 1.3956 1.9266]
PDF version (NAG web site
, 64-bit version, 64-bit version)
© The Numerical Algorithms Group Ltd, Oxford, UK. 2009–2015