D03PEF (PDF version)
D03 Chapter Contents
D03 Chapter Introduction
NAG Library Manual

NAG Library Routine Document

D03PEF

Note:  before using this routine, please read the Users' Note for your implementation to check the interpretation of bold italicised terms and other implementation-dependent details.

+ Contents

    1  Purpose
    7  Accuracy

1  Purpose

D03PEF integrates a system of linear or nonlinear, first-order, time-dependent partial differential equations (PDEs) in one space variable. The spatial discretization is performed using the Keller box scheme and the method of lines is employed to reduce the PDEs to a system of ordinary differential equations (ODEs). The resulting system is solved using a Backward Differentiation Formula (BDF) method.

2  Specification

SUBROUTINE D03PEF ( NPDE, TS, TOUT, PDEDEF, BNDARY, U, NPTS, X, NLEFT, ACC, RSAVE, LRSAVE, ISAVE, LISAVE, ITASK, ITRACE, IND, IFAIL)
INTEGER  NPDE, NPTS, NLEFT, LRSAVE, ISAVE(LISAVE), LISAVE, ITASK, ITRACE, IND, IFAIL
REAL (KIND=nag_wp)  TS, TOUT, U(NPDE,NPTS), X(NPTS), ACC, RSAVE(LRSAVE)
EXTERNAL  PDEDEF, BNDARY

3  Description

D03PEF integrates the system of first-order PDEs
Gix,t,U,Ux,Ut=0,  i=1,2,,NPDE. (1)
In particular the functions Gi must have the general form
Gi=j=1NPDEPi,j Uj t +Qi,  i=1,2,,NPDE,  axb,tt0, (2)
where Pi,j and Qi depend on x, t, U, Ux and the vector U is the set of solution values
U x,t = U 1 x,t ,, U NPDE x,t T , (3)
and the vector Ux is its partial derivative with respect to x. Note that Pi,j and Qi must not depend on U t .
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 mesh should be chosen in accordance with the expected behaviour of the solution.
The PDE system which is defined by the functions Gi must be specified in PDEDEF.
The initial values of the functions Ux,t must be given at t=t0. For a first-order system of PDEs, only one boundary condition is required for each PDE component Ui. The NPDE boundary conditions are separated into na at the left-hand boundary x=a, and nb at the right-hand boundary x=b, such that na+nb=NPDE. The position of the boundary condition for each component should be chosen with care; the general rule is that if the characteristic direction of Ui at the left-hand boundary (say) points into the interior of the solution domain, then the boundary condition for Ui should be specified at the left-hand boundary. Incorrect positioning of boundary conditions generally results in initialization or integration difficulties in the underlying time integration routines.
The boundary conditions have the form:
GiLx,t,U,Ut=0  at ​x=a,  i=1,2,,na (4)
at the left-hand boundary, and
GiRx,t,U,Ut=0  at ​x=b,  i=1,2,,nb (5)
at the right-hand boundary.
Note that the functions GiL and GiR must not depend on Ux, since spatial derivatives are not determined explicitly in the Keller box scheme (see Keller (1970)). If the problem involves derivative (Neumann) boundary conditions then it is generally possible to restate such boundary conditions in terms of permissible variables. Also note that GiL and GiR must be linear with respect to time derivatives, so that the boundary conditions have the general form
j=1NPDEEi,jL Uj t +SiL=0,   i=1,2,,na (6)
at the left-hand boundary, and
j=1NPDEEi,jR Uj t +SiR=0,   i=1,2,,nb (7)
at the right-hand boundary, where Ei,jL, Ei,jR, SiL, and SiR depend on x, t and U only.
The boundary conditions must be specified in BNDARY.
The problem is subject to the following restrictions:
(i) t0<tout, so that integration is in the forward direction;
(ii) Pi,j and Qi must not depend on any time derivatives;
(iii) The evaluation of the function Gi is done at the mid-points of the mesh intervals by calling the PDEDEF for each mid-point in turn. Any discontinuities in the function must therefore be at one or more of the mesh points x1,x2,,xNPTS;
(iv) At least one of the functions Pi,j must be nonzero so that there is a time derivative present in the problem.
In this method of lines approach the Keller box scheme (see Keller (1970)) is applied to each PDE in the space variable only, resulting in a system of ODEs in time for the values of Ui at each mesh point. In total there are NPDE×NPTS ODEs in the time direction. This system is then integrated forwards in time using a BDF method.

