NAG FL Interface
d03raf (dim2_​gen_​order2_​rectangle)

1 Purpose

d03raf integrates a system of linear or nonlinear, time-dependent partial differential equations (PDEs) in two space dimensions on a rectangular domain. The method of lines is employed to reduce the PDEs to a system of ordinary differential equations (ODEs) which are solved using a backward differentiation formula (BDF) method. The resulting system of nonlinear equations is solved using a modified Newton method and a Bi-CGSTAB iterative linear solver with ILU preconditioning. Local uniform grid refinement is used to improve the accuracy of the solution. d03raf originates from the VLUGR2 package (see Blom and Verwer (1993) and Blom et al. (1996)).

2 Specification

Fortran Interface
Subroutine d03raf ( npde, ts, tout, dt, xmin, xmax, ymin, ymax, nx, ny, tols, tolt, pdedef, bndary, pdeiv, monitr, opti, optr, rwk, lenrwk, iwk, leniwk, lwk, lenlwk, itrace, ind, ifail)
Integer, Intent (In) :: npde, nx, ny, opti(4), lenrwk, leniwk, lenlwk, itrace
Integer, Intent (Inout) :: iwk(leniwk), ind, ifail
Real (Kind=nag_wp), Intent (In) :: tout, xmin, xmax, ymin, ymax, tols, tolt, optr(3,npde)
Real (Kind=nag_wp), Intent (Inout) :: ts, dt(3), rwk(lenrwk)
Logical, Intent (Out) :: lwk(lenlwk)
External :: pdedef, bndary, pdeiv, monitr
C Header Interface
#include <nag.h>
void  d03raf_ (const Integer *npde, double *ts, const double *tout, double dt[], const double *xmin, const double *xmax, const double *ymin, const double *ymax, const Integer *nx, const Integer *ny, const double *tols, const double *tolt,
void (NAG_CALL *pdedef)(const Integer *npts, const Integer *npde, const double *t, const double x[], const double y[], const double u[], const double ut[], const double ux[], const double uy[], const double uxx[], const double uxy[], const double uyy[], double res[]),
void (NAG_CALL *bndary)(const Integer *npts, const Integer *npde, const double *t, const double x[], const double y[], const double u[], const double ut[], const double ux[], const double uy[], const Integer *nbpts, const Integer lbnd[], double res[]),
void (NAG_CALL *pdeiv)(const Integer *npts, const Integer *npde, const double *t, const double x[], const double y[], double u[]),
void (NAG_CALL *monitr)(const Integer *npde, const double *t, const double *dt, const double *dtnew, const logical *tlast, const Integer *nlev, const Integer ngpts[], const double xpts[], const double ypts[], const Integer lsol[], const double sol[], Integer *ierr),
const Integer opti[], const double optr[], double rwk[], const Integer *lenrwk, Integer iwk[], const Integer *leniwk, logical lwk[], const Integer *lenlwk, const Integer *itrace, Integer *ind, Integer *ifail)
The routine may be called by the names d03raf or nagf_pde_dim2_gen_order2_rectangle.

3 Description

d03raf integrates the system of PDEs:
Fjt,x,y,u,ut,ux,uy,uxx,uxy,uyy=0,  j=1,2,,npde, (1)
for x and y in the rectangular domain xminxxmax, yminyymax, and time interval t0ttout, where the vector u is the set of solution values
ux,y,t=u1x,y,t,,unpdex,y,tT,  
and ut denotes partial differentiation with respect to t, and similarly for ux etc.
The functions Fj must be supplied by you in pdedef. Similarly the initial values of the functions ux,y,t must be specified at t=t0 in pdeiv.
Note that whilst complete generality is offered by the master equations (1), d03raf is not appropriate for all PDEs. In particular, hyperbolic systems should not be solved using this routine. Also, at least one component of ut must appear in the system of PDEs.
The boundary conditions must be supplied by you in bndary in the form
Gj t,x,y,u,ut,ux,uy = 0 , (2)
for all y when xmin or xmax and for all x when y=ymin or y=ymax and j=1,2,,npde
The domain is covered by a uniform coarse base grid of size nx×ny specified by you, and nested finer uniform subgrids are subsequently created in regions with high spatial activity. The refinement is controlled using a space monitor which is computed from the current solution and a user-supplied space tolerance tols. A number of optional parameters, e.g., the maximum number of grid levels at any time, and some weighting factors, can be specified in the arrays opti and optr. Further details of the refinement strategy can be found in Section 9.
The system of PDEs and the boundary conditions are discretized in space on each grid using a standard second-order finite difference scheme (centred on the internal domain and one-sided at the boundaries), and the resulting system of ODEs is integrated in time using a second-order, two-step, implicit BDF method with variable step size. The time integration is controlled using a time monitor computed at each grid level from the current solution and a user-supplied time tolerance tolt, and some further optional user-specified weighting factors held in optr (see Section 9 for details). The time monitor is used to compute a new step size, subject to restrictions on the size of the change between steps, and (optional) user-specified maximum and minimum step sizes held in dt. The step size is adjusted so that the remaining integration interval is an integer number times Δt. In this way a solution is obtained at t=tout.
A modified Newton method is used to solve the nonlinear equations arising from the time integration. You may specify (in opti) the maximum number of Newton iterations to be attempted. A Jacobian matrix is calculated at the beginning of each time step. If the Newton process diverges or the maximum number of iterations is exceeded, a new Jacobian is calculated using the most recent iterates and the Newton process is restarted. If convergence is not achieved after the (optional) user-specified maximum number of new Jacobian evaluations, the time step is retried with Δt=Δt/4. The linear systems arising from the Newton iteration are solved using a Bi-CGSTAB iterative method, in combination with ILU preconditioning. The maximum number of iterations can be specified by you in opti.
The solution at all grid levels is stored in the workspace arrays, along with other information needed for a restart (i.e., a continuation call). It is not intended that you extract the solution from these arrays, indeed the necessary information regarding these arrays is not included. The user-supplied monitor monitr should be used to obtain the solution at particular levels and times. monitr is called at the end of every time step, with the last step being identified via the input argument tlast.
Within pdeiv, pdedef, bndary and monitr the data structure is as follows. Each point on a particular grid is given an index (ranging from 1 to the total number of points on the grid) and all coordinate or solution information is stored in arrays according to this index, e.g., xi and yi contain the x- and y coordinate of point i, and uij contains the jth solution component uj at point i.
Further details of the underlying algorithm can be found in Section 9 and in Blom and Verwer (1993) and Blom et al. (1996) and the references therein.

