naginterfaces.library.pde.dim1_parab_euler_roe¶
- naginterfaces.library.pde.dim1_parab_euler_roe(uleft, uright, gamma, comm)[source]¶
dim1_parab_euler_roe
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 schemesdim1_parab_convdiff()
,dim1_parab_convdiff_dae()
ordim1_parab_convdiff_remesh()
, but may also be applicable to other conservative upwind schemes requiring numerical flux functions.For full information please refer to the NAG Library document for d03pu
https://support.nag.com/numeric/nl/nagdoc_30.3/flhtml/d03/d03puf.html
- Parameters
- uleftfloat, array-like, shape
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 .
- urightfloat, array-like, shape
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 .
- gammafloat
The ratio of specific heats, .
- commdict, communication object
Communication structure.
On initial entry: need not be set.
- Returns
- fluxfloat, ndarray, shape
contains the numerical flux component , for .
- Raises
- NagValueError
- (errno )
On entry, .
Constraint: .
- (errno )
Right pressure value : .
- (errno )
Left pressure value : .
- (errno )
On entry, .
Constraint: .
- (errno )
On entry, .
Constraint: .
- Notes
dim1_parab_euler_roe
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
dim1_parab_convdiff()
,dim1_parab_convdiff_dae()
anddim1_parab_convdiff_remesh()
, the left and right solution values are derived automatically from the solution values at adjacent spatial points and supplied to the function argument from which you may calldim1_parab_euler_roe
.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