NAG FL Interface
d03pdf  (dim1_parab_coll_old)
d03pda (dim1_parab_coll)

Settings help

FL Name Style:


FL Specification Language:


1 Purpose

d03pdf/​d03pda integrates a system of linear or nonlinear parabolic partial differential equations (PDEs) in one space variable. The spatial discretization is performed using a Chebyshev C0 collocation method, and the method of lines is employed to reduce the PDEs to a system of ordinary differential equations (ODEs). The resulting system is solved using a backward differentiation formula method.
d03pda is a version of d03pdf that has additional arguments in order to make it safe for use in multithreaded applications (see Section 5).

2 Specification

2.1 Specification for d03pdf

Fortran Interface
Subroutine d03pdf ( npde, m, ts, tout, pdedef, bndary, u, nbkpts, xbkpts, npoly, npts, x, uinit, acc, rsave, lrsave, isave, lisave, itask, itrace, ind, ifail)
Integer, Intent (In) :: npde, m, nbkpts, npoly, npts, lrsave, lisave, itask, itrace
Integer, Intent (Inout) :: isave(lisave), ind, ifail
Real (Kind=nag_wp), Intent (In) :: tout, xbkpts(nbkpts), acc
Real (Kind=nag_wp), Intent (Inout) :: ts, u(npde,npts), rsave(lrsave)
Real (Kind=nag_wp), Intent (Out) :: x(npts)
External :: pdedef, bndary, uinit
C Header Interface
#include <nag.h>
void  d03pdf_ (const Integer *npde, const Integer *m, double *ts, const double *tout,
void (NAG_CALL *pdedef)(const Integer *npde, const double *t, const double x[], const Integer *nptl, const double u[], const double ux[], double p[], double q[], double r[], Integer *ires),
void (NAG_CALL *bndary)(const Integer *npde, const double *t, const double u[], const double ux[], const Integer *ibnd, double beta[], double gamma[], Integer *ires),
double u[], const Integer *nbkpts, const double xbkpts[], const Integer *npoly, const Integer *npts, double x[],
void (NAG_CALL *uinit)(const Integer *npde, const Integer *npts, const double x[], double u[]),
const double *acc, double rsave[], const Integer *lrsave, Integer isave[], const Integer *lisave, const Integer *itask, const Integer *itrace, Integer *ind, Integer *ifail)

2.2 Specification for d03pda

Fortran Interface
Integer, Intent (In) :: npde, m, nbkpts, npoly, npts, lrsave, lisave, itask, itrace
Integer, Intent (Inout) :: isave(lisave), ind, iuser(*), iwsav(505), ifail
Real (Kind=nag_wp), Intent (In) :: tout, xbkpts(nbkpts), acc
Real (Kind=nag_wp), Intent (Inout) :: ts, u(npde,npts), rsave(lrsave), ruser(*), rwsav(1100)
Real (Kind=nag_wp), Intent (Out) :: x(npts)
Logical, Intent (Inout) :: lwsav(100)
Character (80), Intent (InOut) :: cwsav(10)
External :: pdedef, bndary, uinit
C Header Interface
#include <nag.h>
void  d03pda_ (const Integer *npde, const Integer *m, double *ts, const double *tout,
void (NAG_CALL *pdedef)(const Integer *npde, const double *t, const double x[], const Integer *nptl, const double u[], const double ux[], double p[], double q[], double r[], Integer *ires, Integer iuser[], double ruser[]),
void (NAG_CALL *bndary)(const Integer *npde, const double *t, const double u[], const double ux[], const Integer *ibnd, double beta[], double gamma[], Integer *ires, Integer iuser[], double ruser[]),
double u[], const Integer *nbkpts, const double xbkpts[], const Integer *npoly, const Integer *npts, double x[],
void (NAG_CALL *uinit)(const Integer *npde, const Integer *npts, const double x[], double u[], Integer iuser[], double ruser[]),
const double *acc, double rsave[], const Integer *lrsave, Integer isave[], const Integer *lisave, const Integer *itask, const Integer *itrace, Integer *ind, Integer iuser[], double ruser[], char cwsav[], logical lwsav[], Integer iwsav[], double rwsav[], Integer *ifail, const Charlen length_cwsav)

3 Description