4 References

Adjerid S and Flaherty J E (1988) A local refinement finite element method for two-dimensional parabolic systems SIAM J. Sci. Statist. Comput. 9 792–811
Blom J G, Trompert R A and Verwer J G (1996) Algorithm 758. VLUGR2: A vectorizable adaptive grid solver for PDEs in 2D Trans. Math. Software 22 302–328
Blom J G and Verwer J G (1993) VLUGR2: A vectorized local uniform grid refinement code for PDEs in 2D Report NM-R9306 CWI, Amsterdam
Brown P N, Hindmarsh A C and Petzold L R (1994) Using Krylov methods in the solution of large scale differential-algebraic systems SIAM J. Sci. Statist. Comput. 15 1467–1488
Trompert R A (1993) Local uniform grid refinement and systems of coupled partial differential equations Appl. Numer. Maths 12 331–355
Trompert R A and Verwer J G (1993) Analysis of the implicit Euler local uniform grid refinement method SIAM J. Sci. Comput. 14 259–278

5 Arguments

1: npde Integer Input
On entry: the number of PDEs in the system.
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 which has been reached. 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: dt3 Real (Kind=nag_wp) array Input/Output
On entry: the initial, minimum and maximum time step sizes respectively.
dt1
Specifies the initial time step size to be used on the first entry, i.e., when ind=0. If dt1=0.0 then the default value dt1=0.01×tout-ts is used. On subsequent entries (ind=1), the value of dt1 is not referenced.
dt2
Specifies the minimum time step size to be attempted by the integrator. If dt2=0.0 the default value dt2=10.0×machine precision is used.
dt3
Specifies the maximum time step size to be attempted by the integrator. If dt3=0.0 the default value dt3=tout-ts is used.
On exit: dt1 contains the time step size for the next time step. dt2 and dt3 are unchanged or set to their default values if zero on entry.
Constraints:
  • if ind=0, dt10.0;
  • if ind=0 and dt1>0.0, 10.0×machine precision×maxts,toutdt1tout-ts and dt2dt1dt3, where the values of dt2 and dt3 will have been reset to their default values if zero on entry;
  • 0 dt2 dt3 .
