NAG Library Routine Document

d03pff (dim1_parab_convdiff)

 Contents

    1  Purpose
    7  Accuracy

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 nagmk26.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)

3
Description

d03pff integrates the system of convection-diffusion equations in conservative form:
j=1npdePi,j Uj t + Fi x =Ci Di x +Si, (1)
or the hyperbolic convection-only system:
Ui t + Fi x =0, (2)
for i=1,2,,npde,  axb,  tt0, where the vector U is the set of solution values
U x,t = U 1 x,t ,, U npde x,t T .  
The functions Pi,j, Fi, Ci and Si depend on x, t and U; and Di depends on x, t, U and Ux, where Ux is the spatial derivative of U. Note that Pi,j, Fi, Ci and Si must not depend on any space derivatives; and none of the functions may depend on time derivatives. In terms of conservation laws, Fi, CiDi x  and Si are the convective flux, diffusion and source terms respectively.
The integration in time is from t0 to tout, over the space interval axb, where a=x1 and b=xnpts are the leftmost and rightmost points of a user-defined mesh x1,x2,,xnpts. The initial values of the functions Ux,t must be given at t=t0.
The PDEs are approximated by a system of ODEs in time for the values of Ui 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 Fi 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, F^i say, must be calculated by you in terms of the left and right values of the solution vector U (denoted by UL and UR respectively), at each mid-point of the mesh xj-1/2=xj-1+xj/2, for j=2,3,,npts. 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 F^i is derived from the solution of the Riemann problem given by
Ui t + Fi y =0, (3)
where y=x-xj-1/2, i.e., y=0 corresponds to x=xj-1/2, with discontinuous initial values U=UL for y<0 and U=UR for y>0, using an approximate Riemann solver. This applies for either of the systems (1) or (2); the numerical flux is independent of the functions Pi,j, Ci, Di and Si. 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 Ut+Fx=0 or equivalently Ut+AUx=0. Provided the system is linear in U, i.e., the Jacobian matrix A does not depend on U, the numerical flux F^ is given by
F^=12 FL+FR-12k=1npdeαkλkek, (4)
where FL (FR) is the flux F calculated at the left (right) value of U, denoted by UL (UR); the λk are the eigenvalues of A; the ek are the right eigenvectors of A; and the αk are defined by
UR-UL=k=1npdeαkek. (5)
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 Pi,j, Ci, Di and Si (but not Fi) must be specified in a pdedef. The numerical flux F^i 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 and sets the matrix with entries Pi,j to the identity matrix, and the functions Ci, Di and Si to zero.
The boundary condition specification has sufficient flexibility to allow for different types of problems. For second-order problems, i.e., Di depending on Ux, 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
GiLx,t,U=0  at ​x=a,  i=1,2,,npde, (6)
at the left-hand boundary, and
GiRx,t,U=0  at ​x=b,  i=1,2,,npde, (7)
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 U 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) Pi,j, Fi, Ci and Si must not depend on any space derivatives;
(ii) Pi,j, Fi, Ci, Di and Si must not depend on any time derivatives;
(iii) t0<tout, so that integration is in the forward direction;
(iv) The evaluation of the terms Pi,j, Ci, Di and Si 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 x1,x2,,xnpts;
(v) At least one of the functions Pi,j must be nonzero so that there is a time derivative present in the PDE problem.
In total there are npde×npts 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:     npde – IntegerInput
On entry: the number of PDEs to be solved.
Constraint: npde1.
2:     ts – Real (Kind=nag_wp)Input/Output
On entry: the initial value of the independent variable t.
On exit: the value of t corresponding to the solution values in u. Normally ts=tout.
Constraint: ts<tout.
3:     tout – Real (Kind=nag_wp)Input
On entry: the final value of t to which the integration is to be carried out.
4:     pdedef – Subroutine, supplied by the NAG Library or the user.External Procedure
pdedef must evaluate the functions Pi,j, Ci, Di and Si which partially define the system of PDEs. Pi,j, Ci and Si may depend on x, t and U; Di may depend on x, t, U and Ux. 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
Subroutine pdedef ( npde, t, x, u, ux, p, c, d, s, ires)
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
#include nagmk26.h
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:     npde – IntegerInput
On entry: the number of PDEs in the system.
2:     t – Real (Kind=nag_wp)Input
On entry: the current value of the independent variable t.
3:     x – Real (Kind=nag_wp)Input
On entry: the current value of the space variable x.
4:     unpde – Real (Kind=nag_wp) arrayInput
On entry: ui contains the value of the component Uix,t, for i=1,2,,npde.
5:     uxnpde – Real (Kind=nag_wp) arrayInput
On entry: uxi contains the value of the component Uix,t x , for i=1,2,,npde.
6:     pnpdenpde – Real (Kind=nag_wp) arrayOutput
On exit: pij must be set to the value of Pi,jx,t,U, for i=1,2,,npde and j=1,2,,npde.
7:     cnpde – Real (Kind=nag_wp) arrayOutput
On exit: ci must be set to the value of Cix,t,U, for i=1,2,,npde.
8:     dnpde – Real (Kind=nag_wp) arrayOutput
On exit: di must be set to the value of Dix,t,U,Ux, for i=1,2,,npde.
9:     snpde – Real (Kind=nag_wp) arrayOutput
On exit: si must be set to the value of Six,t,U, for i=1,2,,npde.
10:   ires – IntegerInput/Output
On entry: set to -1​ or ​1.
On exit: should usually remain unchanged. However, you may set ires to force the integration routine to take certain actions as described below:
ires=2
Indicates to the integrator that control should be passed back immediately to the calling subroutine with the error indicator set to ifail=6.
ires=3
Indicates to the integrator that the current time step should be abandoned and a smaller time step used instead. You may wish to set ires=3 when a physically meaningless input or output value has been generated. If you consecutively set ires=3, d03pff returns to the calling subroutine with the error indicator set to ifail=4.
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:     numflx – 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 u. numflx is called approximately midway between each pair of mesh points in turn by d03pff.
The specification of numflx is:
Fortran Interface
Subroutine numflx ( npde, t, x, uleft, uright, flux, ires)
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
#include nagmk26.h
void  numflx (const Integer *npde, const double *t, const double *x, const double uleft[], const double uright[], double flux[], Integer *ires)
1:     npde – IntegerInput
On entry: the number of PDEs in the system.
2:     t – Real (Kind=nag_wp)Input
On entry: the current value of the independent variable t.
3:     x – Real (Kind=nag_wp)Input
On entry: the current value of the space variable x.
4:     uleftnpde – Real (Kind=nag_wp) arrayInput
On entry: ulefti contains the left value of the component Uix, for i=1,2,,npde.
5:     urightnpde – Real (Kind=nag_wp) arrayInput
On entry: urighti contains the right value of the component Uix, for i=1,2,,npde.
6:     fluxnpde – Real (Kind=nag_wp) arrayOutput
On exit: fluxi must be set to the numerical flux F^i, for i=1,2,,npde.
7:     ires – IntegerInput/Output
On entry: set to -1​ or ​1.
On exit: should usually remain unchanged. However, you may set ires to force the integration routine to take certain actions as described below:
ires=2
Indicates to the integrator that control should be passed back immediately to the calling subroutine with the error indicator set to ifail=6.
ires=3
Indicates to the integrator that the current time step should be abandoned and a smaller time step used instead. You may wish to set ires=3 when a physically meaningless input or output value has been generated. If you consecutively set ires=3, d03pff returns to the calling subroutine with the error indicator set to ifail=4.
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:     bndary – Subroutine, supplied by the user.External Procedure
bndary must evaluate the functions GiL and GiR which describe the physical and numerical boundary conditions, as given by (6) and (7).
The specification of bndary is:
Fortran Interface
Subroutine bndary ( npde, npts, t, x, u, ibnd, g, ires)
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
#include nagmk26.h
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:     npde – IntegerInput
On entry: the number of PDEs in the system.
2:     npts – IntegerInput
On entry: the number of mesh points in the interval a,b.
3:     t – Real (Kind=nag_wp)Input
On entry: the current value of the independent variable t.
4:     xnpts – Real (Kind=nag_wp) arrayInput
On entry: the mesh points in the spatial direction. x1 corresponds to the left-hand boundary, a, and xnpts corresponds to the right-hand boundary, b.
5:     unpde3 – Real (Kind=nag_wp) arrayInput
On entry: contains the value of solution components in the boundary region.
If ibnd=0, uij contains the value of the component Uix ,t at x=xj, for i=1,2,,npde and j=1,2,3.
If ibnd0, uij contains the value of the component Uix,t at x=xnpts-j+1, for i=1,2,,npde and j=1,2,3.
6:     ibnd – IntegerInput
On entry: specifies which boundary conditions are to be evaluated.
ibnd=0
bndary must evaluate the left-hand boundary condition at x=a.
ibnd0
bndary must evaluate the right-hand boundary condition at x=b.
7:     gnpde – Real (Kind=nag_wp) arrayOutput
On exit: gi must contain the ith component of either gL or gR in (6) and (7), depending on the value of ibnd, for i=1,2,,npde.
8:     ires – IntegerInput/Output
On entry: set to -1​ or ​1.
On exit: should usually remain unchanged. However, you may set ires to force the integration routine to take certain actions as described below:
ires=2
Indicates to the integrator that control should be passed back immediately to the calling subroutine with the error indicator set to ifail=6.
ires=3
Indicates to the integrator that the current time step should be abandoned and a smaller time step used instead. You may wish to set ires=3 when a physically meaningless input or output value has been generated. If you consecutively set ires=3, d03pff returns to the calling subroutine with the error indicator set to ifail=4.
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:     unpdenpts – Real (Kind=nag_wp) arrayInput/Output
On entry: uij must contain the initial value of Uix,t at x=xj and t=ts, for i=1,2,,npde and j=1,2,,npts.
On exit: uij will contain the computed solution Uix,t at x=xj and t=ts, for i=1,2,,npde and j=1,2,,npts.
8:     npts – IntegerInput
On entry: the number of mesh points in the interval a,b.
Constraint: npts3.
9:     xnpts – Real (Kind=nag_wp) arrayInput
On entry: the mesh points in the space direction. x1 must specify the left-hand boundary, a, and xnpts must specify the right-hand boundary, b.
Constraint: x1<x2<<xnpts.
10:   acc2 – Real (Kind=nag_wp) arrayInput
On entry: the components of acc contain the relative and absolute error tolerances used in the local error test in the time integration.
If Ei,j is the estimated error for Ui at the jth mesh point, the error test is
Ei,j=acc1×uij+acc2.  
Constraint: acc1 and acc20.0 (but not both zero).
11:   tsmax – Real (Kind=nag_wp)Input
On entry: the maximum absolute step size to be allowed in the time integration. If tsmax=0.0 then no maximum is imposed.
Constraint: tsmax0.0.
12:   rsavelrsave – Real (Kind=nag_wp) arrayCommunication Array
If ind=0, rsave need not be set on entry.
If ind=1, rsave must be unchanged from the previous call to the routine because it contains required information about the iteration.
13:   lrsave – IntegerInput
On entry: the dimension of the array rsave as declared in the (sub)program from which d03pff is called.
Constraint: lrsave11+9×npde×npde×npts+32+3×npde×npde+7×npts+54.
14:   isavelisave – Integer arrayCommunication Array
If ind=0, isave need not be set on entry.
If ind=1, isave must be unchanged from the previous call to the routine because it contains required information about the iteration. In particular:
isave1
Contains the number of steps taken in time.
isave2
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.
isave3
Contains the number of Jacobian evaluations performed by the time integrator.
isave4
Contains the order of the last backward differentiation formula method used.
isave5
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 LU decomposition of the Jacobian matrix.
15:   lisave – IntegerInput
On entry: the dimension of the array isave as declared in the (sub)program from which d03pff is called.
Constraint: lisavenpde×npts+24.
16:   itask – IntegerInput
On entry: the task to be performed by the ODE integrator.
itask=1
Normal computation of output values u at t=tout (by overshooting and interpolating).
itask=2
Take one step in the time direction and return.
itask=3
Stop at first internal integration point at or beyond t=tout.
Constraint: itask=1, 2 or 3.
17:   itrace – IntegerInput
On entry: the level of trace information required from d03pff and the underlying ODE solver. itrace may take the value -1, 0, 1, 2 or 3.
itrace=-1
No output is generated.
itrace=0
Only warning messages from the PDE solver are printed on the current error message unit (see x04aaf).
itrace>0
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 itrace<-1, -1 is assumed and similarly if itrace>3, 3 is assumed.
The advisory messages are given in greater detail as itrace increases. You are advised to set itrace=0, unless you are experienced with Sub-chapter D02M–N.
18:   ind – IntegerInput/Output
On entry: indicates whether this is a continuation call or a new integration.
ind=0
Starts or restarts the integration in time.
ind=1
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: ind=0 or 1.
On exit: ind=1.
19:   ifail – IntegerInput/Output
On entry: ifail must be set to 0, -1​ or ​1. If you are unfamiliar with this argument you should refer to Section 3.4 in How to Use the NAG Library and its Documentation for details.
For environments where it might be inappropriate to halt program execution when an error is detected, the value -1​ or ​1 is recommended. If the output of error messages is undesirable, then the value 1 is recommended. Otherwise, if you are not familiar with this argument, the recommended value is 0. When the value -1​ or ​1 is used it is essential to test the value of ifail on exit.
On exit: ifail=0 unless the routine detects an error or a warning has been flagged (see Section 6).