d03pdf/​d03pda integrates the system of parabolic equations:
j=1npdePi,j Uj t +Qi=x-m x (xmRi),  i=1,2,,npde,  axb,tt0, (1)
where Pi,j, Qi and Ri depend on x, t, U, Ux and the vector U is the set of solution values
U (x,t) = [ U 1 (x,t) ,, U npde (x,t) ] T , (2)
and the vector Ux is its partial derivative with respect to x. Note that Pi,j, Qi and Ri must not depend on U t .
The integration in time is from t0 to tout, over the space interval axb, where a=x1 and b=xnbkpts are the leftmost and rightmost of a user-defined set of break-points x1,x2,,xnbkpts. The coordinate system in space is defined by the value of m; m=0 for Cartesian coordinates, m=1 for cylindrical polar coordinates and m=2 for spherical polar coordinates.
The system is defined by the functions Pi,j, Qi and Ri which must be specified in pdedef.
The initial values of the functions U(x,t) must be given at t=t0, and must be specified in uinit.
The functions Ri, for i=1,2,,npde, which may be thought of as fluxes, are also used in the definition of the boundary conditions for each equation. The boundary conditions must have the form
βi(x,t)Ri(x,t,U,Ux)=γi(x,t,U,Ux),  i=1,2,,npde, (3)
where x=a or x=b.
The boundary conditions must be specified in bndary. Thus, the problem is subject to the following restrictions:
  1. (i)t0<tout, so that integration is in the forward direction;
  2. (ii)Pi,j, Qi and the flux Ri must not depend on any time derivatives;
  3. (iii)the evaluation of the functions Pi,j, Qi and Ri is done at both the break-points and internally selected points for each element in turn, that is Pi,j, Qi and Ri are evaluated twice at each break-point. Any discontinuities in these functions must, therefore, be at one or more of the break-points x1,x2,,xnbkpts;
  4. (iv)at least one of the functions Pi,j must be nonzero so that there is a time derivative present in the problem;
  5. (v)if m>0 and x1=0.0, which is the left boundary point, then it must be ensured that the PDE solution is bounded at this point. This can be done by either specifying the solution at x=0.0 or by specifying a zero flux there, that is βi=1.0 and γi=0.0. See also Section 9.
The parabolic equations are approximated by a system of ODEs in time for the values of Ui at the mesh points. This ODE system is obtained by approximating the PDE solution between each pair of break-points by a Chebyshev polynomial of degree npoly. The interval between each pair of break-points is treated by d03pdf/​d03pda as an element, and on this element, a polynomial and its space and time derivatives are made to satisfy the system of PDEs at npoly-1 spatial points, which are chosen internally by the code and the break-points. In the case of just one element, the break-points are the boundaries. The user-defined break-points and the internally selected points together define the mesh. The smallest value that npoly can take is one, in which case, the solution is approximated by piecewise linear polynomials between consecutive break-points and the method is similar to an ordinary finite element method.
In total there are (nbkpts-1)×npoly+1 mesh points in the spatial direction, and npde×((nbkpts-1)×npoly+1) ODEs in the time direction; one ODE at each break-point for each PDE component and (npoly-1) ODEs for each PDE component between each pair of break-points. The system is then integrated forwards in time using a backward differentiation formula method.

4 References

Berzins M (1990) Developments in the NAG Library software for parabolic equations Scientific Software Systems (eds J C Mason and M G Cox) 59–72 Chapman and Hall
Berzins M and Dew P M (1991) Algorithm 690: Chebyshev polynomial software for elliptic-parabolic systems of PDEs ACM Trans. Math. Software 17 178–206
Zaturska N B, Drazin P G and Banks W H H (1988) On the flow of a viscous fluid driven along a channel by a suction at porous walls Fluid Dynamics Research 4

5 Arguments