5: xmin Real (Kind=nag_wp) Input
6: xmax Real (Kind=nag_wp) Input
On entry: the extents of the rectangular domain in the x-direction, i.e., the x coordinates of the left and right boundaries respectively.
Constraint: xmin<xmax and xmax must be sufficiently distinguishable from xmin for the precision of the machine being used.
7: ymin Real (Kind=nag_wp) Input
8: ymax Real (Kind=nag_wp) Input
On entry: the extents of the rectangular domain in the y-direction, i.e., the y coordinates of the lower and upper boundaries respectively.
Constraint: ymin<ymax and ymax must be sufficiently distinguishable from ymin for the precision of the machine being used.
9: nx Integer Input
On entry: the number of grid points in the x-direction (including the boundary points).
Constraint: nx4.
10: ny Integer Input
On entry: the number of grid points in the y-direction (including the boundary points).
Constraint: ny4.
11: tols Real (Kind=nag_wp) Input
On entry: the space tolerance used in the grid refinement strategy (σ in equation (4)). See Section 9.2.
Constraint: tols>0.0.
12: tolt Real (Kind=nag_wp) Input
On entry: the time tolerance used to determine the time step size (τ in equation (7)). See Section 9.3.
Constraint: tolt>0.0.
13: pdedef Subroutine, supplied by the user. External Procedure
pdedef must evaluate the functions Fj, for j=1,2,,npde, in equation (1) which define the system of PDEs (i.e., the residuals of the resulting ODE system) at all interior points of the domain. Values at points on the boundaries of the domain are ignored and will be overwritten by bndary. pdedef is called for each subgrid in turn.
The specification of pdedef is:
Fortran Interface
Subroutine pdedef ( npts, npde, t, x, y, u, ut, ux, uy, uxx, uxy, uyy, res)
Integer, Intent (In) :: npts, npde
Real (Kind=nag_wp), Intent (In) :: t, x(npts), y(npts), u(npts,npde), ut(npts,npde), ux(npts,npde), uy(npts,npde), uxx(npts,npde), uxy(npts,npde), uyy(npts,npde)
Real (Kind=nag_wp), Intent (Out) :: res(npts,npde)
C Header Interface
void  pdedef_ (const Integer *npts, const Integer *npde, const double *t, const double x[], const double y[], const double u[], const double ut[], const double ux[], const double uy[], const double uxx[], const double uxy[], const double uyy[], double res[])
1: npts Integer Input
On entry: the number of grid points in the current grid.
2: npde Integer Input
On entry: the number of PDEs in the system.
3: t Real (Kind=nag_wp) Input
On entry: the current value of the independent variable t.
4: xnpts Real (Kind=nag_wp) array Input
On entry: xi contains the x coordinate of the ith grid point, for i=1,2,,npts.
5: ynpts Real (Kind=nag_wp) array Input
On entry: yi contains the y coordinate of the ith grid point, for i=1,2,,npts.
6: unptsnpde Real (Kind=nag_wp) array Input
On entry: uij contains the value of the jth PDE component at the ith grid point, for i=1,2,,npts and j=1,2,,npde.
7: utnptsnpde Real (Kind=nag_wp) array Input
On entry: utij contains the value of u t for the jth PDE component at the ith grid point, for i=1,2,,npts and j=1,2,,npde.
8: uxnptsnpde Real (Kind=nag_wp) array Input
On entry: uxij contains the value of u x for the jth PDE component at the ith grid point, for i=1,2,,npts and j=1,2,,npde.
9: uynptsnpde Real (Kind=nag_wp) array Input
On entry: uyij contains the value of u y for the jth PDE component at the ith grid point, for i=1,2,,npts and j=1,2,,npde.
10: uxxnptsnpde Real (Kind=nag_wp) array Input
On entry: uxxij contains the value of 2u x2 for the jth PDE component at the ith grid point, for i=1,2,,npts and j=1,2,,npde.
11: uxynptsnpde Real (Kind=nag_wp) array Input
On entry: uxyij contains the value of 2u xy for the jth PDE component at the ith grid point, for i=1,2,,npts and j=1,2,,npde.
12: uyynptsnpde Real (Kind=nag_wp) array Input
On entry: uyyij contains the value of 2u y2 for the jth PDE component at the ith grid point, for i=1,2,,npts and j=1,2,,npde.
13: resnptsnpde Real (Kind=nag_wp) array Output
On exit: resij must contain the value of Fj, for j=1,2,,npde, at the ith grid point, for i=1,2,,npts, although the residuals at boundary points will be ignored (and overwritten later on) and so they need not be specified here.
pdedef must either be a module subprogram USEd by, or declared as EXTERNAL in, the (sub)program from which d03raf 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 d03raf. If your code inadvertently does return any NaNs or infinities, d03raf is likely to produce unexpected results.
14: bndary Subroutine, supplied by the user. External Procedure
bndary must evaluate the functions Gj, for j=1,2,,npde, in equation (2) which define the boundary conditions at all boundary points of the domain. Residuals at interior points must not be altered by this subroutine.
The specification of bndary is:
Fortran Interface
Subroutine bndary ( npts, npde, t, x, y, u, ut, ux, uy, nbpts, lbnd, res)
Integer, Intent (In) :: npts, npde, nbpts, lbnd(nbpts)
Real (Kind=nag_wp), Intent (In) :: t, x(npts), y(npts), u(npts,npde), ut(npts,npde), ux(npts,npde), uy(npts,npde)
Real (Kind=nag_wp), Intent (Inout) :: res(npts,npde)
C Header Interface
void  bndary_ (const Integer *npts, const Integer *npde, const double *t, const double x[], const double y[], const double u[], const double ut[], const double ux[], const double uy[], const Integer *nbpts, const Integer lbnd[], double res[])
1: npts Integer Input
On entry: the number of grid points in the current grid.
2: npde Integer Input
On entry: the number of PDEs in the system.
3: t Real (Kind=nag_wp) Input
On entry: the current value of the independent variable t.
4: xnpts Real (Kind=nag_wp) array Input
On entry: xi contains the x coordinate of the ith grid point, for i=1,2,,npts.
5: ynpts Real (Kind=nag_wp) array Input
On entry: yi contains the y coordinate of the ith grid point, for i=1,2,,npts.
6: unptsnpde Real (Kind=nag_wp) array Input
On entry: uij contains the value of the jth PDE component at the ith grid point, for i=1,2,,npts and j=1,2,,npde.
7: utnptsnpde Real (Kind=nag_wp) array Input
On entry: utij contains the value of u t for the jth PDE component at the ith grid point, for i=1,2,,npts and j=1,2,,npde.
8: uxnptsnpde Real (Kind=nag_wp) array Input
On entry: uxij contains the value of u x for the jth PDE component at the ith grid point, for i=1,2,,npts and j=1,2,,npde.
9: uynptsnpde Real (Kind=nag_wp) array Input
On entry: uyij contains the value of u y for the jth PDE component at the ith grid point, for i=1,2,,npts and j=1,2,,npde.
10: nbpts Integer Input
On entry: the number of boundary points in the grid.
11: lbndnbpts Integer array Input
On entry: lbndi contains the grid index for the ith boundary point, for i=1,2,,nbpts. Hence the ith boundary point has coordinates xlbndi and ylbndi, and the corresponding solution values are ulbndinpde, etc.
12: resnptsnpde Real (Kind=nag_wp) array Input/Output
On entry: resij contains the value of Fj, for i=1,2,,npde, at the ith grid point, for i=1,2,,npts, as returned by pdedef. The residuals at the boundary points will be overwritten and so need not have been set by pdedef.
On exit: reslbndij must contain the value of Gj, for j=1,2,,npde, at the ith boundary point, for i=1,2,,nbpts.
Note: elements of res corresponding to interior points must not be altered.
bndary must either be a module subprogram USEd by, or declared as EXTERNAL in, the (sub)program from which d03raf 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 d03raf. If your code inadvertently does return any NaNs or infinities, d03raf is likely to produce unexpected results.
15: pdeiv Subroutine, supplied by the user. External Procedure
pdeiv must specify the initial values of the PDE components u at all points in the grid. pdeiv is not referenced if, on entry, ind=1.
The specification of pdeiv is:
Fortran Interface
Subroutine pdeiv ( npts, npde, t, x, y, u)
Integer, Intent (In) :: npts, npde
Real (Kind=nag_wp), Intent (In) :: t, x(npts), y(npts)
Real (Kind=nag_wp), Intent (Out) :: u(npts,npde)
C Header Interface
void  pdeiv_ (const Integer *npts, const Integer *npde, const double *t, const double x[], const double y[], double u[])
1: npts Integer Input
On entry: the number of grid points in the grid.
2: npde Integer Input
On entry: the number of PDEs in the system.
3: t Real (Kind=nag_wp) Input
On entry: the (initial) value of the independent variable t.
4: xnpts Real (Kind=nag_wp) array Input
On entry: xi contains the x coordinate of the ith grid point, for i=1,2,,npts.
5: ynpts Real (Kind=nag_wp) array Input
On entry: yi contains the y coordinate of the ith grid point, for i=1,2,,npts.
6: unptsnpde Real (Kind=nag_wp) array Output
On exit: uij must contain the value of the jth PDE component at the ith grid point, for i=1,2,,npts and j=1,2,,npde.
pdeiv must either be a module subprogram USEd by, or declared as EXTERNAL in, the (sub)program from which d03raf is called. Arguments denoted as Input must not be changed by this procedure.
Note: pdeiv should not return floating-point NaN (Not a Number) or infinity values, since these are not handled by d03raf. If your code inadvertently does return any NaNs or infinities, d03raf is likely to produce unexpected results.
16: monitr Subroutine, supplied by the user. External Procedure
monitr is called by d03raf at the end of every successful time step, and may be used to examine or print the solution or perform other tasks such as error calculations, particularly at the final time step, indicated by the argument tlast. The input arguments contain information about the grid and solution at all grid levels used.
monitr can also be used to force an immediate tidy termination of the solution process and return to the calling program.
The specification of monitr is:
Fortran Interface
Subroutine monitr ( npde, t, dt, dtnew, tlast, nlev, ngpts, xpts, ypts, lsol, sol, ierr)
Integer, Intent (In) :: npde, nlev, ngpts(nlev), lsol(nlev)
Integer, Intent (Inout) :: ierr
Real (Kind=nag_wp), Intent (In) :: t, dt, dtnew, xpts(*), ypts(*), sol(*)
Logical, Intent (In) :: tlast
C Header Interface
void  monitr_ (const Integer *npde, const double *t, const double *dt, const double *dtnew, const logical *tlast, const Integer *nlev, const Integer ngpts[], const double xpts[], const double ypts[], const Integer lsol[], const double sol[], Integer *ierr)
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, i.e., the time at the end of the integration step just completed.
3: dt Real (Kind=nag_wp) Input
On entry: the current time step size Δt, i.e., the time step size used for the integration step just completed.
4: dtnew Real (Kind=nag_wp) Input
On entry: the step size that will be used for the next time step.
5: tlast Logical Input
On entry: indicates if intermediate or final time step. tlast=.FALSE. for an intermediate step, tlast=.TRUE. for the last call to monitr before returning to your program.
6: nlev Integer Input
On entry: the number of grid levels used at time t.
7: ngptsnlev Integer array Input
On entry: ngptsl contains the number of grid points at level l, for l=1,2,,nlev.
8: xpts* Real (Kind=nag_wp) array Input
On entry: contains the x coordinates of the grid points in each level in turn, i.e., xi, for i=1,2,,ngptsl and l=1,2,,nlev.
So for level l, xi=xptsk+i, where k=ngpts1+ngpts2++ngptsl-1, for i=1,2,,ngptsl and l=1,2,,nlev.
9: ypts* Real (Kind=nag_wp) array Input
On entry: contains the y coordinates of the grid points in each level in turn, i.e., yi, for i=1,2,,ngptsl and l=1,2,,nlev.
So for level l, yi=yptsk+i, where k=ngpts1+ngpts2++ngptsl-1, for i=1,2,,ngptsl and l=1,2,,nlev.
10: lsolnlev Integer array Input
On entry: lsoll contains the pointer to the solution in sol at grid level l and time t. (lsoll actually contains the array index immediately preceding the start of the solution in sol.)
11: sol* Real (Kind=nag_wp) array Input
On entry: contains the solution ungptslnpde at time t for each grid level l in turn, positioned according to lsol, i.e., for level l, uij=sollsoll+j-1×ngptsl+i, for i=1,2,,ngptsl, j=1,2,,npde and l=1,2,,nlev.
12: ierr Integer Input/Output
On entry: will be set to 0.
On exit: should be set to 1 to force a tidy termination and an immediate return to the calling program with ifail=4. ierr should remain unchanged otherwise.
monitr must either be a module subprogram USEd by, or declared as EXTERNAL in, the (sub)program from which d03raf is called. Arguments denoted as Input must not be changed by this procedure.
17: opti4 Integer array Input
On entry: may be set to control various options available in the integrator.
opti1=0
All the default options are employed.
opti1>0
The default value of optii, for i=2,3,4, can be obtained by setting optii=0.
opti1
Specifies the maximum number of grid levels allowed (including the base grid). opti10. The default value is opti1=3.
opti2
Specifies the maximum number of Jacobian evaluations allowed during each nonlinear equations solution. opti20. The default value is opti2=2.
opti3
Specifies the maximum number of Newton iterations in each nonlinear equations solution. opti30. The default value is opti3=10.
opti4
Specifies the maximum number of iterations in each linear equations solution. opti40. The default value is opti4=100.
Constraint: opti10 and if opti1>0, optii0, for i=2,3,4.
18: optr3npde Real (Kind=nag_wp) array Input
On entry: may be used to specify the optional vectors umax, ws and wt in the space and time monitors (see Section 9).
If an optional vector is not required then all its components should be set to 1.0.
optr1j, for j=1,2,,npde, specifies ujmax, the approximate maximum absolute value of the jth component of u, as used in (4) and (7). optr1j>0.0, for j=1,2,,npde.
optr2j, for j=1,2,,npde, specifies wjs, the weighting factors used in the space monitor (see (4)) to indicate the relative importance of the jth component of u on the space monitor. optr2j0.0, for j=1,2,,npde.
optr3j, for j=1,2,,npde, specifies wjt, the weighting factors used in the time monitor (see (6)) to indicate the relative importance of the jth component of u on the time monitor. optr3j0.0, for j=1,2,,npde.
Constraints:
  • optr1j>0.0, for j=1,2,,npde;
  • optrij0.0, for i=2,3 and j=1,2,,npde.