6
Error Indicators and Warnings

If on entry ifail=0 or -1, explanatory error messages are output on the current error message unit (as defined by x04aaf).
Errors or warnings detected by the routine:
ifail=1
On entry,tstout,
ortout-ts is too small,
oritask=1, 2 or 3,
ornpts<3,
ornpde<1,
orind0 or 1,
orincorrect user-defined mesh, i.e., xixi+1 for some i=1,2,,npts-1,
orlrsave or lisave are too small,
orind=1 on initial entry to d03pff,
oracc1 or acc2<0.0,
oracc1 or acc2 are both zero,
ortsmax<0.0.
ifail=2
The underlying ODE solver cannot make any further progress, with the values of acc, across the integration range from the current point t=ts. The components of u contain the computed values at the current point t=ts.
ifail=3
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 t=ts. The problem may have a singularity, or the error requirement may be inappropriate. Incorrect specification of boundary conditions may also result in this error.
ifail=4
In setting up the ODE system, the internal initialization routine was unable to initialize the derivative of the ODE system. This could be due to the fact that ires was repeatedly set to 3 in one of pdedef, numflx or bndary when the residual in the underlying ODE solver was being evaluated. Incorrect specification of boundary conditions may also result in this error.
ifail=5
In solving the ODE system, a singular Jacobian has been encountered. Check the problem formulation.
ifail=6
When evaluating the residual in solving the ODE system, ires was set to 2 in at least one of pdedef, numflx or bndary. Integration was successful as far as t=ts.
ifail=7
The values of acc1 and acc2 are so small that the routine is unable to start the integration in time.
ifail=8
In either, pdedef, numflx or bndary, ires was set to an invalid value.
ifail=9 (d02nnf)
A serious error has occurred in an internal call to the specified routine. Check the problem specification and all arguments and array dimensions. Setting itrace=1 may provide more information. If the problem persists, contact NAG.
ifail=10
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 itask2.)
ifail=11
An error occurred during Jacobian formulation of the ODE system (a more detailed error description may be directed to the current advisory message unit when itrace1).
ifail=12
Not applicable.
ifail=13
Not applicable.
ifail=14
One or more of the functions Pi,j, Di or Ci was detected as depending on time derivatives, which is not permissible.
ifail=-99
An unexpected error has been triggered by this routine. Please contact NAG.
See Section 3.9 in How to Use the NAG Library and its Documentation for further information.
ifail=-399
Your licence key may have expired or may not have been installed correctly.
See Section 3.8 in How to Use the NAG Library and its Documentation for further information.
ifail=-999
Dynamic memory allocation failed.
See Section 3.7 in How to Use the NAG Library and its Documentation 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 3.12.1 in How to Use the NAG Library and its Documentation 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.