1: npde Integer Input
On entry: the number of PDEs in the system to be solved.
Constraint: npde1.
2: m Integer Input
On entry: the coordinate system used:
m=0
Indicates Cartesian coordinates.
m=1
Indicates cylindrical polar coordinates.
m=2
Indicates spherical polar coordinates.
Constraint: m=0, 1 or 2.
3: 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.
4: tout Real (Kind=nag_wp) Input
On entry: the final value of t to which the integration is to be carried out.
5: pdedef Subroutine, supplied by the user. External Procedure
pdedef must compute the values of the functions Pi,j, Qi and Ri which define the system of PDEs. The functions may depend on x, t, U and Ux and must be evaluated at a set of points.
The specification of pdedef for d03pdf is:
Fortran Interface
Subroutine pdedef ( npde, t, x, nptl, u, ux, p, q, r, ires)
Integer, Intent (In) :: npde, nptl
Integer, Intent (Inout) :: ires
Real (Kind=nag_wp), Intent (In) :: t, x(nptl), u(npde,nptl), ux(npde,nptl)
Real (Kind=nag_wp), Intent (Out) :: p(npde,npde,nptl), q(npde,nptl), r(npde,nptl)
C Header Interface
void  pdedef (const Integer *npde, const double *t, const double x[], const Integer *nptl, const double u[], const double ux[], double p[], double q[], double r[], Integer *ires)
The specification of pdedef for d03pda is:
Fortran Interface
Subroutine pdedef ( npde, t, x, nptl, u, ux, p, q, r, ires, iuser, ruser)
Integer, Intent (In) :: npde, nptl
Integer, Intent (Inout) :: ires, iuser(*)
Real (Kind=nag_wp), Intent (In) :: t, x(nptl), u(npde,nptl), ux(npde,nptl)
Real (Kind=nag_wp), Intent (Inout) :: ruser(*)
Real (Kind=nag_wp), Intent (Out) :: p(npde,npde,nptl), q(npde,nptl), r(npde,nptl)
C Header Interface
void  pdedef (const Integer *npde, const double *t, const double x[], const Integer *nptl, const double u[], const double ux[], double p[], double q[], double r[], Integer *ires, Integer iuser[], double ruser[])
1: npde Integer Input
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(nptl) Real (Kind=nag_wp) array Input
On entry: contains a set of mesh points at which Pi,j, Qi and Ri are to be evaluated. x(1) and x(nptl) contain successive user-supplied break-points and the elements of the array will satisfy x(1)<x(2)<<x(nptl).
4: nptl Integer Input
On entry: the number of points at which evaluations are required (the value of npoly+1).
5: u(npde,nptl) Real (Kind=nag_wp) array Input
On entry: u(i,j) contains the value of the component Ui(x,t) where x=x(j), for i=1,2,,npde and j=1,2,,nptl.
6: ux(npde,nptl) Real (Kind=nag_wp) array Input
On entry: ux(i,j) contains the value of the component Ui(x,t) x where x=x(j), for i=1,2,,npde and j=1,2,,nptl.
7: p(npde,npde,nptl) Real (Kind=nag_wp) array Output
On exit: p(i,j,k) must be set to the value of Pi,j(x,t,U,Ux) where x=x(k), for i=1,2,,npde, j=1,2,,npde and k=1,2,,nptl.
8: q(npde,nptl) Real (Kind=nag_wp) array Output
On exit: q(i,j) must be set to the value of Qi(x,t,U,Ux) where x=x(j), for i=1,2,,npde and j=1,2,,nptl.
9: r(npde,nptl) Real (Kind=nag_wp) array Output
On exit: r(i,j) must be set to the value of Ri(x,t,U,Ux) where x=x(j), for i=1,2,,npde and j=1,2,,nptl.
10: ires Integer Input/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, d03pdf/​d03pda returns to the calling subroutine with the error indicator set to ifail=4.
Note: the following are additional arguments for specific use with d03pda. Users of d03pdf therefore need not read the remainder of this description.
11: iuser(*) Integer array User Workspace
12: ruser(*) Real (Kind=nag_wp) array User Workspace
pdedef is called with the arguments iuser and ruser as supplied to d03pdf/​d03pda. You should use the arrays iuser and ruser to supply information to pdedef.
pdedef must either be a module subprogram USEd by, or declared as EXTERNAL in, the (sub)program from which d03pdf/​d03pda is called. Arguments denoted as Input must not be changed by this procedure.
Note: pdedef should not return floating-point NaN (Not a Number) or infinity values, since these are not handled by d03pdf/​d03pda. If your code inadvertently does return any NaNs or infinities, d03pdf/​d03pda is likely to produce unexpected results.
6: bndary Subroutine, supplied by the user. External Procedure
bndary must compute the functions βi and γi which define the boundary conditions as in equation (3).
The specification of bndary for d03pdf is:
Fortran Interface
Subroutine bndary ( npde, t, u, ux, ibnd, beta, gamma, ires)
Integer, Intent (In) :: npde, ibnd
Integer, Intent (Inout) :: ires
Real (Kind=nag_wp), Intent (In) :: t, u(npde), ux(npde)
Real (Kind=nag_wp), Intent (Out) :: beta(npde), gamma(npde)
C Header Interface
void  bndary (const Integer *npde, const double *t, const double u[], const double ux[], const Integer *ibnd, double beta[], double gamma[], Integer *ires)
The specification of bndary for d03pda is:
Fortran Interface
Subroutine bndary ( npde, t, u, ux, ibnd, beta, gamma, ires, iuser, ruser)
Integer, Intent (In) :: npde, ibnd
Integer, Intent (Inout) :: ires, iuser(*)
Real (Kind=nag_wp), Intent (In) :: t, u(npde), ux(npde)
Real (Kind=nag_wp), Intent (Inout) :: ruser(*)
Real (Kind=nag_wp), Intent (Out) :: beta(npde), gamma(npde)
C Header Interface
void  bndary (const Integer *npde, const double *t, const double u[], const double ux[], const Integer *ibnd, double beta[], double gamma[], Integer *ires, Integer iuser[], double ruser[])
1: npde Integer Input
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: u(npde) Real (Kind=nag_wp) array Input
On entry: u(i) contains the value of the component Ui(x,t) at the boundary specified by ibnd, for i=1,2,,npde.
4: ux(npde) Real (Kind=nag_wp) array Input
On entry: ux(i) contains the value of the component Ui(x,t) x at the boundary specified by ibnd, for i=1,2,,npde.
5: ibnd Integer Input
On entry: specifies which boundary conditions are to be evaluated.
ibnd=0
bndary must set up the coefficients of the left-hand boundary, x=a.
ibnd0
bndary must set up the coefficients of the right-hand boundary, x=b.
6: beta(npde) Real (Kind=nag_wp) array Output
On exit: beta(i) must be set to the value of βi(x,t) at the boundary specified by ibnd, for i=1,2,,npde.
7: gamma(npde) Real (Kind=nag_wp) array Output
On exit: gamma(i) must be set to the value of γi(x,t,U,Ux) at the boundary specified by ibnd, for i=1,2,,npde.
8: ires Integer Input/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, d03pdf/​d03pda returns to the calling subroutine with the error indicator set to ifail=4.
Note: the following are additional arguments for specific use with d03pda. Users of d03pdf therefore need not read the remainder of this description.
9: iuser(*) Integer array User Workspace
10: ruser(*) Real (Kind=nag_wp) array User Workspace
bndary is called with the arguments iuser and ruser as supplied to d03pdf/​d03pda. You should use the arrays iuser and ruser to supply information to bndary.
bndary must either be a module subprogram USEd by, or declared as EXTERNAL in, the (sub)program from which d03pdf/​d03pda is called. Arguments denoted as Input must not be changed by this procedure.
Note: bndary should not return floating-point NaN (Not a Number) or infinity values, since these are not handled by d03pdf/​d03pda. If your code inadvertently does return any NaNs or infinities, d03pdf/​d03pda is likely to produce unexpected results.
7: u(npde,npts) Real (Kind=nag_wp) array Input/Output
On entry: if ind=1 the value of u must be unchanged from the previous call.
On exit: u(i,j) will contain the computed solution at t=ts.
8: nbkpts Integer Input
On entry: the number of break-points in the interval [a,b].
Constraint: nbkpts2.
9: xbkpts(nbkpts) Real (Kind=nag_wp) array Input
On entry: the values of the break-points in the space direction. xbkpts(1) must specify the left-hand boundary, a, and xbkpts(nbkpts) must specify the right-hand boundary, b.
Constraint: xbkpts(1)<xbkpts(2)<<xbkpts(nbkpts).
10: npoly Integer Input
On entry: the degree of the Chebyshev polynomial to be used in approximating the PDE solution between each pair of break-points.
Constraint: 1npoly49.
11: npts Integer Input
On entry: the number of mesh points in the interval [a,b].
Constraint: npts=(nbkpts-1)×npoly+1.
12: x(npts) Real (Kind=nag_wp) array Output
On exit: the mesh points chosen by d03pdf/​d03pda in the spatial direction. The values of x will satisfy x(1)<x(2)<<x(npts).
13: uinit Subroutine, supplied by the user. External Procedure
uinit must compute the initial values of the PDE components Ui(xj,t0), for i=1,2,,npde and j=1,2,,npts.
The specification of uinit for d03pdf is:
Fortran Interface
Subroutine uinit ( npde, npts, x, u)
Integer, Intent (In) :: npde, npts
Real (Kind=nag_wp), Intent (In) :: x(npts)
Real (Kind=nag_wp), Intent (Out) :: u(npde,npts)
C Header Interface
void  uinit (const Integer *npde, const Integer *npts, const double x[], double u[])
The specification of uinit for d03pda is:
Fortran Interface
Subroutine uinit ( npde, npts, x, u, iuser, ruser)
Integer, Intent (In) :: npde, npts
Integer, Intent (Inout) :: iuser(*)
Real (Kind=nag_wp), Intent (In) :: x(npts)
Real (Kind=nag_wp), Intent (Inout) :: ruser(*)
Real (Kind=nag_wp), Intent (Out) :: u(npde,npts)
C Header Interface
void  uinit (const Integer *npde, const Integer *npts, const double x[], double u[], Integer iuser[], double ruser[])
1: npde Integer Input
On entry: the number of PDEs in the system.
2: npts Integer Input
On entry: the number of mesh points in the interval [a,b].
3: x(npts) Real (Kind=nag_wp) array Input
On entry: x(j), contains the values of the jth mesh point, for j=1,2,,npts.
4: u(npde,npts) Real (Kind=nag_wp) array Output
On exit: u(i,j) must be set to the initial value Ui(xj,t0), for i=1,2,,npde and j=1,2,,npts.
Note: the following are additional arguments for specific use with d03pda. Users of d03pdf therefore need not read the remainder of this description.
5: iuser(*) Integer array User Workspace
6: ruser(*) Real (Kind=nag_wp) array User Workspace
uinit is called with the arguments iuser and ruser as supplied to d03pdf/​d03pda. You should use the arrays iuser and ruser to supply information to uinit.
uinit must either be a module subprogram USEd by, or declared as EXTERNAL in, the (sub)program from which d03pdf/​d03pda is called. Arguments denoted as Input must not be changed by this procedure.
Note: uinit should not return floating-point NaN (Not a Number) or infinity values, since these are not handled by d03pdf/​d03pda. If your code inadvertently does return any NaNs or infinities, d03pdf/​d03pda is likely to produce unexpected results.
14: acc Real (Kind=nag_wp) Input
On entry: a positive quantity for controlling the local error estimate in the time integration. If E(i,j) is the estimated error for Ui at the jth mesh point, the error test is:
|E(i,j)|=acc×(1.0+|u(i,j)|).  
Constraint: acc>0.0.
15: rsave(lrsave) Real (Kind=nag_wp) array Communication 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.
16: lrsave Integer Input
On entry: the dimension of the array rsave as declared in the (sub)program from which d03pdf/​d03pda is called.
Constraint: lrsave11×npde×npts+50+nwkres+lenode.
17: isave(lisave) Integer array Communication Array
If ind=0, isave need not be set on entry.
If ind=1, isave must be unchanged from the previous call to the routine because it contains required information about the iteration. In particular:
isave(1)
Contains the number of steps taken in time.
isave(2)
Contains the number of residual evaluations of the resulting ODE system used. One such evaluation involves computing the PDE functions at all the mesh points, as well as one evaluation of the functions in the boundary conditions.
isave(3)
Contains the number of Jacobian evaluations performed by the time integrator.
isave(4)
Contains the order of the last backward differentiation formula method used.
isave(5)
Contains the number of Newton iterations performed by the time integrator. Each iteration involves an ODE residual evaluation followed by a back-substitution using the LU decomposition of the Jacobian matrix.
18: lisave Integer Input
On entry: the dimension of the array isave as declared in the (sub)program from which d03pdf/​d03pda is called.
Constraint: lisavenpde×npts+24.
19: itask Integer Input
On entry: specifies the task to be performed by the ODE integrator.
itask=1
Normal computation of output values u at t=tout.
itask=2
One step and return.
itask=3
Stop at first internal integration point at or beyond t=tout.
Constraint: itask=1, 2 or 3.
20: itrace Integer Input
On entry: the level of trace information required from d03pdf/​d03pda and the underlying ODE solver. itrace may take the value −1, 0, 1, 2 or 3.
itrace=−1
No output is generated.
itrace=0
Only warning messages from the PDE solver are printed on the current error message unit (see x04aaf).
itrace>0
Output from the underlying ODE solver is printed on the current advisory message unit (see x04abf). This output contains details of Jacobian entries, the nonlinear iteration and the time integration during the computation of the ODE system.
If itrace<−1, −1 is assumed and similarly if itrace>3, 3 is assumed.
The advisory messages are given in greater detail as itrace increases. You are advised to set itrace=0, unless you are experienced with Sub-chapter D02M–N.
21: ind Integer Input/Output
On entry: indicates whether this is a continuation call or a new integration.
ind=0
Starts or restarts the integration in time.
ind=1
Continues the integration after an earlier exit from the routine. In this case, only the arguments tout and ifail should be reset between calls to d03pdf/​d03pda.
Constraint: ind=0 or 1.
On exit: ind=1.
22: ifail Integer Input/Output
Note: for d03pda, ifail does not occur in this position in the argument list. See the additional arguments described below.
On entry: ifail must be set to 0, −1 or 1 to set behaviour on detection of an error; these values have no effect when no error is detected.
A value of 0 causes the printing of an error message and program execution will be halted; otherwise program execution continues. A value of −1 means that an error message is printed while a value of 1 means that it is not.
If halting is not appropriate, the value −1 or 1 is recommended. If message printing is undesirable, then the value 1 is recommended. Otherwise, the value 0 is recommended. 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).
Note: the following are additional arguments for specific use with d03pda. Users of d03pdf therefore need not read the remainder of this description.
22: iuser(*) Integer array User Workspace
23: ruser(*) Real (Kind=nag_wp) array User Workspace
iuser and ruser are not used by d03pdf/​d03pda, but are passed directly to pdedef, bndary and uinit and may be used to pass information to these routines.
24: cwsav(10) Character(80) array Communication Array
25: lwsav(100) Logical array Communication Array
26: iwsav(505) Integer array Communication Array
27: rwsav(1100) Real (Kind=nag_wp) array Communication Array
If ind=0, cwsav, lwsav, iwsav and rwsav need not be set on entry.
If ind=1, cwsav, lwsav, iwsav and rwsav must be unchanged from the previous call to d03pdf/​d03pda.
28: ifail Integer Input/Output
Note: see the argument description for ifail above.

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, acc=value.
Constraint: acc>0.0.
On entry, i=value, xbkpts(i)=value, j=value and xbkpts(j)=value.
Constraint: xbkpts(1)<xbkpts(2)<<xbkpts(nbkpts).
On entry, ind=value.
Constraint: ind=0 or 1.
On entry, itask=value.
Constraint: itask=1, 2 or 3.
On entry, lisave=value.
Constraint: lisavevalue.
On entry, lrsave=value.
Constraint: lrsavevalue.
On entry, m=value.
Constraint: m=0, 1 or 2.
On entry, m=value and xbkpts(1)=value.
Constraint: m0 or xbkpts(1)0.0
On entry, nbkpts=value.
Constraint: nbkpts2.
On entry, npde=value.
Constraint: npde1.
On entry, npoly=value.
Constraint: npoly49.
On entry, npoly=value.
Constraint: npoly1.
On entry, npts=value, nbkpts=value and npoly=value.
Constraint: npts=(nbkpts-1)×npoly+1.
On entry, on initial entry ind=1.
Constraint: on initial entry ind=0.
On entry, tout=value and ts=value.
Constraint: tout>ts.
On entry, tout-ts is too small: tout=value and ts=value.
ifail=2
Underlying ODE solver cannot make further progress from the point ts with the supplied value of acc. ts=value, acc=value.
ifail=3
Repeated errors in an attempted step of underlying ODE solver. Integration was successful as far as ts: ts=value.
ifail=4
In setting up the ODE system an internal auxiliary was unable to initialize the derivative. This could be due to your setting ires=3 in pdedef or bndary.
ifail=5
Singular Jacobian of ODE system. Check problem formulation.
ifail=6
In evaluating residual of ODE system, ires=2 has been set in pdedef or bndary. Integration is successful as far as ts: ts=value.
ifail=7
acc was too small to start integration: acc=value.
ifail=8
ires set to an invalid value in call to pdedef or bndary.
ifail=9
Serious error in internal call to an auxiliary. Increase itrace for further details.
ifail=10
Integration completed, but a small change in acc is unlikely to result in a changed solution. acc=value.
ifail=11
Error during Jacobian formulation for ODE system. Increase itrace for further details.
ifail=14
Flux function appears to depend on time derivatives.
ifail=-99
An unexpected error has been triggered by this routine. Please contact NAG.
See Section 7 in the Introduction to the NAG Library FL Interface for further information.
ifail=-399
Your licence key may have expired or may not have been installed correctly.
See Section 8 in the Introduction to the NAG Library FL Interface for further information.
ifail=-999
Dynamic memory allocation failed.
See Section 9 in the Introduction to the NAG Library FL Interface for further information.

