NAG FL Interface
d03pff (dim1_parab_convdiff)
1
Purpose
d03pff integrates a system of linear or nonlinear convection-diffusion equations in one space dimension, with optional source terms. The system must be posed in conservative form. Convection terms are discretized using a sophisticated upwind scheme involving a user-supplied numerical flux function based on the solution of a Riemann problem at each mesh point. The method of lines is employed to reduce the PDEs to a system of ordinary differential equations (ODEs), and the resulting system is solved using a backward differentiation formula (BDF) method.
2
Specification
Fortran Interface
Subroutine d03pff ( |
npde, ts, tout, pdedef, numflx, bndary, u, npts, x, acc, tsmax, rsave, lrsave, isave, lisave, itask, itrace, ind, ifail) |
Integer, Intent (In) |
:: |
npde, npts, lrsave, lisave, itask, itrace |
Integer, Intent (Inout) |
:: |
isave(lisave), ind, ifail |
Real (Kind=nag_wp), Intent (In) |
:: |
tout, x(npts), acc(2), tsmax |
Real (Kind=nag_wp), Intent (Inout) |
:: |
ts, u(npde,npts), rsave(lrsave) |
External |
:: |
pdedef, numflx, bndary |
|
C Header Interface
#include <nag.h>
void |
d03pff_ (const Integer *npde, double *ts, const double *tout, void (NAG_CALL *pdedef)(const Integer *npde, const double *t, const double *x, const double u[], const double ux[], double p[], double c[], double d[], double s[], Integer *ires), void (NAG_CALL *numflx)(const Integer *npde, const double *t, const double *x, const double uleft[], const double uright[], double flux[], Integer *ires), void (NAG_CALL *bndary)(const Integer *npde, const Integer *npts, const double *t, const double x[], const double u[], const Integer *ibnd, double g[], Integer *ires), double u[], const Integer *npts, const double x[], const double acc[], const double *tsmax, double rsave[], const Integer *lrsave, Integer isave[], const Integer *lisave, const Integer *itask, const Integer *itrace, Integer *ind, Integer *ifail) |
|
C++ Header Interface
#include <nag.h> extern "C" {
void |
d03pff_ (const Integer &npde, double &ts, const double &tout, void (NAG_CALL *pdedef)(const Integer &npde, const double &t, const double &x, const double u[], const double ux[], double p[], double c[], double d[], double s[], Integer &ires), void (NAG_CALL *numflx)(const Integer &npde, const double &t, const double &x, const double uleft[], const double uright[], double flux[], Integer &ires), void (NAG_CALL *bndary)(const Integer &npde, const Integer &npts, const double &t, const double x[], const double u[], const Integer &ibnd, double g[], Integer &ires), double u[], const Integer &npts, const double x[], const double acc[], const double &tsmax, double rsave[], const Integer &lrsave, Integer isave[], const Integer &lisave, const Integer &itask, const Integer &itrace, Integer &ind, Integer &ifail) |
}
|
The routine may be called by the names d03pff or nagf_pde_dim1_parab_convdiff.
3
Description
d03pff integrates the system of convection-diffusion equations in conservative form:
or the hyperbolic convection-only system:
for
, where the vector
is the set of solution values
The functions
,
,
and
depend on
,
and
; and
depends on
,
,
and
, where
is the spatial derivative of
. Note that
,
,
and
must not depend on any space derivatives; and none of the functions may depend on time derivatives. In terms of conservation laws,
,
and
are the convective flux, diffusion and source terms respectively.
The integration in time is from to , over the space interval , where and are the leftmost and rightmost points of a user-defined mesh . The initial values of the functions must be given at .
The PDEs are approximated by a system of ODEs in time for the values of
at mesh points using a spatial discretization method similar to the central-difference scheme used in
d03pcf/d03pca,
d03phf/d03pha and
d03ppf/d03ppa, but with the flux
replaced by a
numerical flux, which is a representation of the flux taking into account the direction of the flow of information at that point (i.e., the direction of the characteristics). Simple central differencing of the numerical flux then becomes a sophisticated upwind scheme in which the correct direction of upwinding is automatically achieved.
The numerical flux vector,
say, must be calculated by you in terms of the
left and
right values of the solution vector
(denoted by
and
respectively), at each mid-point of the mesh
, for
. The left and right values are calculated by
d03pff from two adjacent mesh points using a standard upwind technique combined with a Van Leer slope-limiter (see
LeVeque (1990)). The physically correct value for
is derived from the solution of the Riemann problem given by
where
, i.e.,
corresponds to
, with discontinuous initial values
for
and
for
, using an
approximate Riemann solver. This applies for either of the systems
(1) or
(2); the numerical flux is independent of the functions
,
,
and
. A description of several approximate Riemann solvers can be found in
LeVeque (1990) and
Berzins et al. (1989). Roe's scheme (see
Roe (1981)) is perhaps the easiest to understand and use, and a brief summary follows. Consider the system of PDEs
or equivalently
. Provided the system is linear in
, i.e., the Jacobian matrix
does not depend on
, the numerical flux
is given by
where
(
) is the flux
calculated at the left (right) value of
, denoted by
(
); the
are the eigenvalues of
; the
are the right eigenvectors of
; and the
are defined by
An example is given in
Section 10.
If the system is nonlinear, Roe's scheme requires that a linearized Jacobian is found (see
Roe (1981)).
The functions
,
,
and
(but
not
) must be specified in a
pdedef. The numerical flux
must be supplied in a separate
numflx. For problems in the form
(2), the actual argument
d03pfp may be used for
pdedef.
d03pfp is included in the NAG Library. This sets the matrix with entries
to the identity matrix, and the functions
,
and
to zero.
The boundary condition specification has sufficient flexibility to allow for different types of problems. For second-order problems, i.e.,
depending on
, a boundary condition is required for each PDE at both boundaries for the problem to be well-posed. If there are no second-order terms present, then the continuous PDE problem generally requires exactly one boundary condition for each PDE, that is
npde boundary conditions in total. However, in common with most discretization schemes for first-order problems, a
numerical boundary condition is required at the other boundary for each PDE. In order to be consistent with the characteristic directions of the PDE system, the numerical boundary conditions must be derived from the solution inside the domain in some manner (see below). You must supply both types of boundary conditions, i.e., a total of
npde conditions at each boundary point.
The position of each boundary condition should be chosen with care. In simple terms, if information is flowing into the domain then a physical boundary condition is required at that boundary, and a numerical boundary condition is required at the other boundary. In many cases the boundary conditions are simple, e.g., for the linear advection equation. In general you should calculate the characteristics of the PDE system and specify a physical boundary condition for each of the characteristic variables associated with incoming characteristics, and a numerical boundary condition for each outgoing characteristic.
A common way of providing numerical boundary conditions is to extrapolate the characteristic variables from the inside of the domain. Note that only linear extrapolation is allowed in this routine (for greater flexibility the routine
d03plf should be used). For problems in which the solution is known to be uniform (in space) towards a boundary during the period of integration then extrapolation is unnecessary; the numerical boundary condition can be supplied as the known solution at the boundary. Examples can be found in
Section 10.
The boundary conditions must be specified in
bndary in the form
at the left-hand boundary, and
at the right-hand boundary.
Note that spatial derivatives at the boundary are not passed explicitly to
bndary, but they can be calculated using values of
at and adjacent to the boundaries if required. However, it should be noted that instabilities may occur if such one-sided differencing opposes the characteristic direction at the boundary.
The problem is subject to the following restrictions:
-
(i), , and must not depend on any space derivatives;
-
(ii), , , and must not depend on any time derivatives;
-
(iii), so that integration is in the forward direction;
-
(iv)The evaluation of the terms , , and is done by calling the pdedef at a point approximately midway between each pair of mesh points in turn. Any discontinuities in these functions must therefore be at one or more of the mesh points ;
-
(v)At least one of the functions must be nonzero so that there is a time derivative present in the PDE problem.
In total there are ODEs in the time direction. This system is then integrated forwards in time using a BDF method.
For further details of the algorithm, see
Pennington and Berzins (1994) and the references therein.
4
References
Berzins M, Dew P M and Furzeland R M (1989) Developing software for time-dependent problems using the method of lines and differential-algebraic integrators Appl. Numer. Math. 5 375–397
Hirsch C (1990) Numerical Computation of Internal and External Flows, Volume 2: Computational Methods for Inviscid and Viscous Flows John Wiley
LeVeque R J (1990) Numerical Methods for Conservation Laws Birkhäuser Verlag
Pennington S V and Berzins M (1994) New NAG Library software for first-order partial differential equations ACM Trans. Math. Softw. 20 63–99
Roe P L (1981) Approximate Riemann solvers, parameter vectors, and difference schemes J. Comput. Phys. 43 357–372
5
Arguments
-
1:
– Integer
Input
-
On entry: the number of PDEs to be solved.
Constraint:
.
-
2:
– Real (Kind=nag_wp)
Input/Output
-
On entry: the initial value of the independent variable .
On exit: the value of
corresponding to the solution values in
u. Normally
.
Constraint:
.
-
3:
– Real (Kind=nag_wp)
Input
-
On entry: the final value of to which the integration is to be carried out.
-
4:
– Subroutine, supplied by the NAG Library or the user.
External Procedure
-
pdedef must evaluate the functions
,
,
and
which partially define the system of PDEs.
,
and
may depend on
,
and
;
may depend on
,
,
and
.
pdedef is called approximately midway between each pair of mesh points in turn by
d03pff. The actual argument
d03pfp may be used for
pdedef for problems in the form
(2). (
d03pfp is included in the NAG Library.)
The specification of
pdedef is:
Fortran Interface
Integer, Intent (In) |
:: |
npde |
Integer, Intent (Inout) |
:: |
ires |
Real (Kind=nag_wp), Intent (In) |
:: |
t, x, u(npde), ux(npde) |
Real (Kind=nag_wp), Intent (Out) |
:: |
p(npde,npde), c(npde), d(npde), s(npde) |
|
C Header Interface
void |
pdedef_ (const Integer *npde, const double *t, const double *x, const double u[], const double ux[], double p[], double c[], double d[], double s[], Integer *ires) |
|
C++ Header Interface
#include <nag.h> extern "C" {
void |
pdedef_ (const Integer &npde, const double &t, const double &x, const double u[], const double ux[], double p[], double c[], double d[], double s[], Integer &ires) |
}
|
-
1:
– Integer
Input
-
On entry: the number of PDEs in the system.
-
2:
– Real (Kind=nag_wp)
Input
-
On entry: the current value of the independent variable .
-
3:
– Real (Kind=nag_wp)
Input
-
On entry: the current value of the space variable .
-
4:
– Real (Kind=nag_wp) array
Input
-
On entry: contains the value of the component , for .
-
5:
– Real (Kind=nag_wp) array
Input
-
On entry: contains the value of the component , for .
-
6:
– Real (Kind=nag_wp) array
Output
-
On exit: must be set to the value of , for and .
-
7:
– Real (Kind=nag_wp) array
Output
-
On exit: must be set to the value of , for .
-
8:
– Real (Kind=nag_wp) array
Output
-
On exit: must be set to the value of , for .
-
9:
– Real (Kind=nag_wp) array
Output
-
On exit: must be set to the value of , for .
-
10:
– Integer
Input/Output
-
On entry: set to .
On exit: should usually remain unchanged. However, you may set
ires to force the integration routine to take certain actions as described below:
- Indicates to the integrator that control should be passed back immediately to the calling subroutine with the error indicator set to .
- Indicates to the integrator that the current time step should be abandoned and a smaller time step used instead. You may wish to set when a physically meaningless input or output value has been generated. If you consecutively set , d03pff returns to the calling subroutine with the error indicator set to .
pdedef must either be a module subprogram USEd by, or declared as EXTERNAL in, the (sub)program from which
d03pff is called. Arguments denoted as
Input must
not be changed by this procedure.
Note: pdedef should not return floating-point NaN (Not a Number) or infinity values, since these are not handled by
d03pff. If your code inadvertently
does return any NaNs or infinities,
d03pff is likely to produce unexpected results.
-
5:
– Subroutine, supplied by the user.
External Procedure
-
numflx must supply the numerical flux for each PDE given the
left and
right values of the solution vector
.
numflx is called approximately midway between each pair of mesh points in turn by
d03pff.
The specification of
numflx is:
Fortran Interface
Integer, Intent (In) |
:: |
npde |
Integer, Intent (Inout) |
:: |
ires |
Real (Kind=nag_wp), Intent (In) |
:: |
t, x, uleft(npde), uright(npde) |
Real (Kind=nag_wp), Intent (Out) |
:: |
flux(npde) |
|
C Header Interface
void |
numflx_ (const Integer *npde, const double *t, const double *x, const double uleft[], const double uright[], double flux[], Integer *ires) |
|
C++ Header Interface
#include <nag.h> extern "C" {
void |
numflx_ (const Integer &npde, const double &t, const double &x, const double uleft[], const double uright[], double flux[], Integer &ires) |
}
|
-
1:
– Integer
Input
-
On entry: the number of PDEs in the system.
-
2:
– Real (Kind=nag_wp)
Input
-
On entry: the current value of the independent variable .
-
3:
– Real (Kind=nag_wp)
Input
-
On entry: the current value of the space variable .
-
4:
– Real (Kind=nag_wp) array
Input
-
On entry: contains the left value of the component , for .
-
5:
– Real (Kind=nag_wp) array
Input
-
On entry: contains the right value of the component , for .
-
6:
– Real (Kind=nag_wp) array
Output
-
On exit: must be set to the numerical flux , for .
-
7:
– Integer
Input/Output
-
On entry: set to .
On exit: should usually remain unchanged. However, you may set
ires to force the integration routine to take certain actions as described below:
- Indicates to the integrator that control should be passed back immediately to the calling subroutine with the error indicator set to .
- Indicates to the integrator that the current time step should be abandoned and a smaller time step used instead. You may wish to set when a physically meaningless input or output value has been generated. If you consecutively set , d03pff returns to the calling subroutine with the error indicator set to .
numflx must either be a module subprogram USEd by, or declared as EXTERNAL in, the (sub)program from which
d03pff is called. Arguments denoted as
Input must
not be changed by this procedure.
Note: numflx should not return floating-point NaN (Not a Number) or infinity values, since these are not handled by
d03pff. If your code inadvertently
does return any NaNs or infinities,
d03pff is likely to produce unexpected results.
-
6:
– Subroutine, supplied by the user.
External Procedure
-
bndary must evaluate the functions
and
which describe the physical and numerical boundary conditions, as given by
(6) and
(7).
The specification of
bndary is:
Fortran Interface
Integer, Intent (In) |
:: |
npde, npts, ibnd |
Integer, Intent (Inout) |
:: |
ires |
Real (Kind=nag_wp), Intent (In) |
:: |
t, x(npts), u(npde,3) |
Real (Kind=nag_wp), Intent (Out) |
:: |
g(npde) |
|
C Header Interface
void |
bndary_ (const Integer *npde, const Integer *npts, const double *t, const double x[], const double u[], const Integer *ibnd, double g[], Integer *ires) |
|
C++ Header Interface
#include <nag.h> extern "C" {
void |
bndary_ (const Integer &npde, const Integer &npts, const double &t, const double x[], const double u[], const Integer &ibnd, double g[], Integer &ires) |
}
|
-
1:
– Integer
Input
-
On entry: the number of PDEs in the system.
-
2:
– Integer
Input
-
On entry: the number of mesh points in the interval .
-
3:
– Real (Kind=nag_wp)
Input
-
On entry: the current value of the independent variable .
-
4:
– Real (Kind=nag_wp) array
Input
-
On entry: the mesh points in the spatial direction. corresponds to the left-hand boundary, , and corresponds to the right-hand boundary, .
-
5:
– Real (Kind=nag_wp) array
Input
-
On entry: contains the value of solution components in the boundary region.
If , contains the value of the component at , for and .
If , contains the value of the component at , for and .
-
6:
– Integer
Input
-
On entry: specifies which boundary conditions are to be evaluated.
- bndary must evaluate the left-hand boundary condition at .
- bndary must evaluate the right-hand boundary condition at .
-
7:
– Real (Kind=nag_wp) array
Output
-
On exit:
must contain the
th component of either
or
in
(6) and
(7), depending on the value of
ibnd, for
.
-
8:
– Integer
Input/Output
-
On entry: set to .
On exit: should usually remain unchanged. However, you may set
ires to force the integration routine to take certain actions as described below:
- Indicates to the integrator that control should be passed back immediately to the calling subroutine with the error indicator set to .
- Indicates to the integrator that the current time step should be abandoned and a smaller time step used instead. You may wish to set when a physically meaningless input or output value has been generated. If you consecutively set , d03pff returns to the calling subroutine with the error indicator set to .
bndary must either be a module subprogram USEd by, or declared as EXTERNAL in, the (sub)program from which
d03pff is called. Arguments denoted as
Input must
not be changed by this procedure.
Note: bndary should not return floating-point NaN (Not a Number) or infinity values, since these are not handled by
d03pff. If your code inadvertently
does return any NaNs or infinities,
d03pff is likely to produce unexpected results.
-
7:
– Real (Kind=nag_wp) array
Input/Output
-
On entry: must contain the initial value of at and , for and .
On exit: will contain the computed solution at and , for and .
-
8:
– Integer
Input
-
On entry: the number of mesh points in the interval .
Constraint:
.
-
9:
– Real (Kind=nag_wp) array
Input
-
On entry: the mesh points in the space direction. must specify the left-hand boundary, , and must specify the right-hand boundary, .
Constraint:
.
-
10:
– Real (Kind=nag_wp) array
Input
-
On entry: the components of
acc contain the relative and absolute error tolerances used in the local error test in the time integration.
If
is the estimated error for
at the
th mesh point, the error test is
Constraint:
and (but not both zero).
-
11:
– Real (Kind=nag_wp)
Input
-
On entry: the maximum absolute step size to be allowed in the time integration. If then no maximum is imposed.
Constraint:
.
-
12:
– Real (Kind=nag_wp) array
Communication Array
-
If
,
rsave need not be set on entry.
If
,
rsave must be unchanged from the previous call to the routine because it contains required information about the iteration.
-
13:
– Integer
Input
-
On entry: the dimension of the array
rsave as declared in the (sub)program from which
d03pff is called.
Constraint:
.
-
14:
– Integer array
Communication Array
-
If
,
isave need not be set on entry.
If
,
isave must be unchanged from the previous call to the routine because it contains required information about the iteration. In particular:
- Contains the number of steps taken in time.
- Contains the number of residual evaluations of the resulting ODE system used. One such evaluation involves computing the PDE functions at all the mesh points, as well as one evaluation of the functions in the boundary conditions.
- Contains the number of Jacobian evaluations performed by the time integrator.
- Contains the order of the last backward differentiation formula method used.
- Contains the number of Newton iterations performed by the time integrator. Each iteration involves an ODE residual evaluation followed by a back-substitution using the decomposition of the Jacobian matrix.
-
15:
– Integer
Input
-
On entry: the dimension of the array
isave as declared in the (sub)program from which
d03pff is called.
Constraint:
.
-
16:
– Integer
Input
-
On entry: the task to be performed by the ODE integrator.
- Normal computation of output values at (by overshooting and interpolating).
- Take one step in the time direction and return.
- Stop at first internal integration point at or beyond .
Constraint:
, or .
-
17:
– Integer
Input
-
On entry: the level of trace information required from
d03pff and the underlying ODE solver.
itrace may take the value
,
,
,
or
.
- No output is generated.
- Only warning messages from the PDE solver are printed on the current error message unit (see x04aaf).
- Output from the underlying ODE solver is printed on the current advisory message unit (see x04abf). This output contains details of Jacobian entries, the nonlinear iteration and the time integration during the computation of the ODE system.
If , is assumed and similarly if , is assumed.
The advisory messages are given in greater detail as
itrace increases. You are advised to set
, unless you are experienced with
Sub-chapter D02MN.
-
18:
– Integer
Input/Output
-
On entry: indicates whether this is a continuation call or a new integration.
- Starts or restarts the integration in time.
- Continues the integration after an earlier exit from the routine. In this case, only the
arguments tout and ifail
should be reset between calls to d03pff.
Constraint:
or .
On exit: .
-
19:
– Integer
Input/Output
-
On entry:
ifail must be set to
,
. If you are unfamiliar with this argument you should refer to
Section 4 in the Introduction to the NAG Library FL Interface for details.
For environments where it might be inappropriate to halt program execution when an error is detected, the value
is recommended. If the output of error messages is undesirable, then the value
is recommended. Otherwise, if you are not familiar with this argument, the recommended value is
.
When the value is used it is essential to test the value of ifail on exit.
On exit:
unless the routine detects an error or a warning has been flagged (see
Section 6).
6
Error Indicators and Warnings
If on entry
or
, explanatory error messages are output on the current error message unit (as defined by
x04aaf).
Errors or warnings detected by the routine:
-
On entry, .
Constraint: .
On entry, and are both zero.
On entry, .
Constraint: .
On entry, , , and .
Constraint: .
On entry, .
Constraint: or .
On entry, .
Constraint: , or .
On entry, .
Constraint: .
On entry, .
Constraint: .
On entry, .
Constraint: .
On entry, .
Constraint: .
On entry, on initial entry .
Constraint: on initial entry .
On entry, and .
Constraint: .
On entry, is too small:
and .
On entry, .
Constraint: .
-
Underlying ODE solver cannot make further progress from the point
ts with the supplied values of
acc.
,
,
.
-
Repeated errors in an attempted step of underlying ODE solver. Integration was successful as far as
ts:
.
In the underlying ODE solver, there were repeated error test failures on an attempted step, before completing the requested task, but the integration was successful as far as . The problem may have a singularity, or the error requirement may be inappropriate. Incorrect specification of boundary conditions may also result in this error.
-
In setting up the ODE system an internal auxiliary was unable to initialize the derivative. This could be due to your setting
in
pdedef,
numflx, or
bndary.
-
Singular Jacobian of ODE system. Check problem formulation.
-
In evaluating residual of ODE system,
has been set in
pdedef,
numflx, or
bndary. Integration is successful as far as
ts:
.
-
Values in
acc are too small to start integration:
,
.
-
ires set to an invalid value in call to
pdedef,
numflx, or
bndary.
-
Serious error in internal call to an auxiliary. Increase
itrace for further details.
-
Integration completed, but small changes in
acc are unlikely to result in a changed solution.
,
.
The required task has been completed, but it is estimated that a small change in the values of
acc is unlikely to produce any change in the computed solution. (Only applies when you are not operating in one step mode, that is when
.)
-
Error during Jacobian formulation for ODE system. Increase
itrace for further details.
-
The functions , , or appear to depend on time derivatives.
An unexpected error has been triggered by this routine. Please
contact
NAG.
See
Section 7 in the Introduction to the NAG Library FL Interface for further information.
Your licence key may have expired or may not have been installed correctly.
See
Section 8 in the Introduction to the NAG Library FL Interface for further information.
Dynamic memory allocation failed.
See
Section 9 in the Introduction to the NAG Library FL Interface for further information.
7
Accuracy
d03pff controls the accuracy of the integration in the time direction but not the accuracy of the approximation in space. The spatial accuracy depends on both the number of mesh points and on their distribution in space. In the time integration only the local error over a single step is controlled and so the accuracy over a number of steps cannot be guaranteed. You should therefore test the effect of varying the components of the accuracy argument,
acc.
8
Parallelism and Performance
d03pff is not thread safe and should not be called from a multithreaded user program. Please see
Section 1 in FL Interface Multithreading for more information on thread safety.
d03pff is threaded by NAG for parallel execution in multithreaded implementations of the NAG Library.
d03pff makes calls to BLAS and/or LAPACK routines, which may be threaded within the vendor library used by this implementation. Consult the documentation for the vendor library for further information.
Please consult the
X06 Chapter Introduction for information on how to control and interrogate the OpenMP environment used within this routine. Please also consult the
Users' Note for your implementation for any additional implementation-specific information.
d03pff is designed to solve systems of PDEs in conservative form, with optional source terms which are independent of space derivatives, and optional second-order diffusion terms. The use of the routine to solve systems which are not naturally in this form is discouraged, and you are advised to use one of the central-difference schemes for such problems.
You should be aware of the stability limitations for hyperbolic PDEs. For most problems with small error tolerances the ODE integrator does not attempt unstable time steps, but in some cases a maximum time step should be imposed using
tsmax. It is worth experimenting with this argument, particularly if the integration appears to progress unrealistically fast (with large time steps). Setting the maximum time step to the minimum mesh size is a safe measure, although in some cases this may be too restrictive.
Problems with source terms should be treated with caution, as it is known that for large source terms stable and reasonable looking solutions can be obtained which are in fact incorrect, exhibiting non-physical speeds of propagation of discontinuities (typically one spatial mesh point per time step). It is essential to employ a very fine mesh for problems with source terms and discontinuities, and to check for non-physical propagation speeds by comparing results for different mesh sizes. Further details and an example can be found in
Pennington and Berzins (1994).
The time taken depends on the complexity of the system and on the accuracy requested.
10
Example
For this routine two examples are presented. There is a single example program for d03pff, with a main program and the code to solve the two example problems given in Example 1 (EX1) and Example 2 (EX2).
Example 1 (EX1)
This example is a simple first-order system which illustrates the calculation of the numerical flux using Roe's approximate Riemann solver, and the specification of numerical boundary conditions using extrapolated characteristic variables. The PDEs are
for
and
. The PDEs have an exact solution given by
The initial conditions are given by the exact solution. The characteristic variables are
and
corresponding to the characteristics given by
and
respectively. Hence a physical boundary condition is required for
at the left-hand boundary, and for
at the right-hand boundary (corresponding to the incoming characteristics); and a numerical boundary condition is required for
at the left-hand boundary, and for
at the right-hand boundary (outgoing characteristics). The physical boundary conditions are obtained from the exact solution, and the numerical boundary conditions are calculated by linear extrapolation of the appropriate characteristic variable. The numerical flux is calculated using Roe's approximate Riemann solver: Using the notation in
Section 3, the flux vector
and the Jacobian matrix
are
and the eigenvalues of
are
and
with right eigenvectors
and
respectively. Using equation
(4) the
are given by
that is
is given by
and similarly for
. From equation
(4), the numerical flux vector is
that is
Example 2 (EX2)
This example is an advection-diffusion equation in which the flux term depends explicitly on
:
for
and
. The argument
is taken to be
. The two physical boundary conditions are
and
and the initial condition is
. The integration is run to steady state at which the solution is known to be
across the domain with a narrow boundary layer at both boundaries. In order to write the PDE in conservative form, a source term must be introduced, i.e.,
As in Example 1, the numerical flux is calculated using the Roe approximate Riemann solver. The Riemann problem to solve locally is
The
in the flux term is assumed to be constant at a local level, and so using the notation in
Section 3,
and
. The eigenvalue is
and the eigenvector (a scalar in this case) is
. The numerical flux is therefore
10.1
Program Text
10.2
Program Data
10.3
Program Results