9
Further Comments

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
U1 t + U1 x + U2 x = 0, U2 t +4 U1 x + U2 x = 0,  
for x 0,1  and t 0 . The PDEs have an exact solution given by
U1 x,t = 12 exp x + t + expx-3t + 14 sin 2 π x-3t 2 - sin 2 π x+t 2 + 2 t2 - 2 x t , U2 x,t = expx-3t - expx+t + 12 sin 2 π x-3t 2 + sin 2 π x - 3 t 2 + x2 + 5 t2 - 2 x t .  
The initial conditions are given by the exact solution. The characteristic variables are 2U1+U2 and 2U1-U2 corresponding to the characteristics given by dx/dt=3 and dx/dt=-1 respectively. Hence a physical boundary condition is required for 2U1+U2 at the left-hand boundary, and for 2U1-U2 at the right-hand boundary (corresponding to the incoming characteristics); and a numerical boundary condition is required for 2U1-U2 at the left-hand boundary, and for 2U1+U2 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 F and the Jacobian matrix A are
F= U1+U2 4U1+U2   and   A= 1 1 4 1 ,  
and the eigenvalues of A are 3 and -1 with right eigenvectors 12T and -12T respectively. Using equation (4) the αk are given by
U1R-U1L U2R-U2L =α1 1 2 +α2 -1 2 ,  
that is
α1 = 14 2 U1R - 2 U1L + U2R - U2L   and   α2 = 14 -2 U1R + 2 U1L + U2R - U2L .  
FL is given by
FL = U1L+U2L 4U1L+U2L ,  
and similarly for FR. From equation (4), the numerical flux vector is
F^ = 12 U1L+U2L+0U1R+U2R 4U1L+U2L+4U1R+U2R - 12 α1 3 1 2 - 12 α2 -1 -1 2 ,  
that is
F^ = 12 3U1L-0U1R+32U2L+12 U2R 6U1L+ 2U1R+ 3U2L-0U2R .  
Example 2 (EX2)
This example is an advection-diffusion equation in which the flux term depends explicitly on x:
U t +x U x =ε 2U x2 ,  
for x-1,1 and 0t10. The argument ε is taken to be 0.01. The two physical boundary conditions are U-1,t=3.0 and U1,t=5.0 and the initial condition is Ux,0=x+4. The integration is run to steady state at which the solution is known to be U=4 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.,
U t + xU x =ε 2U x2 +U.  
As in Example 1, the numerical flux is calculated using the Roe approximate Riemann solver. The Riemann problem to solve locally is
U t + x U x =0.  
The x in the flux term is assumed to be constant at a local level, and so using the notation in Section 3, F=xU and A=x. The eigenvalue is x and the eigenvector (a scalar in this case) is 1. The numerical flux is therefore
F^= xUL if ​x0, xUR if ​x<0.  

10.1
Program Text

Program Text (d03pffe.f90)

10.2
Program Data

Program Data (d03pffe.d)

10.3
Program Results

Program Results (d03pffe.r)

© The Numerical Algorithms Group Ltd, Oxford, UK. 2017