19: rwklenrwk Real (Kind=nag_wp) array Communication Array
20: lenrwk Integer Input
On entry: the dimension of the array rwk as declared in the (sub)program from which d03raf is called.
The required value of lenrwk cannot be determined exactly in advance, but a suggested value is
lenrwk=maxpts×npde×5×l+18×npde+9+2×maxpts,  
where l=opti1 if opti10 and l=3 otherwise, and maxpts is the expected maximum number of grid points at any one level. If during the execution the supplied value is found to be too small then the routine returns with ifail=3 and an estimated required size is printed on the current error message unit (see x04aaf).
Constraint: lenrwknx×ny×npde×14+18×npde+2×nx×ny (the required size for the initial grid).
21: iwkleniwk Integer array Communication Array
On entry: if ind=0, iwk need not be set. Otherwise iwk must remain unchanged from a previous call to d03raf.
On exit: the following components of the array iwk concern the efficiency of the integration. Here, m is the maximum number of grid levels allowed (m=opti1 if opti1>1 and m=3 otherwise), and l is a grid level taking the values l=1,2,,nl, where nl is the number of levels used.
iwk1
Contains the number of steps taken in time.
iwk2
Contains the number of rejected time steps.
iwk2+l
Contains the total number of residual evaluations performed (i.e., the number of times pdedef was called) at grid level l.
iwk2+m+l
Contains the total number of Jacobian evaluations performed at grid level l.
iwk2+2×m+l
Contains the total number of Newton iterations performed at grid level l.
iwk2+3×m+l
Contains the total number of linear solver iterations performed at grid level l.
iwk2+4×m+l
Contains the maximum number of Newton iterations performed at any one time step at grid level l.
iwk2+5×m+l
Contains the maximum number of linear solver iterations performed at any one time step at grid level l.
Note: the total and maximum numbers are cumulative over all calls to d03raf. If the specified maximum number of Newton or linear solver iterations is exceeded at any stage, the maximums above are set to the specified maximum plus one.
22: leniwk Integer Input
On entry: the dimension of the array iwk as declared in the (sub)program from which d03raf is called.
The required value of leniwk cannot be determined exactly in advance, but a suggested value is
leniwk=maxpts×14+5×m+7×m+2, 
where maxpts is the expected maximum number of grid points at any one level and m=opti1 if opti1>0 and m=3 otherwise. If during the execution the supplied value is found to be too small then the routine returns with ifail=3 and an estimated required size is printed on the current error message unit (see x04aaf).
Constraint: leniwk19×nx×ny+9 (the required size for the initial grid).
23: lwklenlwk Logical array Workspace
24: lenlwk Integer Input
On entry: the dimension of the array lwk as declared in the (sub)program from which d03raf is called.
The required value of lenlwk cannot be determined exactly in advanced, but a suggested value is
lenlwk=maxpts+1, 
where maxpts is the expected maximum number of grid points at any one level. If during the execution the supplied value is found to be too small then the routine returns with ifail=3 and an estimated required size is printed on the current error message unit (see x04aaf).
Constraint: lenlwknx×ny+1 (the required size for the initial grid).
25: itrace Integer Input
On entry: the level of trace information required from d03raf. itrace may take the value -1, 0, 1, 2 or 3.
itrace=-1
No output is generated.
itrace=0
Only warning messages are printed.
itrace>0
Output from the underlying solver is printed on the current advisory message unit (see x04abf). This output contains details of the time integration, the nonlinear iteration and the linear solver.
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. Setting itrace=1 allows you to monitor the progress of the integration without possibly excessive information.
26: ind Integer Input/Output
On entry: must be set to 0 or 1, alternatively 10 or 11.
ind=0
Starts the integration in time. pdedef is assumed to be serial.
ind=1
Continues the integration after an earlier exit from the routine. In this case, only the following parameters may be reset between calls to d03raf: tout, dt, tols, tolt, opti, optr, itrace and ifail. pdedef is assumed to be serial.
ind=10
Starts the integration in time. pdedef is assumed to have been parallelized by you, as described in Section 8. In all other respects, this is equivalent to ind=0.
ind=11
Continues the integration after an earlier exit from the routine. In this case, only the following parameters may be reset between calls to d03raf: tout, dt, tols, tolt, opti, optr, itrace and ifail. pdedef is assumed to have been parallelized by you, as described in Section 8. In all other respects, this is equivalent to ind=1.
Constraint: 0ind1 or 10ind11.
On exit: ind=1, if ind on input was 0 or 1, or ind=11, if ind on input was 10 or 11.
Note: for users of serial versions of the NAG Library, it is recommended that you only use ind=0 or 1. See Section 8 for more information on the use of ind.
27: ifail Integer Input/Output
On entry: ifail must be set to 0, -1 or 1. If you are unfamiliar with this argument you should refer to Section 4 in the Introduction to the NAG Library FL Interface for details.
For environments where it might be inappropriate to halt program execution when an error is detected, the value -1 or 1 is recommended. If the output of error messages is undesirable, then the value 1 is recommended. Otherwise, if you are not familiar with this argument, the recommended value is 0. When the value -1 or 1 is used it is essential to test the value of ifail on exit.
On exit: ifail=0 unless the routine detects an error or a warning has been flagged (see Section 6).

