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

NAG Library Routine Document

D03PRF

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

D03PRF integrates a system of linear or nonlinear, first-order, time-dependent partial differential equations (PDEs) in one space variable, with scope for coupled ordinary differential equations (ODEs), and automatic adaptive spatial remeshing. The spatial discretization is performed using the Keller box scheme (see Keller (1970)) and the method of lines is employed to reduce the PDEs to a system of ODEs. The resulting system is solved using a Backward Differentiation Formula (BDF) method or a Theta method (switching between Newton's method and functional iteration).

2  Specification

INTEGER  NPDE, NPTS, NLEFT, NCODE, NXI, NEQN, ITOL, NXFIX, NRMESH, IPMINF, LRSAVE, ISAVE(LISAVE), LISAVE, ITASK, ITRACE, IND, IFAIL
REAL (KIND=nag_wp)  TS, TOUT, U(NEQN), X(NPTS), XI(*), RTOL(*), ATOL(*), ALGOPT(30), XFIX(*), DXMESH, TRMESH, XRATIO, CON, RSAVE(LRSAVE)
LOGICAL  REMESH
CHARACTER(1)  NORM, LAOPT
EXTERNAL  PDEDEF, BNDARY, UVINIT, ODEDEF, MONITF

3  Description

D03PRF integrates the system of first-order PDEs and coupled ODEs given by the master equations:
Gi x,t,U,Ux,Ut,V,V. = 0 ,   i=1,2,,NPDE ,   axb,tt0, (1)
Rit,V,V.,ξ,U*,Ux*,Ut*=0,  i=1,2,,NCODE. (2)
In the PDE part of the problem given by (1), the functions Gi must have the general form
Gi=j=1NPDEPi,j Uj t +j=1NCODEQi,jV.j+Si=0,  i=1,2,,NPDE, (3)
where Pi,j, Qi,j and Si depend on x, t, U, Ux and V.
The vector U is the set of PDE solution values
Ux,t=U1x,t,,UNPDEx,tT,
and the vector Ux is the partial derivative with respect to x. The vector V is the set of ODE solution values
Vt=V1t,,VNCODEtT,
and V. denotes its derivative with respect to time.
In the ODE part given by (2), ξ represents a vector of nξ spatial coupling points at which the ODEs are coupled to the PDEs. These points may or may not be equal to some of the PDE spatial mesh points. U*, Ux* and Ut* are the functions U, Ux and Ut evaluated at these coupling points. Each Ri may only depend linearly on time derivatives. Hence equation (2) may be written more precisely as
R=A-BV.-CUt*, (4)
where R=R1,,RNCODET, A is a vector of length NCODE, B is an NCODE by NCODE matrix, C is an NCODE by nξ×NPDE matrix and the entries in A, B and C may depend on t, ξ, U*, Ux* and V. In practice you only need to supply a vector of information to define the ODEs and not the matrices B and C. (See Section 5 for the specification of ODEDEF.)
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 mesh x1,x2,,xNPTS defined initially by you and (possibly) adapted automatically during the integration according to user-specified criteria.
The PDE system which is defined by the functions Gi must be specified in PDEDEF.
The initial t=t0 values of the functions Ux,t and Vt must be specified in UVINIT. Note that UVINIT will be called again following any remeshing, and so Ux,t0 should be specified for all values of x in the interval axb, and not just the initial mesh points.
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 master equation form:
GiL x,t,U,Ut,V,V. = 0   at ​ x = a ,   i=1,2,,na , (5)
at the left-hand boundary, and
GiR x,t,U,Ut,V,V. = 0   at ​ x = b ,   i=1,2,,nb , (6)
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 routines. 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=1 NPDE E i,j L Uj t + j=1 NCODE H i,j L V.j + KiL = 0 ,   i=1,2,,na , (7)
at the left-hand boundary, and
j=1 NPDE E i,j R Uj t + j=1 NCODE H i,j R V.j + KiR = 0 ,   i=1,2,,nb , (8)
at the right-hand boundary, where Ei,jL, Ei,jR, Hi,jL, Hi,jR, KiL and KiR depend on x,t,U and V only.
The boundary conditions must be specified in BNDARY.
The problem is subject to the following restrictions:
(i) Pi,j, Qi,j and Si must not depend on any time derivatives;
(ii) t0<tout, so that integration is in the forward direction;
(iii) The evaluation of the function Gi is done approximately at the mid-points of the mesh Xi, for i=1,2,,NPTS, by calling PDEDEF for each mid-point in turn. Any discontinuities in the function must therefore be at one or more of the fixed mesh points specified by XFIX;
(iv) At least one of the functions Pi,j must be nonzero so that there is a time derivative present in the PDE problem.
The algebraic-differential equation system which is defined by the functions Ri must be specified in ODEDEF. You must also specify the coupling points ξ in the array XI.
The first-order equations are approximated by a system of ODEs in time for the values of Ui at mesh points. In this method of lines approach the Keller box scheme 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+NCODE ODEs in time direction. This system is then integrated forwards in time using a Backward Differentiation Formula (BDF) or a Theta method.
The adaptive space remeshing can be used to generate meshes that automatically follow the changing time-dependent nature of the solution, generally resulting in a more efficient and accurate solution using fewer mesh points than may be necessary with a fixed uniform or non-uniform mesh. Problems with travelling wavefronts or variable-width boundary layers for example will benefit from using a moving adaptive mesh. The discrete time-step method used here (developed by Furzeland (1984)) automatically creates a new mesh based on the current solution profile at certain time-steps, and the solution is then interpolated onto the new mesh and the integration continues.
The method requires you to supply MONITF which specifies in an analytic or numeric form the particular aspect of the solution behaviour you wish to track. This so-called monitor function is used to choose a mesh which equally distributes the integral of the monitor function over the domain. A typical choice of monitor function is the second space derivative of the solution value at each point (or some combination of the second space derivatives if more than one solution component), which results in refinement in regions where the solution gradient is changing most rapidly.
You must specify the frequency of mesh updates along with certain other criteria such as adjacent mesh ratios. Remeshing can be expensive and you are encouraged to experiment with the different options in order to achieve an efficient solution which adequately tracks the desired features of the solution.
Note that unless the monitor function for the initial solution values is zero at all user-specified initial mesh points, a new initial mesh is calculated and adopted according to the user-specified remeshing criteria. UVINIT will then be called again to determine the initial solution values at the new mesh points (there is no interpolation at this stage) and the integration proceeds.

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
Berzins M and Furzeland R M (1992) An adaptive theta method for the solution of stiff and nonstiff differential equations Appl. Numer. Math. 9 1–19
Furzeland R M (1984) The construction of adaptive space meshes TNER.85.022 Thornton Research Centre, Chester
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 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 evaluate the functions Gi which define the system of PDEs. PDEDEF is called approximately midway between each pair of mesh points in turn by D03PRF.
The specification of PDEDEF is:
SUBROUTINE PDEDEF ( NPDE, T, X, U, UDOT, UX, NCODE, V, VDOT, RES, IRES)
INTEGER  NPDE, NCODE, IRES
REAL (KIND=nag_wp)  T, X, U(NPDE), UDOT(NPDE), UX(NPDE), V(NCODE), VDOT(NCODE), 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:     UDOT(NPDE) – REAL (KIND=nag_wp) arrayInput
On entry: UDOTi 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:     NCODE – INTEGERInput
On entry: the number of coupled ODEs in the system.
8:     V(NCODE) – REAL (KIND=nag_wp) arrayInput
On entry: if NCODE>0, Vi contains the value of the component Vit, for i=1,2,,NCODE.
9:     VDOT(NCODE) – REAL (KIND=nag_wp) arrayInput
On entry: if NCODE>0, VDOTi contains the value of component V.it, for i=1,2,,NCODE.
10:   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 +j=1NCODEQi,jV.j, (9)
i.e., only terms depending explicitly on time derivatives, or
Gi=j=1NPDEPi,j Uj t +j=1NCODEQi,jV.j+Si, (10)
i.e., all terms in equation (3).
The definition of G is determined by the input value of IRES.
11:   IRES – INTEGERInput/Output
On entry: the form of Gi that must be returned in the array RES.
IRES=-1
Equation (9) must be used.
IRES=1
Equation (10) 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 D03PRF 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 D03PRF is called. Parameters denoted as Input must not be changed by this procedure.
5:     BNDARY – SUBROUTINE, supplied by the user.External Procedure
BNDARY must evaluate the functions GiL and GiR which describe the boundary conditions, as given in (5) and (6).
The specification of BNDARY is:
SUBROUTINE BNDARY ( NPDE, T, IBND, NOBC, U, UDOT, NCODE, V, VDOT, RES, IRES)
INTEGER  NPDE, IBND, NOBC, NCODE, IRES
REAL (KIND=nag_wp)  T, U(NPDE), UDOT(NPDE), V(NCODE), VDOT(NCODE), 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: specifies which boundary conditions are to be evaluated.
IBND=0
BNDARY must compute the left-hand boundary condition at x=a.
IBND0
BNDARY must compute of the right-hand boundary condition at x=b.
4:     NOBC – INTEGERInput
On entry: specifies the number na 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:     UDOT(NPDE) – REAL (KIND=nag_wp) arrayInput
On entry: UDOTi contains the value of the component Uix,t t , for i=1,2,,NPDE.
7:     NCODE – INTEGERInput
On entry: the number of coupled ODEs in the system.
8:     V(NCODE) – REAL (KIND=nag_wp) arrayInput
On entry: if NCODE>0, Vi contains the value of the component Vit, for i=1,2,,NCODE.
9:     VDOT(NCODE) – REAL (KIND=nag_wp) arrayInput
On entry: if NCODE>0, VDOTi contains the value of component V.it, for i=1,2,,NCODE.
Note: VDOTi, for i=1,2,,NCODE, may only appear linearly as in (11) and (12).
10:   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=1NPDEEi,jL Uj t +j=1NCODEHi,jLV.j, (11)
i.e., only terms depending explicitly on time derivatives, or
GiL=j=1NPDEEi,jL Uj t +j=1NCODEHi,jLV.j+KiL, (12)
i.e., all terms in equation (7), and similarly for GiR.
The definitions of GL and GR are determined by the input value of IRES.
11:   IRES – INTEGERInput/Output
On entry: the form of GiL (or GiR) that must be returned in the array RES.
IRES=-1
Equation (11) must be used.
IRES=1
Equation (12) 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 D03PRF 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 D03PRF is called. Parameters denoted as Input must not be changed by this procedure.
6:     UVINIT – SUBROUTINE, supplied by the user.External Procedure
UVINIT must supply the initial t=t0 values of Ux,t and Vt for all values of x in the interval a,b.
The specification of UVINIT is:
SUBROUTINE UVINIT ( NPDE, NPTS, NXI, X, XI, U, NCODE, V)
INTEGER  NPDE, NPTS, NXI, NCODE
REAL (KIND=nag_wp)  X(NPTS), XI(NXI), U(NPDE,NPTS), V(NCODE)
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:     NXI – INTEGERInput
On entry: the number of ODE/PDE coupling points.
4:     X(NPTS) – REAL (KIND=nag_wp) arrayInput
On entry: the current mesh. Xi contains the value of xi, for i=1,2,,NPTS.
5:     XI(NXI) – REAL (KIND=nag_wp) arrayInput
On entry: if NXI>0, XIi contains the ODE/PDE coupling point, ξi, for i=1,2,,NXI.
6:     U(NPDE,NPTS) – REAL (KIND=nag_wp) arrayOutput
On exit: if NXI>0, Uij contains the value of the component Uixj,t0, for i=1,2,,NPDE and j=1,2,,NPTS.
7:     NCODE – INTEGERInput
On entry: the number of coupled ODEs in the system.
8:     V(NCODE) – REAL (KIND=nag_wp) arrayOutput
On exit: if NCODE>0, Vi must contain the value of component Vit0, for i=1,2,,NCODE.
UVINIT must either be a module subprogram USEd by, or declared as EXTERNAL in, the (sub)program from which D03PRF is called. Parameters denoted as Input must not be changed by this procedure.
7:     U(NEQN) – REAL (KIND=nag_wp) arrayInput/Output
On entry: if IND=1, the value of U must be unchanged from the previous call.
On exit: UNPDE×j-1+i contains the computed solution Uixj,t, for i=1,2,,NPDE and j=1,2,,NPTS, evaluated at t=TS.
8:     NPTS – INTEGERInput
On entry: the number of mesh points in the interval [a,b].
Constraint: NPTS3.
9:     X(NPTS) – REAL (KIND=nag_wp) arrayInput/Output
On entry: the initial 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.
On exit: the final values of the mesh points.
10:   NLEFT – INTEGERInput
On entry: the number na of boundary conditions at the left-hand mesh point X1.
Constraint: 0NLEFTNPDE.
11:   NCODE – INTEGERInput
On entry: the number of coupled ODE components.
Constraint: NCODE0.
12:   ODEDEF – SUBROUTINE, supplied by the NAG Library or the user.External Procedure
ODEDEF must evaluate the functions R, which define the system of ODEs, as given in (4).
If you wish to compute the solution of a system of PDEs only (i.e., NCODE=0), ODEDEF must be the dummy routine D03PEK. (D03PEK is included in the NAG Library.)
The specification of ODEDEF is:
SUBROUTINE ODEDEF ( NPDE, T, NCODE, V, VDOT, NXI, XI, UCP, UCPX, UCPT, R, IRES)
INTEGER  NPDE, NCODE, NXI, IRES
REAL (KIND=nag_wp)  T, V(NCODE), VDOT(NCODE), XI(NXI), UCP(NPDE,*), UCPX(NPDE,*), UCPT(NPDE,*), R(NCODE)
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:     NCODE – INTEGERInput
On entry: the number of coupled ODEs in the system.
4:     V(NCODE) – REAL (KIND=nag_wp) arrayInput
On entry: if NCODE>0, Vi contains the value of the component Vit, for i=1,2,,NCODE.
5:     VDOT(NCODE) – REAL (KIND=nag_wp) arrayInput
On entry: if NCODE>0, VDOTi contains the value of component V.it, for i=1,2,,NCODE.
6:     NXI – INTEGERInput
On entry: the number of ODE/PDE coupling points.
7:     XI(NXI) – REAL (KIND=nag_wp) arrayInput
On entry: if NXI>0, XIi contains the ODE/PDE coupling point, ξi, for i=1,2,,NXI.
8:     UCP(NPDE,*) – REAL (KIND=nag_wp) arrayInput
On entry: if NXI>0, UCPij contains the value of Uix,t at the coupling point x=ξj, for i=1,2,,NPDE and j=1,2,,NXI.
9:     UCPX(NPDE,*) – REAL (KIND=nag_wp) arrayInput
On entry: if NXI>0, UCPXij contains the value of Uix,t x  at the coupling point x=ξj, for i=1,2,,NPDE and j=1,2,,NXI.
10:   UCPT(NPDE,*) – REAL (KIND=nag_wp) arrayInput
On entry: if NXI>0, UCPTij contains the value of Ui t  at the coupling point x=ξj, for i=1,2,,NPDE and j=1,2,,NXI.
11:   R(NCODE) – REAL (KIND=nag_wp) arrayOutput
On exit: if NCODE>0, Ri must contain the ith component of R, for i=1,2,,NCODE, where R is defined as
R=-BV.-CUt*, (13)
i.e., only terms depending explicitly on time derivatives, or
R=A-BV.-CUt*, (14)
i.e., all terms in equation (4). The definition of R is determined by the input value of IRES.
12:   IRES – INTEGERInput/Output
On entry: the form of R that must be returned in the array R.
IRES=-1
Equation (13) must be used.
IRES=1
Equation (14) must be used.
On exit: should usually remain unchanged. However, you may reset 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 D03PRF returns to the calling subroutine with the error indicator set to IFAIL=4.
ODEDEF must either be a module subprogram USEd by, or declared as EXTERNAL in, the (sub)program from which D03PRF is called. Parameters denoted as Input must not be changed by this procedure.
13:   NXI – INTEGERInput
On entry: the number of ODE/PDE coupling points.
Constraints:
  • if NCODE=0, NXI=0;
  • if NCODE>0, NXI0.
14:   XI(*) – REAL (KIND=nag_wp) arrayInput
Note: the dimension of the array XI must be at least max1,NXI.
On entry: XIi, for i=1,2,,NXI, must be set to the ODE/PDE coupling points, ξi.
Constraint: X1XI1<XI2<<XINXIXNPTS.
15:   NEQN – INTEGERInput
On entry: the number of ODEs in the time direction.
Constraint: NEQN=NPDE×NPTS+NCODE.
16:   RTOL(*) – REAL (KIND=nag_wp) arrayInput
Note: the dimension of the array RTOL must be at least 1 if ITOL=1 or 2 and at least NEQN if ITOL=3 or 4.
On entry: the relative local error tolerance.
Constraint: RTOLi0.0 for all relevant i.
17:   ATOL(*) – REAL (KIND=nag_wp) arrayInput
Note: the dimension of the array ATOL must be at least 1 if ITOL=1 or 3 and at least NEQN if ITOL=2 or 4.
On entry: the absolute local error tolerance.
Constraint: ATOLi0.0 for all relevant i.
Note: corresponding elements of RTOL and ATOL cannot both be 0.0.
18:   ITOL – INTEGERInput
A value to indicate the form of the local error test. ITOL indicates to D03PRF whether to interpret either or both of RTOL or ATOL as a vector or scalar. The error test to be satisfied is ei/wi<1.0, where wi is defined as follows:
On entry:
ITOL RTOL ATOL wi
1 scalar scalar RTOL1×Ui+ATOL1
2 scalar vector RTOL1×Ui+ATOLi
3 vector scalar RTOLi×Ui+ATOL1
4 vector vector RTOLi×Ui+ATOLi
In the above, ei denotes the estimated local error for the ith component of the coupled PDE/ODE system in time, Ui, for i=1,2,,NEQN.
The choice of norm used is defined by the parameter NORM.
Constraint: ITOL=1, 2, 3 or 4.
19:   NORM – CHARACTER(1)Input
On entry: the type of norm to be used.
NORM='M'
Maximum norm.
NORM='A'
Averaged L2 norm.
If Unorm denotes the norm of the vector U of length NEQN, then for the averaged L2 norm
Unorm= 1NEQN i=1NEQN Ui/wi 2 ,
while for the maximum norm
Unorm=maxiUi/wi.
See the description of ITOL for the formulation of the weight vector w.
Constraint: NORM='M' or 'A'.
20:   LAOPT – CHARACTER(1)Input
On entry: the type of matrix algebra required.
LAOPT='F'
Full matrix methods to be used.
LAOPT='B'
Banded matrix methods to be used.
LAOPT='S'
Sparse matrix methods to be used.
Constraint: LAOPT='F', 'B' or 'S'.
Note: you are recommended to use the banded option when no coupled ODEs are present (i.e., NCODE=0).
21:   ALGOPT(30) – REAL (KIND=nag_wp) arrayInput
On entry: may be set to control various options available in the integrator. If you wish to employ all the default options, then ALGOPT1 should be set to 0.0. Default values will also be used for any other elements of ALGOPT set to zero. The permissible values, default values, and meanings are as follows:
ALGOPT1
Selects the ODE integration method to be used. If ALGOPT1=1.0, a BDF method is used and if ALGOPT1=2.0, a Theta method is used. The default value is ALGOPT1=1.0.
If ALGOPT1=2.0, then ALGOPTi, for i=2,3,4, are not used.
ALGOPT2
Specifies the maximum order of the BDF integration formula to be used. ALGOPT2 may be 1.0, 2.0, 3.0, 4.0 or 5.0. The default value is ALGOPT2=5.0.
ALGOPT3
Specifies what method is to be used to solve the system of nonlinear equations arising on each step of the BDF method. If ALGOPT3=1.0 a modified Newton iteration is used and if ALGOPT3=2.0 a functional iteration method is used. If functional iteration is selected and the integrator encounters difficulty, then there is an automatic switch to the modified Newton iteration. The default value is ALGOPT3=1.0.
ALGOPT4
Specifies whether or not the Petzold error test is to be employed. The Petzold error test results in extra overhead but is more suitable when algebraic equations are present, such as Pi,j=0.0, for j=1,2,,NPDE, for some i or when there is no V.it dependence in the coupled ODE system. If ALGOPT4=1.0, then the Petzold test is used. If ALGOPT4=2.0, then the Petzold test is not used. The default value is ALGOPT4=1.0.
If ALGOPT1=1.0, then ALGOPTi, for i=5,6,7, are not used.
ALGOPT5
Specifies the value of Theta to be used in the Theta integration method. 0.51ALGOPT50.99. The default value is ALGOPT5=0.55.
ALGOPT6
Specifies what method is to be used to solve the system of nonlinear equations arising on each step of the Theta method. If ALGOPT6=1.0, a modified Newton iteration is used and if ALGOPT6=2.0, a functional iteration method is used. The default value is ALGOPT6=1.0.
ALGOPT7
Specifies whether or not the integrator is allowed to switch automatically between modified Newton and functional iteration methods in order to be more efficient. If ALGOPT7=1.0, then switching is allowed and if ALGOPT7=2.0, then switching is not allowed. The default value is ALGOPT7=1.0.
ALGOPT11
Specifies a point in the time direction, tcrit, beyond which integration must not be attempted. The use of tcrit is described under the parameter ITASK. If ALGOPT10.0, a value of 0.0, for ALGOPT11, say, should be specified even if ITASK subsequently specifies that tcrit will not be used.
ALGOPT12
Specifies the minimum absolute step size to be allowed in the time integration. If this option is not required, ALGOPT12 should be set to 0.0.
ALGOPT13
Specifies the maximum absolute step size to be allowed in the time integration. If this option is not required, ALGOPT13 should be set to 0.0.
ALGOPT14
Specifies the initial step size to be attempted by the integrator. If ALGOPT14=0.0, then the initial step size is calculated internally.
ALGOPT15
Specifies the maximum number of steps to be attempted by the integrator in any one call. If ALGOPT15=0.0, then no limit is imposed.
ALGOPT23
Specifies what method is to be used to solve the nonlinear equations at the initial point to initialize the values of U, Ut, V and V.. If ALGOPT23=1.0, a modified Newton iteration is used and if ALGOPT23=2.0, functional iteration is used. The default value is ALGOPT23=1.0.
ALGOPT29 and ALGOPT30 are used only for the sparse matrix algebra option, i.e., LAOPT='S'.
ALGOPT29
Governs the choice of pivots during the decomposition of the first Jacobian matrix. It should lie in the range 0.0<ALGOPT29<1.0, with smaller values biasing the algorithm towards maintaining sparsity at the expense of numerical stability. If ALGOPT29 lies outside this range then the default value is used. If the routines regard the Jacobian matrix as numerically singular then increasing ALGOPT29 towards 1.0 may help, but at the cost of increased fill-in. The default value is ALGOPT29=0.1.
ALGOPT30
Used as a relative pivot threshold during subsequent Jacobian decompositions (see ALGOPT29) below which an internal error is invoked. ALGOPT30 must be greater than zero, otherwise the default value is used. If ALGOPT30 is greater than 1.0 no check is made on the pivot size, and this may be a necessary option if the Jacobian is found to be numerically singular (see ALGOPT29). The default value is ALGOPT30=0.0001.
22:   REMESH – LOGICALInput
On entry: indicates whether or not spatial remeshing should be performed.
REMESH=.TRUE.
Indicates that spatial remeshing should be performed as specified.
REMESH=.FALSE.
Indicates that spatial remeshing should be suppressed.
Note:  REMESH should not be changed between consecutive calls to D03PRF. Remeshing can be switched off or on at specified times by using appropriate values for the parameters NRMESH and TRMESH at each call.
23:   NXFIX – INTEGERInput
On entry: the number of fixed mesh points.
Constraint: 0NXFIXNPTS-2.
Note: the end points X1 and XNPTS are fixed automatically and hence should not be specified as fixed points.
24:   XFIX(*) – REAL (KIND=nag_wp) arrayInput
Note: the dimension of the array XFIX must be at least max1,NXFIX.
On entry: XFIXi, for i=1,2,,NXFIX, must contain the value of the x coordinate at the ith fixed mesh point.
Constraint: XFIXi<XFIXi+1, for i=1,2,,NXFIX-1, and each fixed mesh point must coincide with a user-supplied initial mesh point, that is XFIXi=Xj for some j, 2jNPTS-1.
Note: the positions of the fixed mesh points in the array X remain fixed during remeshing, and so the number of mesh points between adjacent fixed points (or between fixed points and end points) does not change. You should take this into account when choosing the initial mesh distribution.
25:   NRMESH – INTEGERInput
On entry: indicates the form of meshing to be performed.
NRMESH<0
Indicates that a new mesh is adopted according to the parameter DXMESH. The mesh is tested every NRMESH timesteps.
NRMESH=0
Indicates that remeshing should take place just once at the end of the first time step reached when t>TRMESH.
NRMESH>0
Indicates that remeshing will take place every NRMESH time steps, with no testing using DXMESH.
Note: NRMESH may be changed between consecutive calls to D03PRF to give greater flexibility over the times of remeshing.
26:   DXMESH – REAL (KIND=nag_wp)Input
On entry: determines whether a new mesh is adopted when NRMESH is set less than zero. A possible new mesh is calculated at the end of every NRMESH time steps, but is adopted only if
xinew>xiold+DXMESH×xi+1old-xiold,
or
xinew<xiold-DXMESH×xiold-xi- 1old.
DXMESH thus imposes a lower limit on the difference between one mesh and the next.
Constraint: DXMESH0.0.
27:   TRMESH – REAL (KIND=nag_wp)Input
On entry: specifies when remeshing will take place when NRMESH is set to zero. Remeshing will occur just once at the end of the first time step reached when t is greater than TRMESH.
Note: TRMESH may be changed between consecutive calls to D03PRF to force remeshing at several specified times.
28:   IPMINF – INTEGERInput
On entry: the level of trace information regarding the adaptive remeshing. Details are directed to the current advisory message unit (see X04ABF).
IPMINF=0
No trace information.
IPMINF=1
Brief summary of mesh characteristics.
IPMINF=2
More detailed information, including old and new mesh points, mesh sizes and monitor function values.
Constraint: IPMINF=0, 1 or 2.
29:   XRATIO – REAL (KIND=nag_wp)Input
On entry: input bound on adjacent mesh ratio (greater than 1.0 and typically in the range 1.5 to 3.0). The remeshing routines will attempt to ensure that
xi-xi-1/XRATIO<xi+1-xi<XRATIO×xi-xi-1.
Suggested value: XRATIO=1.5.
Constraint: XRATIO>1.0.
30:   CON – REAL (KIND=nag_wp)Input
On entry: an input bound on the sub-integral of the monitor function Fmonx over each space step. The remeshing routines will attempt to ensure that
x1xi+1FmonxdxCONx1xNPTSFmonxdx,
(see Furzeland (1984)). CON gives you more control over the mesh distribution e.g., decreasing CON allows more clustering. A typical value is 2/NPTS-1, but you are encouraged to experiment with different values. Its value is not critical and the mesh should be qualitatively correct for all values in the range given below.
Suggested value: CON=2.0/NPTS-1.
Constraint: 0.1/NPTS-1CON10.0/NPTS-1.
31:   MONITF – SUBROUTINE, supplied by the NAG Library or the user.External Procedure
MONITF must supply and evaluate a remesh monitor function to indicate the solution behaviour of interest.
If you specify REMESH=.FALSE., i.e., no remeshing, then MONITF will not be called and the dummy routine D03PEL may be used for MONITF. (D03PEL is included in the NAG Library.)
The specification of MONITF is:
SUBROUTINE MONITF ( T, NPTS, NPDE, X, U, FMON)
INTEGER  NPTS, NPDE
REAL (KIND=nag_wp)  T, X(NPTS), U(NPDE,NPTS), FMON(NPTS)
1:     T – REAL (KIND=nag_wp)Input
On entry: the current value of the independent variable t.
2:     NPTS – INTEGERInput
On entry: the number of mesh points in the interval a,b.
3:     NPDE – INTEGERInput
On entry: the number of PDEs in the system.
4:     X(NPTS) – REAL (KIND=nag_wp) arrayInput
On entry: the current mesh. Xi contains the value of xi, for i=1,2,,NPTS.
5:     U(NPDE,NPTS) – REAL (KIND=nag_wp) arrayInput
On entry: Uij contains the value of Uix,t at x=Xj and time t, for i=1,2,,NPDE and j=1,2,,NPTS.
6:     FMON(NPTS) – REAL (KIND=nag_wp) arrayOutput
On exit: FMONi must contain the value of the monitor function Fmonx at mesh point x=Xi.
Constraint: FMONi0.0.
MONITF must either be a module subprogram USEd by, or declared as EXTERNAL in, the (sub)program from which D03PRF is called. Parameters denoted as Input must not be changed by this procedure.
32:   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.
33:   LRSAVE – INTEGERInput
On entry: the dimension of the array RSAVE as declared in the (sub)program from which D03PRF is called. Its size depends on the type of matrix algebra selected.
If LAOPT='F', LRSAVENEQN×NEQN+NEQN+nwkres+lenode.
If LAOPT='B', LRSAVE3×ml+mu+2×NEQN+nwkres+lenode.
If LAOPT='S', LRSAVE4×NEQN+11×NEQN/2+1+nwkres+lenode.
Where
ml and mu are the lower and upper half bandwidths given by ml=NPDE+NLEFT-1 such that
mu=2×NPDE-NLEFT-1, for problems involving PDEs only; or
ml=mu=NEQN-1, for coupled PDE/ODE problems.
nwkres= NPDE×3×NPDE+6×NXI+NPTS+15+NXI+NCODE+7×NPTS+NXFIX+1, when ​NCODE>0​ and ​NXI>0; or NPDE×3×NPDE+NPTS+21+NCODE+7×NPTS+NXFIX+2, when ​NCODE>0​ and ​NXI=0; or NPDE×3×NPDE+NPTS+21+7×NPTS+NXFIX+3, when ​NCODE=0.  
lenode= 6+intALGOPT2×NEQN+50, when the BDF method is used; or 9×NEQN+50, when the Theta method is used.  
Note: when using the sparse option, the value of LRSAVE may be too small when supplied to the integrator. An estimate of the minimum size of LRSAVE is printed on the current error message unit if ITRACE>0 and the routine returns with IFAIL=15.
34:   ISAVE(LISAVE) – INTEGER arrayCommunication Array
If IND=0, ISAVE need not be set.
If IND=1, ISAVE must be unchanged from the previous call to the routine because it contains required information about the iteration. In particular the following components of the array ISAVE concern the efficiency of the integration:
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 evaluating 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 ODE method last used in the time integration.
ISAVE5
Contains the number of Newton iterations performed by the time integrator. Each iteration involves residual evaluation of the resulting ODE system followed by a back-substitution using the LU decomposition of the Jacobian matrix.
The rest of the array is used as workspace.
35:   LISAVE – INTEGERInput
On entry: the dimension of the array ISAVE as declared in the (sub)program from which D03PRF is called. Its size depends on the type of matrix algebra selected:
  • if LAOPT='F', LISAVE25+NXFIX;
  • if LAOPT='B', LISAVENEQN+25+NXFIX;
  • if LAOPT='S', LISAVE25×NEQN+25+NXFIX.
Note: when using the sparse option, the value of LISAVE may be too small when supplied to the integrator. An estimate of the minimum size of LISAVE is printed on the current error message unit if ITRACE>0 and the routine returns with IFAIL=15.
36:   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.
ITASK=4
Normal computation of output values U at t=TOUT but without overshooting t=tcrit where tcrit is described under the parameter ALGOPT.
ITASK=5
Take one step in the time direction and return, without passing tcrit, where tcrit is described under the parameter ALGOPT.
Constraint: ITASK=1, 2, 3, 4 or 5.
37:   ITRACE – INTEGERInput
On entry: the level of trace information required from D03PRF 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
The 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.
38:   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 and the remeshing parameters NRMESH, DXMESH, TRMESH, XRATIO and CON may be reset between calls to D03PRF.
Constraint: IND=0 or 1.
On exit: IND=1.
39:   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,TOUT-TS is too small,
orITASK1, 2, 3, 4 or 5,
orat least one of the coupling points defined in array XI is outside the interval [X1,XNPTS],
orNPTS<3,
orNPDE<1,
orNLEFT not in the range 0 to NPDE,
orNORM'A' or 'M',
orLAOPT'F', 'B' or 'S',
orITOL1, 2, 3 or 4,
orIND0 or 1,
ormesh points Xi are badly ordered,
orLRSAVE is too small,
orLISAVE is too small,
orNCODE and NXI are incorrectly defined,
orIND=1 on initial entry to D03PRF,
oran element of RTOL or ATOL<0.0,
orcorresponding elements of RTOL and ATOL are both 0.0,
orNEQNNPDE×NPTS+NCODE,
orNXFIX not in the range 0 to NPTS-2,
orfixed mesh point(s) do not coincide with any of the user-supplied mesh points,
orDXMESH<0.0,
orIPMINF0, 1 or 2,
orXRATIO1.0,
orCON not in the range 0.1/NPTS-1 to 10/NPTS-1.
IFAIL=2
The underlying ODE solver cannot make any further progress, with the values of ATOL and RTOL, 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 positioning 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, BNDARY or ODEDEF, 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, BNDARY or ODEDEF. Integration was successful as far as t=TS.
IFAIL=7
The values of ATOL and RTOL are so small that the routine is unable to start the integration in time.
IFAIL=8
In either, PDEDEF, BNDARY or ODEDEF, 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 an 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 ATOL and RTOL 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 or 5.)
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). If using the sparse matrix algebra option, the values of ALGOPT29 and ALGOPT30 may be inappropriate.
IFAIL=12
In solving the ODE system, the maximum number of steps specified in ALGOPT15 have been taken.
IFAIL=13
Some error weights wi became zero during the time integration (see the description of ITOL). Pure relative error control ATOLi=0.0 was requested on a variable (the ith) which has become zero. The integration was successful as far as t=TS.
IFAIL=14
Not applicable.
IFAIL=15
When using the sparse option, the value of LISAVE or LRSAVE was insufficient (more detailed information may be directed to the current error message unit).
IFAIL=16
REMESH has been changed between calls to D03PRF.