4  References

Berzins M (1990) Developments in the NAG Library software for parabolic equations Scientific Software Systems (eds J C Mason and M G Cox) 59–72 Chapman and Hall
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
Keller H B (1970) A new difference scheme for parabolic problems Numerical Solutions of Partial Differential Equations (ed J Bramble) 2 327–350 Academic Press
Pennington S V and Berzins M (1994) New NAG Library software for first-order partial differential equations ACM Trans. Math. Softw. 20 63–99

5  Parameters

1:     NPDE – INTEGERInput
On entry: the number of PDEs in the system to be solved.
Constraint: NPDE1.
2:     TS – REAL (KIND=nag_wp)Input/Output
On entry: the initial value of the independent variable t.
Constraint: TS<TOUT.
On exit: the value of t corresponding to the solution values in U. Normally 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 user.External Procedure
PDEDEF must compute the functions Gi which define the system of PDEs. PDEDEF is called approximately midway between each pair of mesh points in turn by D03PEF.
The specification of PDEDEF is:
SUBROUTINE PDEDEF ( NPDE, T, X, U, UT, UX, RES, IRES)
INTEGER  NPDE, IRES
REAL (KIND=nag_wp)  T, X, U(NPDE), UT(NPDE), UX(NPDE), RES(NPDE)
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:     U(NPDE) – REAL (KIND=nag_wp) arrayInput
On entry: Ui contains the value of the component Uix,t, for i=1,2,,NPDE.
5:     UT(NPDE) – REAL (KIND=nag_wp) arrayInput
On entry: UTi contains the value of the component Uix,t t , for i=1,2,,NPDE.
6:     UX(NPDE) – REAL (KIND=nag_wp) arrayInput
On entry: UXi contains the value of the component Uix,t x , for i=1,2,,NPDE.
7:     RES(NPDE) – REAL (KIND=nag_wp) arrayOutput
On exit: RESi must contain the ith component of G, for i=1,2,,NPDE, where G is defined as
Gi=j=1NPDEPi,j Uj t , (8)
i.e., only terms depending explicitly on time derivatives, or
Gi=j=1NPDEPi,j Uj t +Qi, (9)
i.e., all terms in equation (2).
The definition of G is determined by the input value of IRES.
8:     IRES – INTEGERInput/Output
On entry: the form of Gi that must be returned in the array RES.
IRES=-1
Equation (8) must be used.
IRES=1
Equation (9) must be used.
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 (sub)routine 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, then D03PEF 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 D03PEF is called. Parameters denoted as Input must not be changed by this procedure.
5:     BNDARY – SUBROUTINE, supplied by the user.External Procedure
BNDARY must compute the functions GiL and GiR which define the boundary conditions as in equations (4) and (5).
The specification of BNDARY is:
SUBROUTINE BNDARY ( NPDE, T, IBND, NOBC, U, UT, RES, IRES)
INTEGER  NPDE, IBND, NOBC, IRES
REAL (KIND=nag_wp)  T, U(NPDE), UT(NPDE), RES(NOBC)
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:     IBND – INTEGERInput
On entry: determines the position of the boundary conditions.
IBND=0
BNDARY must compute the left-hand boundary condition at x=a.
IBND0
Indicates that BNDARY must compute the right-hand boundary condition at x=b.
4:     NOBC – INTEGERInput
On entry: specifies the number of boundary conditions at the boundary specified by IBND.
5:     U(NPDE) – REAL (KIND=nag_wp) arrayInput
On entry: Ui contains the value of the component Uix,t at the boundary specified by IBND, for i=1,2,,NPDE.
6:     UT(NPDE) – REAL (KIND=nag_wp) arrayInput
On entry: UTi contains the value of the component Uix,t t  at the boundary specified by IBND, for i=1,2,,NPDE.
7:     RES(NOBC) – REAL (KIND=nag_wp) arrayOutput
On exit: RESi must contain the ith component of GL or GR, depending on the value of IBND, for i=1,2,,NOBC, where GL is defined as
GiL=j=1NPDE Ei,jL Uj t , (10)
i.e., only terms depending explicitly on time derivatives, or
GiL=j=1NPDE Ei,jL Uj t +SiL, (11)
i.e., all terms in equation (6), and similarly for GiR.
The definitions of GL and GR are determined by the input value of IRES.
8:     IRES – INTEGERInput/Output
On entry: the form GiL (or GiR) that must be returned in the array RES.
IRES=-1
Equation (10) must be used.
IRES=1
Equation (11) must be used.
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 (sub)routine 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, then D03PEF 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 D03PEF is called. Parameters denoted as Input must not be changed by this procedure.
6:     U(NPDE,NPTS) – REAL (KIND=nag_wp) arrayInput/Output
On entry: the initial values of Ux,t at t=TS and the mesh points Xj, for j=1,2,,NPTS.
On exit: Uij will contain the computed solution at t=TS.
7:     NPTS – INTEGERInput
On entry: the number of mesh points in the interval a,b.
Constraint: NPTS3.
8:     X(NPTS) – REAL (KIND=nag_wp) arrayInput
On entry: the mesh points in the spatial direction. X1 must specify the left-hand boundary, a, and XNPTS must specify the right-hand boundary, b.
Constraint: X1<X2<<XNPTS.
9:     NLEFT – INTEGERInput
On entry: the number na of boundary conditions at the left-hand mesh point X1.
Constraint: 0NLEFTNPDE.
10:   ACC – REAL (KIND=nag_wp)Input
On entry: a positive quantity for controlling the local error estimate in the time integration. If Ei,j is the estimated error for Ui at the jth mesh point, the error test is:
Ei,j=ACC×1.0+Uij.
Constraint: ACC>0.0.
11:   RSAVE(LRSAVE) – 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.
12:   LRSAVE – INTEGERInput
On entry: the dimension of the array RSAVE as declared in the (sub)program from which D03PEF is called.
Constraint: LRSAVE4×NPDE+NLEFT+14×NPDE×NPTS+3×NPDE+21×NPDE+7×NPTS+54.
13:   ISAVE(LISAVE) – 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.
14:   LISAVE – INTEGERInput
On entry: the dimension of the array ISAVE as declared in the (sub)program from which D03PEF is called.
Constraint: LISAVENPDE×NPTS+24.
15:   ITASK – INTEGERInput
On entry: specifies the task to be performed by the ODE integrator.
ITASK=1
Normal computation of output values U at t=TOUT.
ITASK=2
Take one step and return.
ITASK=3
Stop at the first internal integration point at or beyond t=TOUT.
Constraint: ITASK=1, 2 or 3.
16:   ITRACE – INTEGERInput
On entry: the level of trace information required from D03PEF and the underlying ODE solver as follows:
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=1
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.
ITRACE=2
Output from the underlying ODE solver is similar to that produced when ITRACE=1, except that the advisory messages are given in greater detail.
ITRACE3
Output from the underlying ODE solver is similar to that produced when ITRACE=2, except that the advisory messages are given in greater detail.
You are advised to set ITRACE=0, unless you are experienced with sub-chapter D02M–N.
17:   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 parameters TOUT and IFAIL should be reset between calls to D03PEF.
Constraint: IND=0 or 1.
On exit: IND=1.
18:   IFAIL – INTEGERInput/Output
On entry: IFAIL must be set to 0, -1​ or ​1. If you are unfamiliar with this parameter you should refer to Section 3.3 in the Essential Introduction 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 parameter, 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,TOUTTS,
orTOUT-TS is too small,
orITASK1, 2 or 3,
ormesh points ​Xi are not ordered correctly,
orNPTS<3,
orNPDE<1,
orNLEFT is not in the range 0 to NPDE,
orACC0.0,
orIND0 or 1,
orLRSAVE is too small,
orLISAVE is too small,
orD03PEF called initially with IND=1.
IFAIL=2
The underlying ODE solver cannot make any further progress across the integration range from the current point t=TS with the supplied value of ACC. The components of U contain the computed values at the current point t=TS.
IFAIL=3
In the underlying ODE solver, there were repeated errors or corrector convergence test failures on an attempted step, before completing the requested task. The problem may have a singularity or ACC is too small for the integration to continue. Incorrect positioning of boundary conditions may also result in this error. Integration was successful as far as t=TS.
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 the PDEDEF or BNDARY, when the residual in the underlying ODE solver was being evaluated. Incorrect positioning of boundary conditions may also result in this error.
IFAIL=5
In solving the ODE system, a singular Jacobian has been encountered. You should check their problem formulation.
IFAIL=6
When evaluating the residual in solving the ODE system, IRES was set to 2 in one of PDEDEF or BNDARY. Integration was successful as far as t=TS.
IFAIL=7
The value of ACC is so small that the routine is unable to start the integration in time.
IFAIL=8
In either, PDEDEF 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 parameters 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 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).