6 Error Indicators and Warnings

If on entry ifail=0 or -1, explanatory error messages are output on the current error message unit (as defined by x04aaf).
Errors or warnings detected by the routine:
ifail=1
On entry, dt1=value.
Constraint: if ind=0, dt10.0.
On entry, dt1=value and dt2=value. Note that dt2 was reset to default if zero on entry.
Constraint: if ind=0, dt1dt2.
On entry, dt1=value and dt3=value. Note that dt3 was reset to default if zero on entry.
Constraint: if ind=0, dt1dt3.
On entry, dt2=value.
Constraint: dt20.0.
On entry, dt2=value and dt3=value.
Constraint: dt2dt3.
On entry, dt3=value.
Constraint: dt30.0.
On entry, i=value, j=value and optrij=value.
Constraint: optrij0.0.
On entry, ind=0 and dt1 too large: dt1=value and maximum value=value.
On entry, ind=0 and dt1 too small: dt1=value and minimum value=value.
On entry, ind is not equal to 0 or 1: ind=value.
On entry, j=value and optr1j=value.
Constraint: optr1j>0.0.
On entry, leniwk too small for initial grid level: leniwk=value, minimum value =value. Note that subsequent levels will require more. Consult document.
On entry, lenlwk too small for initial grid level: lenlwk=value, minimum value =value. Note that subsequent levels will require more. Consult document.
On entry, lenrwk too small for initial grid level: lenrwk=value, minimum value =value. Note that subsequent levels will require more. Consult document.
On entry, npde=value.
Constraint: npde1.
On entry, nx=value.
Constraint: nx4.
On entry, ny=value.
Constraint: ny4.
On entry, opti1=value, i=value and optii=value.
Constraint: if opti1>0, optii0.
On entry, opti1=value.
Constraint: opti10.
On entry, tols=value.
Constraint: tols>0.0.
On entry, tolt=value.
Constraint: tolt>0.0.
On entry, tout=value and ts=value.
Constraint: tout>ts.
On entry, tout-ts too small: tout-ts=value.
On entry, xmax too close to xmin: xmax=value and xmin=value.
On entry, xmin=value and xmax=value.
Constraint: xmax>xmin.
On entry, ymax too close to ymin: ymax=value and ymin=value.
On entry, ymin=value and ymax=value.
Constraint: ymax>ymin.
ifail=2
Attempted time-step smaller than specified minimum. Check problem formulation in pdedef, bndary and pdeiv. Try increasing itrace for more information.
The time step size to be attempted is less than the specified minimum size. This may occur following time step failures and subsequent step size reductions caused by one or more of the following:
  • the requested accuracy could not be achieved, i.e., tolt is too small,
  • the maximum number of linear solver iterations, Newton iterations or Jacobian evaluations is too small,
  • ILU decomposition of the Jacobian matrix could not be performed, possibly due to singularity of the Jacobian.