7 Accuracy

d03pdf/​d03pda controls the accuracy of the integration in the time direction but not the accuracy of the approximation in space. The spatial accuracy depends on the degree of the polynomial approximation npoly, and on both the number of break-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 argument, acc.

8 Parallelism and Performance

Background information to multithreading can be found in the Multithreading documentation.
d03pdf/​d03pda is threaded by NAG for parallel execution in multithreaded implementations of the NAG Library.
d03pdf/​d03pda 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

d03pdf/​d03pda is designed to solve parabolic systems (possibly including elliptic equations) with second-order derivatives in space. The argument specification allows you to include equations with only first-order derivatives in the space direction but there is no guarantee that the method of integration will be satisfactory for such systems. The position and nature of the boundary conditions in particular are critical in defining a stable problem.
The time taken depends on the complexity of the parabolic system and on the accuracy requested.

10 Example

The problem consists of a fourth-order PDE which can be written as a pair of second-order elliptic-parabolic PDEs for U1(x,t) and U2(x,t),
0= 2U1 x2 -U2 (4)
U2 t = 2U2 x2 +U2 U1 x -U1 U2 x (5)
where −1x1 and t0. The boundary conditions are given by
U1 x =0  and  U1=1  at ​x=−1,   and U1 x =0  and  U1=−1  at ​x=1.  
The initial conditions at t=0 are given by
U1=-sinπx2  and  U2=π24sinπx2.  
The absence of boundary conditions for U2(x,t) does not pose any difficulties provided that the derivative flux boundary conditions are assigned to the first PDE (4) which has the correct flux, U1 x . The conditions on U1(x,t) at the boundaries are assigned to the second PDE by setting β2=0.0 in equation (3) and placing the Dirichlet boundary conditions on U1(x,t) in the function γ2.

10.1 Program Text

Note: the following programs illustrate the use of d03pdf and d03pda.
Program Text (d03pdfe.f90)
Program Text (d03pdae.f90)

10.2 Program Data

Program Data (d03pdfe.d)
Program Data (d03pdae.d)

10.3 Program Results

Program Results (d03pdfe.r)
Program Results (d03pdae.r)
GnuplotProduced by GNUPLOT 5.0 patchlevel 0 Example Program Solution, U(1,x,t), of Elliptic-parabolic Pair using Chebyshev Collocation and BDF U(1,x,t) gnuplot_plot_1 gnuplot_plot_2 0.000010 0.000100 0.001000 0.010000 0.100000 1.000000 Time (logscale) −1 −0.5 0 0.5 1 x −1 −0.5 0 0.5 1
GnuplotProduced by GNUPLOT 5.0 patchlevel 0 Solution, U(2,x,t), of Elliptic-parabolic Pair using Chebyshev Collocation and BDF U(2,x,t) gnuplot_plot_1 gnuplot_plot_2 0.000010 0.000100 0.001000 0.010000 0.100000 1.000000 Time (logscale) −1 −0.5 0 0.5 1 x −3 −2 −1 0 1 2 3