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

NAG Library Routine Document

D03PSF

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

D03PSF integrates a system of linear or nonlinear convection-diffusion equations in one space dimension, with optional source terms and scope for coupled ordinary differential equations (ODEs). The system must be posed in conservative form. This routine also includes the option of automatic adaptive spatial remeshing. 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 partial differential equations (PDEs) to a system of ODEs, and the resulting system is solved using a backward differentiation formula (BDF) method or a Theta method.

2  Specification

INTEGER  NPDE, NPTS, 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(NXI), RTOL(*), ATOL(*), ALGOPT(30), XFIX(*), DXMESH, TRMESH, XRATIO, CON, RSAVE(LRSAVE)
LOGICAL  REMESH
CHARACTER(1)  NORM, LAOPT
EXTERNAL  PDEDEF, NUMFLX, BNDARY, UVINIT, ODEDEF, MONITF

3  Description

D03PSF 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 PDE solution values
Ux,t=U1x,t,,UNPDEx,tT.  
The optional coupled ODEs are of the general form
Rit,V,V.,ξ,U*,Ux*,Ut*=0,  i=1,2,,NCODE, (3)
where the vector V is the set of ODE solution values
Vt=V1t,,VNCODEtT,  
V. denotes its derivative with respect to time, and Ux is the spatial derivative of U.
In (2), Pi,j, Fi and Ci depend on x, t, U and V; Di depends on x, t, U, Ux and V; and Si depends on x, t, U, V and linearly on V.. Note that Pi,j, Fi, Ci and Si must not depend on any space derivatives, and Pi,j, Fi, Ci and Di must not depend on any time derivatives. In terms of conservation laws, Fi, CiDi x  and Si are the convective flux, diffusion and source terms respectively.
In (3), ξ 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 PDE spatial mesh points. U*, Ux* and Ut* are the functions U, Ux and Ut evaluated at these coupling points. Each Ri may depend only linearly on time derivatives. Hence (3) may be written more precisely as
R=L-MV.-NUt*, (4)
where R=R1,,RNCODET, L is a vector of length NCODE, M is an NCODE by NCODE matrix, N is an NCODE by nξ×NPDE matrix and the entries in L, M and N 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 L, M and N. (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 user-defined mesh x1,x2,,xNPTS defined initially by you and (possibly) adapted automatically during the integration according to user-specified criteria.
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 initial remeshing, and so Ux,t0 should be specified for all values of x in the interval axb, and not just the initial mesh points.
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, 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-12=xj-1+xj/2, for j=2,3,,NPTS. The left and right values are calculated by D03PSF 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, (5)
where y=x-xj-12, i.e., y=0 corresponds to x=xj-12, 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, (6)
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. (7)
Examples are given in the documents for D03PFF and D03PLF.
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 PDEDEF. The numerical flux F^i must be supplied in NUMFLX. For problems in the form (2), the actual argument D03PLP may be used for PDEDEF. D03PLP 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.
For second-order problems, i.e., diffusion terms are present, a boundary condition is required for each PDE at both boundaries for the problem to be well-posed. If there are no diffusion 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 when using banded matrix algebra the fixed bandwidth means that only linear extrapolation is allowed, i.e., using information at just two interior points adjacent to the boundary). 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. Another method of supplying numerical boundary conditions involves the solution of the characteristic equations associated with the outgoing characteristics. Examples of both methods can be found in the documents for D03PFF and D03PLF.
The boundary conditions must be specified in BNDARY in the form
GiLx,t,U,V,V.=0  at ​x=a,  i=1,2,,NPDE, (8)
at the left-hand boundary, and
GiRx,t,U,V,V.=0  at ​x=b,  i=1,2,,NPDE, (9)
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 algebraic-differential equation system which is defined by the functions Ri must be specified in ODEDEF. You must also specify the coupling points ξ (if any) in the array XI.
In total there are NPDE×NPTS+NCODE ODEs in the time direction. This system is then integrated forwards in time using a BDF or Theta method, optionally switching between Newton's method and functional iteration (see Berzins et al. (1989) and the references therein).
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 a MONITF which specifies in an analytical or numerical form the particular aspect of the solution behaviour you wish to track. This so-called monitor function is used by the routine 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 there is 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 together 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.
The problem is subject to the following restrictions:
(i) In (1), V.jt, for j=1,2,,NCODE, may only appear linearly in the functions Si, for i=1,2,,NPDE, with a similar restriction for GiL and GiR;
(ii) Pi,j, Fi, Ci and Si must not depend on any space derivatives; and Pi,j, Ci, Di and Fi 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 fixed mesh points specified by XFIX;
(v) At least one of the functions Pi,j must be nonzero so that there is a time derivative present in the PDE problem.
For further details of the scheme, 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
Furzeland R M (1984) The construction of adaptive space meshes TNER.85.022 Thornton Research Centre, Chester
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  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.
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 and Ci may depend on x, t, U and V; Di may depend on x, t, U, Ux and V; and Si may depend on x, t, U, V and linearly on V.. PDEDEF is called approximately midway between each pair of mesh points in turn by D03PSF. The actual argument D03PLP may be used for PDEDEF for problems in the form (2). (D03PLP is included in the NAG Library.)
The specification of PDEDEF is:
SUBROUTINE PDEDEF ( NPDE, T, X, U, UX, NCODE, V, VDOT, P, C, D, S, IRES)
INTEGER  NPDE, NCODE, IRES
REAL (KIND=nag_wp)  T, X, U(NPDE), UX(NPDE), V(NCODE), VDOT(NCODE), P(NPDE,NPDE), C(NPDE), D(NPDE), S(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:     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:     NCODE – INTEGERInput
On entry: the number of coupled ODEs in the system.
7:     VNCODE – REAL (KIND=nag_wp) arrayInput
On entry: if NCODE>0, Vi contains the value of the component Vit, for i=1,2,,NCODE.
8:     VDOTNCODE – REAL (KIND=nag_wp) arrayInput
On entry: if NCODE>0, VDOTi contains the value of component V.it, for i=1,2,,NCODE.
Note: V.it, for i=1,2,,NCODE, may only appear linearly in Sj, for j=1,2,,NPDE.
9:     PNPDENPDE – REAL (KIND=nag_wp) arrayOutput
On exit: Pij must be set to the value of Pi,jx,t,U,V, for i=1,2,,NPDE and j=1,2,,NPDE.
10:   CNPDE – REAL (KIND=nag_wp) arrayOutput
On exit: Ci must be set to the value of Cix,t,U,V, for i=1,2,,NPDE.
11:   DNPDE – REAL (KIND=nag_wp) arrayOutput
On exit: Di must be set to the value of Dix,t,U,Ux,V, for i=1,2,,NPDE.
12:   SNPDE – REAL (KIND=nag_wp) arrayOutput
On exit: Si must be set to the value of Six,t,U,V,V., for i=1,2,,NPDE.
13:   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 (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 D03PSF 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 D03PSF is called. Parameters denoted as Input must not be changed by this procedure.
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 D03PSF.
The specification of NUMFLX is:
SUBROUTINE NUMFLX ( NPDE, T, X, NCODE, V, ULEFT, URIGHT, FLUX, IRES)
INTEGER  NPDE, NCODE, IRES
REAL (KIND=nag_wp)  T, X, V(NCODE), ULEFT(NPDE), URIGHT(NPDE), FLUX(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:     NCODE – INTEGERInput
On entry: the number of coupled ODEs in the system.
5:     VNCODE – REAL (KIND=nag_wp) arrayInput
On entry: if NCODE>0, Vi contains the value of the component Vit, for i=1,2,,NCODE.
6:     ULEFTNPDE – REAL (KIND=nag_wp) arrayInput
On entry: ULEFTi contains the left value of the component Uix, for i=1,2,,NPDE.
7:     URIGHTNPDE – REAL (KIND=nag_wp) arrayInput
On entry: URIGHTi contains the right value of the component Uix, for i=1,2,,NPDE.
8:     FLUXNPDE – REAL (KIND=nag_wp) arrayOutput
On exit: FLUXi must be set to the numerical flux F^i, for i=1,2,,NPDE.
9:     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 (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 D03PSF 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 D03PSF is called. Parameters denoted as Input must not be changed by this procedure.
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 (8) and (9).
The specification of BNDARY is:
SUBROUTINE BNDARY ( NPDE, NPTS, T, X, U, NCODE, V, VDOT, IBND, G, IRES)
INTEGER  NPDE, NPTS, NCODE, IBND, IRES
REAL (KIND=nag_wp)  T, X(NPTS), U(NPDE,NPTS), V(NCODE), VDOT(NCODE), G(NPDE)
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:     UNPDENPTS – REAL (KIND=nag_wp) arrayInput
On entry: Uij contains the value of the component Uix,t at x=Xj, for i=1,2,,NPDE and j=1,2,,NPTS.
Note:  if banded matrix algebra is to be used then the functions GiL and GiR may depend on the value of Uix,t at the boundary point and the two adjacent points only.
6:     NCODE – INTEGERInput
On entry: the number of coupled ODEs in the system.
7:     VNCODE – REAL (KIND=nag_wp) arrayInput
On entry: if NCODE>0, Vi contains the value of the component Vit, for i=1,2,,NCODE.
8:     VDOTNCODE – REAL (KIND=nag_wp) arrayInput
On entry: if NCODE>0, VDOTi contains the value of component V.it, for i=1,2,,NCODE.
Note: V.it, for i=1,2,,NCODE, may only appear linearly in GjL and GjR, for j=1,2,,NPDE.
9:     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.
10:   GNPDE – REAL (KIND=nag_wp) arrayOutput
On exit: Gi must contain the ith component of either GiL or GiR in (8) and (9), depending on the value of IBND, for i=1,2,,NPDE.
11:   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 (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 D03PSF 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 D03PSF is called. Parameters denoted as Input must not be changed by this procedure.
7:     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 axb.
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:     XNPTS – REAL (KIND=nag_wp) arrayInput
On entry: the current mesh. Xi contains the value of xi, for i=1,2,,NPTS.
5:     XINXI – REAL (KIND=nag_wp) arrayInput
On entry: if NXI>0, XIi contains the ODE/PDE coupling point, ξi, for i=1,2,,NXI.
6:     UNPDENPTS – 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:     VNCODE – 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 D03PSF is called. Parameters denoted as Input must not be changed by this procedure.
8:     UNEQN – 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, and UNPTS×NPDE+k contains Vkt, for k=1,2,,NCODE, all evaluated at t=TS.
9:     NPTS – INTEGERInput
On entry: the number of mesh points in the interval a,b.
Constraint: NPTS3.
10:   XNPTS – REAL (KIND=nag_wp) arrayInput/Output
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.
On exit: the final values of the mesh points.
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:     VNCODE – REAL (KIND=nag_wp) arrayInput
On entry: if NCODE>0, Vi contains the value of the component Vit, for i=1,2,,NCODE.
5:     VDOTNCODE – 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:     XINXI – REAL (KIND=nag_wp) arrayInput
On entry: if NXI>0, XIi contains the ODE/PDE coupling point, ξi, for i=1,2,,NXI.
8:     UCPNPDE* – 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:     UCPXNPDE* – 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:   UCPTNPDE* – 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:   RNCODE – REAL (KIND=nag_wp) arrayOutput
On exit: Ri must contain the ith component of R, for i=1,2,,NCODE, where R is defined as
R=L-MV.-NUt*, (10)
or
R=-MV.-NUt*. (11)
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 (10) must be used.
IRES=-1
Equation (11) 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 D03PSF 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 D03PSF 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:   XINXI – REAL (KIND=nag_wp) arrayInput
On entry: if NXI>0, XIi, for i=1,2,,NXI, must be set to the ODE/PDE coupling points.
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
On entry: a value to indicate the form of the local error test. If ei is the estimated local error for Ui, for i=1,2,,NEQN, and , denotes the norm, then the error test to be satisfied is ei<1.0. ITOL indicates to D03PSF whether to interpret either or both of RTOL and ATOL as a vector or scalar in the formation of the weights wi used in the calculation of the norm (see the description of NORM):
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
Constraint: ITOL=1, 2, 3 or 4.
19:   NORM – CHARACTER(1)Input
On entry: the type of norm to be used.
NORM='1'
Averaged L1 norm.
NORM='2'
Averaged L2 norm.
If Unorm denotes the norm of the vector U of length NEQN, then for the averaged L1 norm
Unorm=1NEQNi=1NEQNUi/wi,  
and for the averaged L2 norm
Unorm=1NEQNi=1NEQNUi/wi2,  
See the description of ITOL for the formulation of the weight vector w.
Constraint: NORM='1' or '2'.
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 (NCODE=0). Also, the banded option should not be used if the boundary conditions involve solution components at points other than the boundary and the immediately adjacent two points.
21:   ALGOPT30 – 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 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 the 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 the 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 matrix 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 D03PSF. 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.
Constraints:
  • XFIXi<XFIXi+1, for i=1,2,,NXFIX-1;
  • 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 XNPTS 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: specifies the spatial remeshing frequency and criteria for the calculation and adoption of a new mesh.
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 D03PSF 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 D03PSF 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: an input bound on the 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
xixi+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.0/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:     XNPTS – REAL (KIND=nag_wp) arrayInput
On entry: the current mesh. Xi contains the value of xi, for i=1,2,,NPTS.
5:     UNPDENPTS – 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:     FMONNPTS – 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 D03PSF is called. Parameters denoted as Input must not be changed by this procedure.
32:   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.
33:   LRSAVE – INTEGERInput
On entry: the dimension of the array RSAVE as declared in the (sub)program from which D03PSF is called. Its size depends on the type of matrix algebra selected.
If LAOPT='F', LRSAVENEQN×NEQN+NEQN+nwkres+lenode.
If LAOPT='B', LRSAVE3×mlu+1×NEQN+nwkres+lenode.
If LAOPT='S', LRSAVE4×NEQN+11×NEQN/2+1+nwkres+lenode.
Where
mlu is the lower or upper half bandwidths such that
mlu=3×NPDE-1, for PDE problems only (no coupled ODEs); or
mlu=NEQN-1, for coupled PDE/ODE problems.
nwkres= NPDE×2×NPTS+6×NXI+3×NPDE+26+NXI+NCODE+7×NPTS+NXFIX+1, when ​NCODE>0​ and ​NXI>0; or NPDE×2×NPTS+3×NPDE+32+NCODE+7×NPTS+NXFIX+2, when ​NCODE>0​ and ​NXI=0; or ​ NPDE×2×NPTS+3×NPDE+32+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 LAOPT='S', 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:   ISAVELISAVE – 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 BDF method last used in the time integration, if applicable. When the Theta method is used, ISAVE4 contains no useful information.
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.
35:   LISAVE – INTEGERInput
On entry: the dimension of the array ISAVE as declared in the (sub)program from which D03PSF is called. Its size depends on the type of matrix algebra selected:
  • if LAOPT='F', LISAVE25;
  • if LAOPT='B', LISAVENEQN+NXFIX+25;
  • if LAOPT='S', LISAVE25×NEQN+NXFIX+25.
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 D03PSF 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, then -1 is assumed and similarly if ITRACE>3, then 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.
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, IFAIL, NRMESH and TRMESH may be reset between calls to D03PSF.
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,TSTOUT,
orTOUT-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],
orthe coupling points are not in strictly increasing order,
orNPTS<3,
orNPDE<1,
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 D03PSF,
orNEQNNPDE×NPTS+NCODE,
oran element of RTOL or ATOL<0.0,
orcorresponding elements of RTOL and ATOL are both 0.0,
orNORM'1' or '2',
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.0/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 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, BNDARY or ODEDEF, 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, 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 one of PDEDEF, NUMFLX, 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 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 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 when ITRACE1). 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
One or more of the functions Pi,j, Di or Ci was detected as depending on time derivatives, which is not permissible.
IFAIL=15
When using the sparse option, the value of LISAVE or LRSAVE was not sufficient (more detailed information may be directed to the current error message unit, see X04AAF).
IFAIL=16
REMESH has been changed between calls to D03PSF.
IFAIL=17
FMON is negative at one or more mesh points, or zero mesh spacing has been obtained due to an inappropriate choice of monitor function.
IFAIL=-99
An unexpected error has been triggered by this routine. Please contact NAG.
See Section 3.8 in the Essential Introduction for further information.
IFAIL=-399
Your licence key may have expired or may not have been installed correctly.
See Section 3.7 in the Essential Introduction for further information.
IFAIL=-999
Dynamic memory allocation failed.
See Section 3.6 in the Essential Introduction for further information.

7  Accuracy

D03PSF 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  Parallelism and Performance

D03PSF is threaded by NAG for parallel execution in multithreaded implementations of the NAG Library.
D03PSF 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

D03PSF 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 scheme routines 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 ALGOPT13. It is worth experimenting with this parameter, 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, 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.

10  Example

For this routine two examples are presented, with a main program and two example problems given in Example 1 (EX1) and Example 2 (EX2).
Example 1 (EX1)
This example is a simple model of the advection and diffusion of a cloud of material:
U t + W U x = C 2 U x2 ,  
for x0,1 and t00.3. In this example the constant wind speed W=1 and the diffusion coefficient C=0.002.
The cloud does not reach the boundaries during the time of integration, and so the two (physical) boundary conditions are simply U0,t=U1,t=0.0, and the initial condition is
U x,0 = sin πx-a b-a ,  axb,  
and Ux,0=0 elsewhere, where a=0.2 and b=0.4.
The numerical flux is simply F^=WUL.
The monitor function for remeshing is taken to be the absolute value of the second derivative of U.
Example 2 (EX2)
This example is a linear advection equation with a nonlinear source term and discontinuous initial profile:
u t + u x =-puu-1u-12,  
for 0x1 and t0. The discontinuity is modelled by a ramp function of width 0.01 and gradient 100, so that the exact solution at any time t0 is
ux,t=1.0+maxminδ,0,-1,  
where δ=1000.1-x+t. The initial profile is given by the exact solution. The characteristic points into the domain at x=0 and out of the domain at x=1, and so a physical boundary condition u0,t=1 is imposed at x=0, with a numerical boundary condition at x=1 which can be specified as u1,t=0 since the discontinuity does not reach x=1 during the time of integration.
The numerical flux is simply F^=UL at all times.
The remeshing monitor function (described below) is chosen to create an increasingly fine mesh towards the discontinuity in order to ensure good resolution of the discontinuity, but without loss of efficiency in the surrounding regions. However, refinement must be limited so that the time step required for stability does not become unrealistically small. The region of refinement must also keep up with the discontinuity as it moves across the domain, and hence it cannot be so small that the discontinuity moves out of the refined region between remeshing.
The above requirements mean that the use of the first or second spatial derivative of U for the monitor function is inappropriate; the large relative size of either derivative in the region of the discontinuity leads to extremely small mesh-spacing in a very limited region, and the solution is then far more expensive than for a very fine fixed mesh.
An alternative monitor function based on a cosine function proves very successful. It is only semi-automatic as it requires some knowledge of the solution (for problems without an exact solution an initial approximate solution can be obtained using a coarse fixed mesh). On each call to MONITF the discontinuity is located by finding the maximum spatial derivative of the solution. On the first call the desired width of the region of nonzero monitor function is set (this can be changed at a later time if desired). Then on each call the monitor function is assigned using a cosine function so that it has a value of one at the discontinuity down to zero at the edges of the predetermined region of refinement, and zero outside the region. Thus the monitor function and the subsequent refinement are limited, and the region is large enough to ensure that there is always sufficient refinement at the discontinuity.

10.1  Program Text

Program Text (d03psfe.f90)

10.2  Program Data

Program Data (d03psfe.d)

10.3  Program Results

Program Results (d03psfe.r)

GnuplotProduced by GNUPLOT 4.6 patchlevel 3 Example Program 1 Advection and Diffusion of a Cloud of Material u(x,t) gnuplot_plot_1 gnuplot_plot_2 0 0.05 0.1 0.15 0.2 0.25 0.3 Time 0 0.2 0.4 0.6 0.8 1 x 0 0.1 0.2 0.3 0.4 0.5 0.6 0.7 0.8 0.9 1
GnuplotProduced by GNUPLOT 4.6 patchlevel 3 Example Program 2 Linear Advection Equation with Non-linear Source Term u(x,t) gnuplot_plot_1 gnuplot_plot_2 0 0.05 0.1 0.15 0.2 0.25 0.3 0.35 0.4 Time 0 0.2 0.4 0.6 0.8 1 x −0.2 0 0.2 0.4 0.6 0.8 1 1.2

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

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