Setting itrace to a higher value may provide further information.
In the latter two cases you are advised to check their problem formulation in pdedef and/or bndary, and the initial values in pdeiv if appropriate.
ifail=3
One or more of the workspace arrays are too small. Try increasing itrace for more information.
ifail=4
IERR set to 1 in monitr. Integration completed as far as ts: ts=value.
ifail=5
Integration completed, but maximum number of levels too small for required accuracy.
The integration has been completed but the maximum number of levels specified in opti1 was insufficient at one or more time steps, meaning that the requested space accuracy could not be achieved. To avoid this warning either increase the value of opti1 or decrease the value of tols.
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

There are three sources of error in the algorithm: space and time discretization, and interpolation (linear) between grid levels. The space and time discretization errors are controlled separately using the arguments tols and tolt described in the following section, and you should test the effects of varying these arguments. Interpolation errors are generally implicitly controlled by the refinement criterion since in areas where interpolation errors are potentially large, the space monitor will also be large. It can be shown that the global spatial accuracy is comparable to that which would be obtained on a uniform grid of the finest grid size. A full error analysis can be found in Trompert and Verwer (1993).

8 Parallelism and Performance

d03raf is not thread safe and should not be called from a multithreaded user program. Please see Section 1 in FL Interface Multithreading for more information on thread safety.
d03raf is threaded by NAG for parallel execution in multithreaded implementations of the NAG Library.
d03raf 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.
d03raf requires a user-supplied routine pdedef to evaluate the functions Fj, for j=1,2,,npde. The parallelism within d03raf will be more efficient if pdedef can also be parallelized. This is often the case, but you must add some OpenMP directives to your version of pdedef to implement the parallelism. For example, if the body of code for pdedef is as follows (adapted from the first test case in the document for d03raf):
res(1:npts,1:npde) = ut(1:npts,1:npde) - diffusion*(uxx(1:npts,1:         &
  npde)+uyy(1:npts,1:npde)) - damkohler*(one+heat_release-u(1:npts,       &
  1:npde))*exp(-activ_energy/u(1:npts,1:npde))