7  Accuracy

D03PRF 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 parameters, ATOL and RTOL.

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 in Section 9). In general, a second-order problem can be solved with slightly greater accuracy using the Keller box scheme instead of a finite difference scheme (D03PPF/D03PPA 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 (D03PSF for example), or the addition of a second-order artificial dissipation term.
The time taken depends on the complexity of the system, the accuracy requested, and the frequency of the mesh updates. For a given system with fixed accuracy and mesh-update frequency it is approximately proportional to NEQN.

9  Example

This example is the first-order system
U1 t + U1 x + U2 x = 0, U2 t +4 U1 x + U2 x = 0,
for x0,1 and t0.
The initial conditions are
U1x,0 = ex, U2x,0 = x2+sin2πx2,
and the Dirichlet boundary conditions for U1 at x=0 and U2 at x=1 are given by the exact solution:
U1x,t = 12 ex+t+ex-3t+14 sin2π x-3t 2-sin2π x+t 2+2t2-2xt, U2x,t = ex-3t-ex+t+12 sin2π x-3t 2+sin2π x+t 2+x2+5t2-2xt.

9.1  Program Text

Program Text (d03prfe.f90)

9.2  Program Data

Program Data (d03prfe.d)

9.3  Program Results

Program Results (d03prfe.r)

Produced by GNUPLOT 4.4 patchlevel 0 Example Program Solution of First-order System using Moving Mesh U(1,x,t) U(1,x,t) 0 0.05 0.1 0.15 0.2 0.25 Time 0 0.2 0.4 0.6 0.8 1 x 0.8 1 1.2 1.4 1.6 1.8 2 2.2 2.4 2.6 2.8
Produced by GNUPLOT 4.4 patchlevel 0 Solution of First-order System using Moving Mesh U(2,x,t) U(2,x,t) 0 0.05 0.1 0.15 0.2 0.25 Time 0 0.2 0.4 0.6 0.8 1 x -2 -1.5 -1 -0.5 0 0.5 1 1.5

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

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