NAG Library Function Document
nag_pde_parab_1d_euler_hll (d03pwc)
1 Purpose
nag_pde_parab_1d_euler_hll (d03pwc) calculates a numerical flux function using a modified HLL (Harten–Lax–van Leer) Approximate Riemann Solver for the Euler equations in conservative form. It is designed primarily for use with the upwind discretization schemes
nag_pde_parab_1d_cd (d03pfc),
nag_pde_parab_1d_cd_ode (d03plc) or
nag_pde_parab_1d_cd_ode_remesh (d03psc), but may also be applicable to other conservative upwind schemes requiring numerical flux functions.
2 Specification
#include <nag.h> |
#include <nagd03.h> |
|
3 Description
nag_pde_parab_1d_euler_hll (d03pwc) calculates a numerical flux function at a single spatial point using a modified HLL (Harten–Lax–van Leer) Approximate Riemann Solver (see
Toro (1992),
Toro (1996) and
Toro et al. (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
nag_pde_parab_1d_cd (d03pfc),
nag_pde_parab_1d_cd_ode (d03plc) and
nag_pde_parab_1d_cd_ode_remesh (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 nag_pde_parab_1d_euler_hll (d03pwc).
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 an 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.
4 References
Toro E F (1992) The weighted average flux method applied to the Euler equations Phil. Trans. R. Soc. Lond. A341 499–530
Toro E F (1996) Riemann Solvers and Upwind Methods for Fluid Dynamics Springer–Verlag
Toro E F, Spruce M and Spears W (1994) Restoration of the contact surface in the HLL Riemann solver J. Shock Waves 4 25–34
5 Arguments
- 1:
– const doubleInput
-
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 doubleInput
-
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:
– doubleInput
-
On entry: the ratio of specific heats, .
Constraint:
.
- 4:
– doubleOutput
-
On exit: contains the numerical flux component , for .
- 5:
– Nag_D03_Save *Communication Structure
-
saved may contain data concerning the computation required by nag_pde_parab_1d_euler_hll (d03pwc) as passed through to
numflx from one of the integrator functions
nag_pde_parab_1d_cd (d03pfc),
nag_pde_parab_1d_cd_ode (d03plc) or
nag_pde_parab_1d_cd_ode_remesh (d03psc). You should not change the components of
saved.
- 6:
– NagError *Input/Output
-
The NAG error argument (see
Section 3.6 in the Essential Introduction).
6 Error Indicators and Warnings
- NE_ALLOC_FAIL
-
Dynamic memory allocation failed.
See
Section 3.2.1.2 in the Essential Introduction 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.
An unexpected error has been triggered by this function. Please contact
NAG.
See
Section 3.6.6 in the Essential Introduction for further information.
- NE_NO_LICENCE
-
Your licence key may have expired or may not have been installed correctly.
See
Section 3.6.5 in the Essential Introduction for further information.
- NE_REAL
-
Left pressure value : .
On entry, .
Constraint: .
On entry, .
Constraint: .
On entry, .
Constraint: .
Right pressure value : .
7 Accuracy
nag_pde_parab_1d_euler_hll (d03pwc) performs an exact calculation of the HLL (Harten–Lax–van Leer) numerical flux function, and so the result will be accurate to machine precision.
8 Parallelism and Performance
Not applicable.
nag_pde_parab_1d_euler_hll (d03pwc) 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. The time taken is independent of the input arguments.
10 Example
This example uses
nag_pde_parab_1d_cd_ode (d03plc) and nag_pde_parab_1d_euler_hll (d03pwc) to solve the Euler equations in the domain
for
with initial conditions for the primitive variables
,
and
given by
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 (d03pwce.c)
10.2 Program Data
Program Data (d03pwce.d)
10.3 Program Results
Program Results (d03pwce.r)