This example can be parallelized, as the updating of res for each value in the range 1,,npts is independent of every other value. Thus this should be parallelized in OpenMP (using an explicit loop rather than Fortran array syntax) as follows:
!$OMP DO
   Do i = 1, npts
      res(i,1:npde) = ut(i,1:npde) -diffusion*(uxx(i,1:npde)+uyy(i,1:npde &
       )) - damkohler*(1.0E0_nag_wp+heat_release-u(i,1:npde))*exp(-       &
        activ_energy/u(i,1:npde))
   End Do
!$OMP END DO
Note that the OpenMP PARALLEL directive must not be specified, as the OpenMP DO directive will bind to the PARALLEL region within the d03raf code. Also note that this assumes the default OpenMP behaviour that all variables are SHARED, except for loop indices that are PRIVATE.
To avoid problems for existing Library users, who will not have specified any OpenMP directives in their pdedef routine, the default assumption of d03raf is that pdedef has not been parallelized, and executes calls to pdedef in serial mode. You must indicate that pdedef has been parallelized by setting ind to 10 or 11 as appropriate. See Section 5 for details.
If the code within pdedef cannot be parallelized, you must not add any OpenMP directives to your code, and must not set ind to 10 or 11. If ind is set to 10 or 11 and pdedef has not been parallelized, results on multiple threads will be unpredictable and may give rise to incorrect results and/or program crashes or deadlocks. Please contact NAG for advice if required. Overloading ind in this manner is not entirely satisfactory, consequently it is likely that replacement interfaces for d03raf will be included in a future NAG Library release.

9 Further Comments

9.1 Algorithm Outline

The local uniform grid refinement method is summarised as follows:
  1. 1.Initialize the course base grid, an initial solution and an initial time step.
  2. 2.Solve the system of PDEs on the current grid with the current time step.
  3. 3.If the required accuracy in space and the maximum number of grid levels have not yet been reached:
    1. (a)Determine new finer grid at forward time level.
    2. (b)Get solution values at previous time level(s) on new grid.
    3. (c)Interpolate internal boundary values from old grid at forward time.
    4. (d)Get initial values for the Newton process at forward time.
    5. (e)Go to 2.
  4. 4.Update the coarser grid solution using the finer grid values.
  5. 5.Estimate error in time integration. If time error is acceptable advance time level.
  6. 6.Determine new step size then go to 2 with coarse base as current grid.

9.2 Refinement Strategy