7  Accuracy

D03PEF 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 accuracy parameter, ACC.

8  Further Comments

The Keller box scheme can be used to solve higher-order problems which have been reduced to first-order by the introduction of new variables (see the example problem in D03PKF). In general, a second-order problem can be solved with slightly greater accuracy using the Keller box scheme instead of a finite difference scheme (D03PCF/D03PCA or D03PHF/D03PHA for example), but at the expense of increased CPU time due to the larger number of function evaluations required.
It should be noted that the Keller box scheme, in common with other central-difference schemes, may be unsuitable for some hyperbolic first-order problems such as the apparently simple linear advection equation Ut+aUx=0, where a is a constant, resulting in spurious oscillations due to the lack of dissipation. This type of problem requires a discretization scheme with upwind weighting (D03PFF for example), or the addition of a second-order artificial dissipation term.
The time taken depends on the complexity of the system and on the accuracy requested.

9  Example

This example is the simple first-order system
U1 t + U1 x + U2 x =0, U2 t +4 U1 x + U2 x =0,
for t0,1 and x0,1.
The initial conditions are
U1x,0=expx,  U2x,0=sinx,
and the Dirichlet boundary conditions for U1 at x=0 and U2 at x=1 are given by the exact solution:
U1x,t=12 expx+t+expx-3t+14 sinx-3t-sinx+t , U2x,t=expx-3t-expx+t+12 sinx+t+sinx-3t .

9.1  Program Text

Program Text (d03pefe.f90)

9.2  Program Data

Program Data (d03pefe.d)

9.3  Program Results

Program Results (d03pefe.r)

Produced by GNUPLOT 4.4 patchlevel 0 Example Program Solution, U(1,x,t), of First-order System using Keller, Box and BDF U(1,x,t) 0.1 0.2 0.3 0.4 0.5 0.6 0.7 0.8 0.9 1 Time 0 0.2 0.4 0.6 0.8 1 x 0.5 1 1.5 2 2.5 3 3.5
Produced by GNUPLOT 4.4 patchlevel 0 Solution, U(2,x,t), of First-order System using Keller, Box and BDF U(2,x,t) 0.1 0.2 0.3 0.4 0.5 0.6 0.7 0.8 0.9 1 Time 0 0.2 0.4 0.6 0.8 1 x -8 -7 -6 -5 -4 -3 -2 -1 0

D03PEF (PDF version)
D03 Chapter Contents
D03 Chapter Introduction
NAG Library Manual

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