hide long namesshow long names
hide short namesshow short names
Integer type:  int32  int64  nag_int  show int32  show int32  show int64  show int64  show nag_int  show nag_int

PDF version (NAG web site, 64-bit version, 64-bit version)
Chapter Contents
Chapter Introduction
NAG Toolbox

NAG Toolbox: nag_pde_1d_parab_euler_hll (d03pw)

Purpose

nag_pde_1d_parab_euler_hll (d03pw) 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_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

[flux, ifail] = d03pw(uleft, uright, gamma)
[flux, ifail] = nag_pde_1d_parab_euler_hll(uleft, uright, gamma)

Description

nag_pde_1d_parab_euler_hll (d03pw) 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_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_hll (d03pw).
The Euler equations for a perfect gas in conservative form are:
(U)/(t) + (F)/(x) = 0,
U t + F x =0,
(1)
with
U =
[ ρ m e ]
  and  F =
[ m (m2)/ρ + (γ − 1) (e − (m2)/(2 ρ)) (me)/ρ + m/ρ(γ − 1) (e − (m2)/(2ρ)) ]
,
U=[ ρ m e ]   and  F=[ m m2ρ+(γ-1) (e-m22 ρ ) meρ+mρ(γ-1) (e-m22ρ ) ] ,
(2)
where ρρ is the density, mm is the momentum, ee is the specific total energy and γγ is the (constant) ratio of specific heats. The pressure pp is given by
p = (γ1) (e(ρu2)/2) ,
p=(γ-1) (e-ρu22) ,
(3)
where u = m / ρu=m/ρ is the velocity.
The function calculates an approximation to the numerical flux function F(UL,UR) = F(U*(UL,UR))F(UL,UR)=F(U*(UL,UR)), where U = ULU=UL and U = URU=UR are the left and right solution values, and U*(UL,UR)U*(UL,UR) is the intermediate state ω(0)ω(0) arising from the similarity solution U(y,t) = ω(y / t)U(y,t)=ω(y/t) of the Riemann problem defined by
(U)/(t) + (F)/(y) = 0,
U t + F y =0,
(4)
with UU and FF as in (2), and initial piecewise constant values U = ULU=UL for y < 0y<0 and U = URU=UR for y > 0y>0. The spatial domain is < y < -<y<, where y = 0y=0 is the point at which the numerical flux is required.

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

Parameters

Compulsory Input Parameters

1:     uleft(33) – double array
uleft(i)ulefti must contain the left value of the component UiUi, for i = 1,2,3i=1,2,3. That is, uleft(1)uleft1 must contain the left value of ρρ, uleft(2)uleft2 must contain the left value of mm and uleft(3)uleft3 must contain the left value of ee.
Constraints:
  • uleft(1)0.0uleft10.0;
  • Left pressure, pl0.0pl0.0, where plpl is calculated using (3).
2:     uright(33) – double array
uright(i)urighti must contain the right value of the component UiUi, for i = 1,2,3i=1,2,3. That is, uright(1)uright1 must contain the right value of ρρ, uright(2)uright2 must contain the right value of mm and uright(3)uright3 must contain the right value of ee.
Constraints:
  • uright(1)0.0uright10.0;
  • Right pressure, pr0.0pr0.0, where prpr is calculated using (3).
3:     gamma – double scalar
The ratio of specific heats, γγ.
Constraint: gamma > 0.0gamma>0.0.

Optional Input Parameters

None.

Input Parameters Omitted from the MATLAB Interface

None.

Output Parameters

1:     flux(33) – double array
flux(i)fluxi contains the numerical flux component iF^i, for i = 1,2,3i=1,2,3.
2:     ifail – int64int32nag_int scalar
ifail = 0ifail=0 unless the function detects an error (see [Error Indicators and Warnings]).

Error Indicators and Warnings

Errors or warnings detected by the function:
  ifail = 1ifail=1
On entry,gamma0.0gamma0.0.
  ifail = 2ifail=2
On entry,the left and/or right density or derived pressure value is less than 0.00.0.

Accuracy

nag_pde_1d_parab_euler_hll (d03pw) performs an exact calculation of the HLL (Harten–Lax–van Leer) numerical flux function, and so the result will be accurate to machine precision.

Further Comments

nag_pde_1d_parab_euler_hll (d03pw) must only be used to calculate the numerical flux for the Euler equations in exactly the form given by (2), with uleft(i)ulefti and uright(i)urighti containing the left and right values of ρ,mρ,m and ee, for i = 1,2,3i=1,2,3, respectively. The time taken is independent of the input parameters.

Example

function nag_pde_1d_parab_euler_hll_example
uleft = [5.99924;
     117.5701059;
     2304.275075187626];
uright = [5.99924;
     117.5701059;
     2304.275075187626];
gamma = 1.4;
[flux, ifail] = nag_pde_1d_parab_euler_hll(uleft, uright, gamma)
 

flux =

   1.0e+04 *

    0.0118
    0.2765
    5.4190


ifail =

                    0


function d03pw_example
uleft = [5.99924;
     117.5701059;
     2304.275075187626];
uright = [5.99924;
     117.5701059;
     2304.275075187626];
gamma = 1.4;
[flux, ifail] = d03pw(uleft, uright, gamma)
 

flux =

   1.0e+04 *

    0.0118
    0.2765
    5.4190


ifail =

                    0



PDF version (NAG web site, 64-bit version, 64-bit version)
Chapter Contents
Chapter Introduction
NAG Toolbox

© The Numerical Algorithms Group Ltd, Oxford, UK. 2009–2013