For each grid point i a space monitor μis is determined by
μ i s = max j=1,npde γ j Δ x 2 2 x 2 u j x i , y i ,t + Δ y 2 2 y 2 u j x i , y i ,t , (3)
where Δx and Δy are the grid widths in the x and y directions; and xi, yi are the x and y coordinates at grid point i. The argument γj is obtained from
γj=wjs ujmax σ , (4)
where σ is the user-supplied space tolerance; wjs is a weighting factor for the relative importance of the jth PDE component on the space monitor; and ujmax is the approximate maximum absolute value of the jth component. A value for σ must be supplied by you. Values for wjs and ujmax must also be supplied but may be set to the value 1.0 if little information about the solution is known.
A new level of refinement is created if
maxiμis>0.9  or  1.0, (5)
depending on the grid level at the previous step in order to avoid fluctuations in the number of grid levels between time steps. If (5) is satisfied then all grid points for which μis>0.25 are flagged and surrounding cells are quartered in size.
No derefinement takes place as such, since at each time step the solution on the base grid is computed first and new finer grids are then created based on the new solution. Hence derefinement occurs implicitly. See Section 9.1.

9.3 Time Integration

The time integration is controlled using a time monitor calculated at each level l up to the maximum level used, given by
μlt=1Nj=1npdewjti=1 ngptsl Δtαij utxi,yi,t 2 (6)
where ngptsl is the total number of points on grid level l; N=ngptsl×npde; Δt is the current time step; ut is the time derivative of u which is approximated by first-order finite differences; wjt is the time equivalent of the space weighting factor wjs; and αij is given by
αij=τ ujmax100+uxi,yi,t (7)
where ujmax is as before, and τ is the user-specified time tolerance.
An integration step is rejected and retried at all levels if
maxlμlt>1.0. (8)

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 stems from combustion theory and is a model for a single, one-step reaction of a mixture of two chemicals (see Adjerid and Flaherty (1988)). The PDE for the temperature of the mixture u is
u t =d 2u x2 + 2u y2 +D1+α-u exp-δu  
for 0x,y1 and t0, with initial conditions ux,y,0=1 for 0x,y1, and boundary conditions
ux 0,y,t = 0, u 1,y,t = 1   for ​ 0y1,  
uy x,0,t = 0, u x,1,t = 1   for ​ 0x 1.  
The heat release argument α=1, the Damkohler number D=Rexpδ/αδ, the activation energy δ=20, the reaction rate R=5, and the diffusion argument d=0.1.
For small times the temperature gradually increases in a circular region about the origin, and at about t=0.24 ‘ignition’ occurs causing the temperature to suddenly jump from near unity to 1+α, and a reaction front forms and propagates outwards, becoming steeper. Thus during the solution, just one grid level is used up to the ignition point, then two levels, and then three as the reaction front steepens.
Example 2 (EX2)
This example is taken from a multispecies food web model, in which predator-prey relationships in a spatial domain are simulated (see Brown et al. (1994)). In this example there is just one species each of prey and predator, and the two PDEs for the concentrations c1 and c2 of the prey and the predator respectively are
c1 t =c1b1+a11c1+a12c2+d1 2c1 x2 + 2c1 y2 ,  
0=c2 b2+a21c1+a22c2+d2 2 c2 x2 + 2 c2 y2 ,  
with
a11=a22=-1, a12=-0.5×10-6, and a21=104, and b1=1+αxy+βsin4πxsin4πy,  
where α=50 and β=300, and b2=-b1.
The initial conditions are taken to be simple peaked functions which satisfy the boundary conditions and very nearly satisfy the PDEs:
c1=10+16x1-xy1-y2,  
c2=b2+a21c1,  
and the boundary conditions are of Neumann type, i.e., zero normal derivatives everywhere.
During the solution a number of peaks and troughs develop across the domain, and so the number of levels required increases with time. Since the solution varies rapidly in space across the whole of the domain, refinement at intermediate levels tends to occur at all points of the domain.

10.1 Program Text

Program Text (d03rafe.f90)

10.2 Program Data

Program Data (d03rafe.d)

10.3 Program Results

Program Results (d03rafe.r)
GnuplotProduced by GNUPLOT 4.6 patchlevel 3 Example Program 1 Model for a Single, One-step Reaction of a Mixture of Two Chemicals U(x,y,t) gnuplot_plot_1 gnuplot_plot_2 t = 0.001 gnuplot_plot_3 t = 0.228 gnuplot_plot_4 t = 0.240 gnuplot_plot_5 t = 0.25 0 0.2 0.4 0.6 0.8 1 y 0 0.2 0.4 0.6 0.8 1 x 1 1.2 1.4 1.6 1.8 2
GnuplotProduced by GNUPLOT 4.6 patchlevel 3 Example Program 2 Multispecies Food Web Model Concentrations of Predator Predator Concentration gnuplot_plot_1 gnuplot_plot_2 t = 0.0005 gnuplot_plot_3 t = 0.0076 gnuplot_plot_4 t = 0.0173 gnuplot_plot_5 t = 0.025 0 0.2 0.4 0.6 0.8 1 y 0 0.2 0.4 0.6 0.8 1 x 0 20 40 60 80 100 120 140
GnuplotProduced by GNUPLOT 4.6 patchlevel 3 Multispecies Food Web Model Concentrations of Prey Prey Concentration gnuplot_plot_1 gnuplot_plot_2 t = 0.0005 gnuplot_plot_3 t = 0.0076 gnuplot_plot_4 t = 0.0173 gnuplot_plot_5 t = 0.025 0 0.2 0.4 0.6 0.8 1 y 0 0.2 0.4 0.6 0.8 1 x 0 200000 400000 600000 800000 1e+06 1.2e+06 1.4e+06