hide long namesshow long names
hide short namesshow short names
Integer type:  int32  int64  nag_int  show int32  show int32  show int64  show int64  show nag_int  show nag_int

PDF version (NAG web site, 64-bit version, 64-bit version)
Chapter Contents
Chapter Introduction
NAG Toolbox

NAG Toolbox: nag_pde_2d_gen_order2_rectilinear (d03rb)

 Contents

    1  Purpose
    2  Syntax
    7  Accuracy
    9  Example

Purpose

nag_pde_2d_gen_order2_rectilinear (d03rb) integrates a system of linear or nonlinear, time-dependent partial differential equations (PDEs) in two space dimensions on a rectilinear domain. The method of lines is employed to reduce the NPDEs 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. nag_pde_2d_gen_order2_rectilinear (d03rb) originates from the VLUGR2 package (see Blom and Verwer (1993) and Blom et al. (1996)).

Syntax

[ts, dt, rwk, iwk, ind, ifail] = d03rb(ts, tout, dt, tols, tolt, inidom, pdedef, bndary, pdeiv, monitr, opti, optr, rwk, iwk, lenlwk, itrace, ind, 'npde', npde, 'lenrwk', lenrwk, 'leniwk', leniwk)
[ts, dt, rwk, iwk, ind, ifail] = nag_pde_2d_gen_order2_rectilinear(ts, tout, dt, tols, tolt, inidom, pdedef, bndary, pdeiv, monitr, opti, optr, rwk, iwk, lenlwk, itrace, ind, 'npde', npde, 'lenrwk', lenrwk, 'leniwk', leniwk)

Description

nag_pde_2d_gen_order2_rectilinear (d03rb) integrates the system of PDEs:
Fj t,x,y,u,ut,ux,uy,uxx,uxy,uyy = 0 ,   j=1,2,,npde ,   x,y Ω ,   t0 t tout , (1)
where Ω is an arbitrary rectilinear domain, i.e., a domain bounded by perpendicular straight lines. If the domain is rectangular then it is recommended that nag_pde_2d_gen_order2_rectangle (d03ra) is used.
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 for x,yΩ must be specified at t=t0 in pdeiv.
Note that whilst complete generality is offered by the master equations (1), nag_pde_2d_gen_order2_rectilinear (d03rb) is not appropriate for all PDEs. In particular, hyperbolic systems should not be solved using this function. 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 ,   j=1,2,,npde ,   x,yΩ ,   t0ttout . (2)
The domain is covered by a uniform coarse base grid 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 Further Comments.
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 Further Comments 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.
In order to define the base grid you must first specify a virtual uniform rectangular grid which contains the entire base grid. The position of the virtual grid in physical x,y space is given by the x,y coordinates of its boundaries. The number of points nx and ny in the x and y directions must also be given, corresponding to the number of columns and rows respectively. This is sufficient to determine precisely the x,y coordinates of all virtual grid points. Each virtual grid point is then referred to by integer coordinates vx,vy, where 0,0 corresponds to the lower-left corner and nx-1,ny-1 corresponds to the upper-right corner. vx and vy are also referred to as the virtual column and row indices respectively.
The base grid is then specified with respect to the virtual grid, with each base grid point coinciding with a virtual grid point. Each base grid point must be given an index, starting from 1, and incrementing row-wise from the leftmost point of the lowest row. Also, each base grid row must be numbered consecutively from the lowest row in the grid, so that row 1 contains grid point 1.
As an example, consider the domain consisting of the two separate squares shown in Figure 1. The left-hand diagram shows the virtual grid and its integer coordinates (i.e., its column and row indices), and the right-hand diagram shows the base grid point indices and the base row indices (in brackets).
Figure 1
Figure 1
Hence the base grid point with index 6 say is in base row 2, virtual column 4, and virtual row 1, i.e., virtual grid integer coordinates 4,1; and the base grid point with index 19 say is in base row 5, virtual column 2, and virtual row 5, i.e., virtual grid integer coordinates 2,5.
The base grid must then be defined in inidom by specifying the number of base grid rows, the number of base grid points, the number of boundaries, the number of boundary points, and the following integer arrays:
Finally, ilbnd contains the types of the boundaries and corners, as follows:
Boundaries:
External corners (90°):
Internal corners (270°):
Figure 2 shows the boundary types of a domain with a hole. Notice the logic behind the labelling of the corners: each one includes the types of the two adjacent boundary edges, in a clockwise fashion (outside the domain).
Figure 2
Figure 2
As an example, consider the domain shown in Figure 3. The left-hand diagram shows the physical domain and the right-hand diagram shows the base and virtual grids. The numbers outside the base grid are the indices of the left and rightmost base grid points, and the numbers inside the base grid are the boundary or corner numbers, indicating the order in which the boundaries are stored in lbnd.
Figure 3
Figure 3
The function nag_pde_2d_gen_order2_checkgrid (d03ry) can be called from inidom to obtain a simple graphical representation of the base grid, and to verify the data that you have specified in inidom.
Subgrids are stored internally using the same data structure, and solution information is communicated to you in pdeiv, pdedef and bndary in arrays according to the grid index on the particular level, e.g., xi and yi contain the x,y coordinates of grid point i, and uij contains the jth solution component uj at grid point i.
The grid data and the solutions at all grid levels are 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 provided. 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. The function nag_pde_2d_gen_order2_rectilinear_extractgrid (d03rz) should be called from monitr to obtain grid information at a particular level.
Further details of the underlying algorithm can be found in Further Comments and in Blom and Verwer (1993) and Blom et al. (1996) and the references therein.

References

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
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

Parameters

Compulsory Input Parameters

1:     ts – double scalar
The initial value of the independent variable t.
Constraint: ts<tout.
2:     tout – double scalar
The final value of t to which the integration is to be carried out.
3:     dt3 – double array
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.
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 .
4:     tols – double scalar
The space tolerance used in the grid refinement strategy (σ in equation (4)). See Refinement Strategy.
Constraint: tols>0.0.
5:     tolt – double scalar
The time tolerance used to determine the time step size (τ in equation (7)). See Time Integration .
Constraint: tolt>0.0.
6:     inidom – function handle or string containing name of m-file
inidom must specify the base grid in terms of the data structure described in Description. inidom is not referenced if, on entry, ind=1. nag_pde_2d_gen_order2_checkgrid (d03ry) can be called from inidom to obtain a simple graphical representation of the base grid, and to verify the data that you have specified in inidom. nag_pde_2d_gen_order2_rectilinear (d03rb) also checks the validity of the data, but you are strongly advised to call nag_pde_2d_gen_order2_checkgrid (d03ry) to ensure that the base grid is exactly as required.
Note: the boundaries of the base grid should consist of as many points as are necessary to employ second-order space discretization, i.e., a boundary enclosing the internal part of the domain must include at least 3 grid points including the corners. If Neumann boundary conditions are to be applied the minimum is 4.
[xmin, xmax, ymin, ymax, nx, ny, npts, nrows, nbnds, nbpts, lrow, irow, icol, llbnd, ilbnd, lbnd, ierr] = inidom(maxpts, ierr)

Input Parameters

1:     maxpts int64int32nag_int scalar
The maximum number of base grid points allowed by the available workspace.
2:     ierr int64int32nag_int scalar
Will be initialized by nag_pde_2d_gen_order2_rectilinear (d03rb) to some value prior to internal calls to inidom.

Output Parameters

1:     xmin – double scalar
2:     xmax – double scalar
The extents of the virtual grid in the x-direction, i.e., the x coordinates of the left and right boundaries respectively.
3:     ymin – double scalar
4:     ymax – double scalar
The extents of the virtual grid in the y-direction, i.e., the y coordinates of the left and right boundaries respectively.
5:     nx int64int32nag_int scalar
6:     ny int64int32nag_int scalar
The number of virtual grid points in the x- and y-direction respectively (including the boundary points).
7:     npts int64int32nag_int scalar
The total number of points in the base grid. If the required number of points is greater than maxpts then inidom must be exited immediately with ierr set to -1 to avoid overwriting memory.
8:     nrows int64int32nag_int scalar
The total number of rows of the virtual grid that contain base grid points. This is the maximum base row index.
9:     nbnds int64int32nag_int scalar
The total number of physical boundaries and corners in the base grid.
10:   nbpts int64int32nag_int scalar
The total number of boundary points in the base grid.
11:   lrow: int64int32nag_int array
lrowi, for i=1,2,,nrows, must contain the base grid index of the first grid point in base grid row i.
12:   irow: int64int32nag_int array
irowi, for i=1,2,,nrows, must contain the virtual row number vy that corresponds to base grid row i.
13:   icol: int64int32nag_int array
icoli, for i=1,2,,npts, must contain the virtual column number vx that contains base grid point i.
14:   llbnd: int64int32nag_int array
llbndi, for i=1,2,,nbnds, must contain the element of lbnd corresponding to the start of the ith boundary or corner.
Note: the order of the boundaries and corners in llbnd must be first all the boundaries and then all the corners. The end points of a boundary (i.e., the adjacent corner points) must not be included in the list of points on that boundary. Also, if a corner is shared by two pairs of physical boundaries then it has two types and must therefore be treated as two corners.
15:   ilbnd: int64int32nag_int array
ilbndi, for i=1,2,,nbnds, must contain the type of the ith boundary (or corner), as given in Description.
16:   lbnd: int64int32nag_int array
lbndi, for i=1,2,,nbpts, must contain the grid index of the ith boundary point. The order of the boundaries is as specified in llbnd, but within this restriction the order of the points in lbnd is arbitrary.
17:   ierr int64int32nag_int scalar
If the required number of grid points is larger than maxpts, ierr must be set to -1 to force a termination of the integration and an immediate return to the calling program with ifail=3. Otherwise, ierr should remain unchanged.
7:     pdedef – function handle or string containing name of m-file
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.
[res] = pdedef(npts, npde, t, x, y, u, ut, ux, uy, uxx, uxy, uyy)

Input Parameters

1:     npts int64int32nag_int scalar
The number of grid points in the current grid.
2:     npde int64int32nag_int scalar
The number of PDEs in the system.
3:     t – double scalar
The current value of the independent variable t.
4:     xnpts – double array
xi contains the x coordinate of the ith grid point, for i=1,2,,npts.
5:     ynpts – double array
yi contains the y coordinate of the ith grid point, for i=1,2,,npts.
6:     unptsnpde – double array
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 – double array
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 – double array
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 – double array
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 – double array
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 – double array
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 – double array
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.

Output Parameters

1:     resnptsnpde – double array
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.
8:     bndary – function handle or string containing name of m-file
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 function.
[res] = bndary(npts, npde, t, x, y, u, ut, ux, uy, nbnds, nbpts, llbnd, ilbnd, lbnd, res)

Input Parameters

1:     npts int64int32nag_int scalar
The number of grid points in the current grid.
2:     npde int64int32nag_int scalar
The number of PDEs in the system.
3:     t – double scalar
The current value of the independent variable t.
4:     xnpts – double array
xi contains the x coordinate of the ith grid point, for i=1,2,,npts.
5:     ynpts – double array
yi contains the y coordinate of the ith grid point, for i=1,2,,npts.
6:     unptsnpde – double array
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 – double array
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 – double array
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 – double array
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:   nbnds int64int32nag_int scalar
The total number of physical boundaries and corners in the grid.
11:   nbpts int64int32nag_int scalar
The total number of boundary points in the grid.
12:   llbndnbnds int64int32nag_int array
llbndi, for i=1,2,,nbnds, contains the element of lbnd corresponding to the start of the ith boundary (or corner).
13:   ilbndnbnds int64int32nag_int array
ilbndi, for i=1,2,,nbnds, contains the type of the ith boundary, as given in Description.
14:   lbndnbpts int64int32nag_int array
lbndi, contains the grid index of the ith boundary point, where the order of the boundaries is as specified in llbnd. Hence the ith boundary point has coordinates xlbndi and ylbndi, and the corresponding solution values are ulbndij, for i=1,2,,nbpts and j=1,2,,npde.
15:   resnptsnpde – double array
Contains function values returned by pdedef.

Output Parameters

1:     resnptsnpde – double array
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, i.e., points not included in lbnd, must not be altered.
9:     pdeiv – function handle or string containing name of m-file
pdeiv must specify the initial values of the PDE components u at all points in the base grid. pdeiv is not referenced if, on entry, ind=1.
[u] = pdeiv(npts, npde, t, x, y)

Input Parameters

1:     npts int64int32nag_int scalar
The number of grid points in the base grid.
2:     npde int64int32nag_int scalar
The number of PDEs in the system.
3:     t – double scalar
The (initial) value of the independent variable t.
4:     xnpts – double array
xi contains the x coordinate of the ith grid point, for i=1,2,,npts.
5:     ynpts – double array
yi contains the y coordinate of the ith grid point, for i=1,2,,npts.

Output Parameters

1:     unptsnpde – double array
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.
10:   monitr – function handle or string containing name of m-file
monitr is called by nag_pde_2d_gen_order2_rectilinear (d03rb) 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. nag_pde_2d_gen_order2_rectilinear_extractgrid (d03rz) should be called from monitr in order to extract the number of points and their x,y coordinates on a particular grid.
monitr can also be used to force an immediate tidy termination of the solution process and return to the calling program.
[ierr] = monitr(npde, t, dt, dtnew, tlast, nlev, xmin, ymin, dxb, dyb, lgrid, istruc, lsol, sol, ierr)

Input Parameters

1:     npde int64int32nag_int scalar
The number of PDEs in the system.
2:     t – double scalar
The current value of the independent variable t, i.e., the time at the end of the integration step just completed.
3:     dt – double scalar
The current time step size Δt, i.e., the time step size used for the integration step just completed.
4:     dtnew – double scalar
The time step size that will be used for the next time step.
5:     tlast – logical scalar
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 int64int32nag_int scalar
The number of grid levels used at time t.
7:     xmin – double scalar
8:     ymin – double scalar
The x,y coordinates of the lower-left corner of the virtual grid.
9:     dxb – double scalar
10:   dyb – double scalar
The sizes of the base grid spacing in the x- and y-direction respectively.
11:   lgrid: int64int32nag_int array
Contains pointers to the start of the grid structures in istruc, and must be passed unchanged to nag_pde_2d_gen_order2_rectilinear_extractgrid (d03rz) in order to extract the grid information.
12:   istruc: int64int32nag_int array
Contains the grid structures for each grid level and must be passed unchanged to nag_pde_2d_gen_order2_rectilinear_extractgrid (d03rz) in order to extract the grid information.
13:   lsolnlev int64int32nag_int array
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.)
14:   sollenrwk-6×npde+1 – double array
Contains the solution u at time t for each grid level l in turn, positioned according to lsol. More precisely
uij= sol lsoll+ j-1 × nl+i  
represents the jth component of the solution at the ith grid point in the lth level, for i=1,2,,nl, j=1,2,,npde and l=1,2,,nlev, where nl is the number of grid points at level l (obtainable by a call to nag_pde_2d_gen_order2_rectilinear_extractgrid (d03rz)).
15:   ierr int64int32nag_int scalar
Will be initialized by nag_pde_2d_gen_order2_rectilinear (d03rb) to some value prior to internal calls to ierr.

Output Parameters

1:     ierr int64int32nag_int scalar
Should be set to 1 to force a termination of the integration and an immediate return to the calling program with ifail=4. ierr should remain unchanged otherwise.
11:   opti4 int64int32nag_int array
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.
12:   optr3npde – double array
May be used to specify the optional vectors umax, ws and wt in the space and time monitors (see Further Comments).
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.
13:   rwklenrwk – double array
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 function returns with ifail=3 and an estimated required size is printed on the current error message unit (see nag_file_set_unit_error (x04aa)).
Note:  the size of lenrwk cannot be checked upon initial entry to nag_pde_2d_gen_order2_rectilinear (d03rb) since the number of grid points on the base grid is not known.
14:   iwkleniwk int64int32nag_int array
If ind=0, iwk need not be set. Otherwise iwk must remain unchanged from a previous call to nag_pde_2d_gen_order2_rectilinear (d03rb).
15:   lenlwk int64int32nag_int scalar
The dimension of the array lwk.
The required value of lenlwk cannot be determined exactly in advance, 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 function returns with ifail=3 and an estimated required size is printed on the current error message unit (see nag_file_set_unit_error (x04aa)).
Note:  the size of lenlwk cannot be checked upon initial entry to nag_pde_2d_gen_order2_rectilinear (d03rb) since the number of grid points on the base grid is not known.
16:   itrace int64int32nag_int scalar
The level of trace information required from nag_pde_2d_gen_order2_rectilinear (d03rb). 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 nag_file_set_unit_advisory (x04ab)). This output contains details of the time integration, the nonlinear iteration and the linear solver.
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. Setting itrace=1 allows you to monitor the progress of the integration without possibly excessive information.
17:   ind int64int32nag_int scalar
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 function. In this case, only the following parameters may be reset between calls to nag_pde_2d_gen_order2_rectilinear (d03rb): tout, dt, tols, tolt, opti, optr, itrace and ifail. pdedef is assumed to be serial.
ind=10
Equivalent to ind=0. This option is included only for compatibility with other NAG library products.
ind=11
Equivalent to ind=1. This option is included only for compatibility with other NAG library products.
Constraint: 0ind1 or 10ind11.

Optional Input Parameters

1:     npde int64int32nag_int scalar
Default: the dimension of the array optr.
The number of PDEs in the system.
Constraint: npde1.
2:     lenrwk int64int32nag_int scalar
Default: the dimension of the array rwk.
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 function returns with ifail=3 and an estimated required size is printed on the current error message unit (see nag_file_set_unit_error (x04aa)).
Note:  the size of lenrwk cannot be checked upon initial entry to nag_pde_2d_gen_order2_rectilinear (d03rb) since the number of grid points on the base grid is not known.
3:     leniwk int64int32nag_int scalar
Default: the dimension of the array iwk.
The dimension of the array iwk.
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 function returns with ifail=3 and an estimated required size is printed on the current error message unit (see nag_file_set_unit_error (x04aa)).
Note:  the size of leniwk cannot be checked upon initial entry to nag_pde_2d_gen_order2_rectilinear (d03rb) since the number of grid points on the base grid is not known.

Output Parameters

1:     ts – double scalar
The value of t which has been reached. Normally ts=tout.
2:     dt3 – double array
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.
3:     rwklenrwk – double array
Communication array, used to store information between calls to nag_pde_2d_gen_order2_rectilinear (d03rb).
4:     iwkleniwk int64int32nag_int array
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 nag_pde_2d_gen_order2_rectilinear (d03rb). If the specified maximum number of Newton or linear solver iterations is exceeded at any stage, then the maximums above are set to the specified maximum plus one.
5:     ind int64int32nag_int scalar
ind=1, if ind on input was 0 or 1, or ind=11, if ind on input was 10 or 11.
6:     ifail int64int32nag_int scalar
ifail=0 unless the function detects an error (see Error Indicators and Warnings).

Error Indicators and Warnings

Errors or warnings detected by the function:

Cases prefixed with W are classified as warnings and do not generate an error of type NAG:error_n. See nag_issue_warnings.

   ifail=1
On entry,npde<1,
ortoutts,
ortout is too close to ts,
orind=0 and dt1<0.0,
or dti<0.0, for i=2​ or ​3,
ordt2>dt3,
orind=0 and 0.0<dt1<10×machine precision×maxts,tout,
orind=0 and dt1>tout-ts,
orind=0 and dt1<dt2​ or ​dt1>dt3,
ortols or tolt0.0,
oropti1<0,
or opti1>0 and optij<0, for j=2, 3 or 4,
oroptr1j0.0, for some j=1,2,,npde,
oroptr2j<0.0, for some j=1,2,,npde,
oroptr3j<0.0, for some j=1,2,,npde,
orind0 or 1,
orind=1 on initial entry to nag_pde_2d_gen_order2_rectilinear (d03rb).
   ifail=2
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 is too small for the required number of grid points. At the initial time step this error may result because you set ierr to -1 in inidom or the internal check on the number of grid points following the call to inidom. An estimate of the required sizes for the current stage is output, but more space may be required at a later stage.
W  ifail=4
ierr was set to 1 in monitr, forcing control to be passed back to calling program. Integration was successful as far as t=ts.
W  ifail=5
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=6
One or more of the output arguments of inidom was incorrectly specified, i.e.,
xminxmax,
orxmax too close to xmin,
oryminymax,
orymax too close to ymin,
ornx or ny<4,
ornrows<4,
ornrows>ny,
ornpts>nx×ny,
ornbnds<8,
ornbpts<12,
ornbptsnpts,
orlrowi<1 or lrowi>npts, for some i=1,2,,nrows,
orlrowilrowi-1, for some i=2,3,,nrows,
orirowi<0 or irowi>ny, for some i=1,2,,nrows,
orirowiirowi-1, for some i=2,3,,nrows,
oricoli<0 or icoli>nx, for some i=1,2,,npts,
orllbndi<1 or llbndi>nbpts, for some i=1,2,,nbnds,
orllbndillbndi-1, for some i=2,3,,nbnds,
orilbndi1, 2, 3, 4, 12, 23, 34, 41, 21, 32, 43 or 14, for some i=1,2,,nbnds,
orlbndi<1 or lbndi>npts, for some i=1,2,,nbpts.
   ifail=-99
An unexpected error has been triggered by this routine. Please contact NAG.
   ifail=-399
Your licence key may have expired or may not have been installed correctly.
   ifail=-999
Dynamic memory allocation failed.

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 Further Comments, 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).

Further Comments

Algorithm Outline

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

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 uj 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,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 values 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 Algorithm Outline.

Time Integration

The time integration is controlled using a time monitor calculated at each level l up to the maximum level used, given by
μ l t = 1 N j = 1 npde w j t i = 1 ngpts l Δ t α i j u t x i , y i ,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
α i j = τ u j max 100 + u x i , y i ,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)

Example

This example is taken from Blom and Verwer (1993) and is the two-dimensional Burgers' system
u t =-u u x -v u y +ε 2 u x2 + 2 u y2 ,  
v t =-u v x -v v y +ε 2v x2 + 2v y2 ,  
with ε=10-3 on the domain given in Figure 3. Dirichlet boundary conditions are used on all boundaries using the exact solution
u=34-141+exp-4x+4y-t/32ε ,  
v = 3 4 + 1 4 1 + exp - 4 x + 4 y - t / 32 ε .  
The solution contains a wave front at y=x+0.25t which propagates in a direction perpendicular to the front with speed 2/8.
function d03rb_example


fprintf('d03rb example results\n\n');

global iout xsol ysol usol vsol inds;

ts    = 0;
twant = [0.25; 1];
dt    = [0.001; 1e-07; 0];
tols  = 0.1;
tolt  = 0.05;
opti  = int64([5;0;0;0]);
optr = [1, 1; 1, 1; 1, 1];
rwk  = zeros(426000, 1);
iwk  = zeros(117037, 1, 'int64');
lenlwk = int64(6000);
itrace = int64(0);
ind  = int64(0);
inds  = int64(0);

for iout = 1:2
  tout = twant(iout);
  [ts, dt, rwk, iwk, ind, ifail] = ...
    d03rb(...
          ts, tout, dt, tols, tolt, @inidom, @pdedef, ...
          @bndary, @pdeiv, @monitr, opti, optr, ...
          rwk, iwk, lenlwk, itrace, ind);
  fprintf('\nStatistics\n');
  fprintf('Time = %8.4f\n', ts);
  fprintf('Total number of accepted timesteps = %d\n', iwk(1));
  fprintf('Total number of rejected timesteps = %d\n', iwk(2));
  fprintf('\n             T o t a l  n u m b e r  o f    \n');
  fprintf('          Residual  Jacobian    Newton  Lin sys\n');
  fprintf('             evals     evals     iters    iters\n');
  fprintf(' At level \n');
  maxlev = opti(1);
  for j = 1:maxlev
    if (iwk(j+2) ~= 0)
      fprintf('%6d %10d %10d %10d %10d\n', j, iwk(j+2), iwk(j+2+maxlev), ...
                                 iwk(j+2+2*maxlev), iwk(j+2+3*maxlev) );
    end
  end
  fprintf('\n             M a x i m u m  n u m b e r  o f\n');
  fprintf('              Newton iters    Lin sys iters \n');
  fprintf(' At level \n');
  for j = 1:maxlev
    if (iwk(j+2) ~= 0)
      fprintf('%6d%14d%14d\n', j, iwk(j+2+4*maxlev), iwk(j+2+5*maxlev) );
    end
  end
end

% Plot solutions u and v at t=1.0
fig1 =  figure;
scatter3(xsol,ysol,usol,15,vsol,'o','fill');
title('2D Burgers'' equation at t=1 on disjoint domain: u');
xlabel('x');
ylabel('y');
zlabel('u(x,y;t=1)');
fig2 =  figure;
scatter3(xsol,ysol,vsol,15,vsol,'o','fill');
view(-54,32);
xlabel('x');
ylabel('y');
zlabel('v(x,y;t=1)');
title('2D Burgers'' equation at t=1 on disjoint domain: v');


function [res] = bndary(npts, npde, t, x, y, u, ut, ux, uy, nbnds, ...
                        nbpts, llbnd, ilbnd, lbnd, res)

  epsilon = 1e-3;

  for k = llbnd(1):nbpts
    i = lbnd(k);
    a = (-4*x(i)+4*y(i)-t)/(32*epsilon);
    if (a <= 0)
      res(i,1) = u(i,1) - (0.75-0.25/(1+exp(a)));
      res(i,2) = u(i,2) - (0.75+0.25/(1+exp(a)));
    else
      res(i,1) = u(i,1) - (0.75-0.25*exp(-a)/(exp(-a)+1));
      res(i,2) = u(i,2) - (0.75+0.25*exp(-a)/(exp(-a)+1));
    end
  end
  

function [xmin, xmax, ymin, ymax, nx, ny, npts, nrows, nbnds, nbpts, ...
          lrow, irow, icol, llbnd, ilbnd, lbnd, ierr] = inidom(maxpts, ierr)

  nrows = int64(11);
  npts  = int64(105);
  nbnds = int64(28);
  nbpts = int64(72);

  lrow  = zeros(nrows, 1, 'int64');
  irow  = zeros(nrows, 1, 'int64');
  icol  = zeros(npts, 1, 'int64');
  llbnd = zeros(nbnds, 1, 'int64');
  ilbnd = zeros(nbnds, 1, 'int64');
  lbnd  = zeros(nbpts, 1, 'int64');



icold = int64([0,1,2,0,1,2,3,4,5,6,7,8,9,10,0,1,2,3,4,5,6,7,8,9, ...
                  10,0,1,2,3,4,5,6,7,8,9,10,0,1,2,3,4,5,8,9,10,0,1,2,3,4,5, ...
                  6,7,8,9,10,0,1,2,3,4,5,6,7,8,9,10,0,1,2,3,4,5,6,7,8,9,10, ...
                  0,1,2,3,4,5,6,7,8,0,1,2,3,4,5,6,7,8,0,1,2,3,4,5,6,7,8]);

ilbndd = int64([1,2,3,4,1,4,1,2,3,4,3,4,1,2,12,23,34,41,14,41,12, ...
                   23,34,41,43,14,21,32]);

irowd  = int64([0,1,2,3,4,5,6,7,8,9,10]);

lbndd  = int64([2,4,15,26,37,46,57,68,79,88,98,99,100,101,102,103,104,96,...
                   86,85,84,83,82,70,59,48,39,28,17,6,8,9,10,11,12,13,18,...
                   29,40,49,60,72,73,74,75,76,77,67,56,45,36,25,33,32,42, ...
                   52,53,43,1,97,105,87,81,3,7,71,78,14,31,51,54,34]);

llbndd = int64([1,2,11,18,19,24,31,37,42,48,53,55,56,58,59,60,61,62, ...
                   63,64,65,66,67,68,69,70,71,72]);

lrowd  = int64([1,4,15,26,37,46,57,68,79,88,97]);

nx = int64(11);
ny = int64(11);

% check maxpts against rough estimate of npts
  if (maxpts < nx*ny)
     ierr = -1;
     return;
  end

  xmin = 0;
  ymin = 0;
  xmax = 1;
  ymax = 1;

lrow(1:nrows)  = lrowd(1:nrows);
irow(1:nrows)  = irowd(1:nrows);
llbnd(1:nbnds) = llbndd(1:nbnds);
ilbnd(1:nbnds) = ilbndd(1:nbnds);
lbnd(1:nbpts)  = lbndd(1:nbpts);
icol(1:npts)   = icold(1:npts);

function [res] = pdedef(npts, npde, t, x, y, u, ut, ux, uy, uxx, uxy, uyy)
  res = zeros(npts, npde);

  epsilon = 1e-3;
  uxxyy = epsilon*(uxx(1:npts,1:2)+uyy(1:npts,1:2));
  uuxuy(1:npts,1) = -u(1:npts,1).*ux(1:npts,1)-u(1:npts,2).*uy(1:npts,1);
  uuxuy(1:npts,2) = -u(1:npts,1).*ux(1:npts,2)-u(1:npts,2).*uy(1:npts,2);
  res(1:npts,1:2) = ut(1:npts,1:2)-(uuxuy+uxxyy);

function [u] = pdeiv(npts, npde, t, x, y)
  u = zeros(npts, npde);

  epsilon = 1e-3;

  for i = 1:npts
    a = (-4*x(i)+4*y(i)-t)/(32*epsilon);
    if (a <= 0)
       u(i,1) = 0.75 - 0.25/(1+exp(a));
       u(i,2) = 0.75 + 0.25/(1+exp(a));
    else
       u(i,1) = 0.75 - 0.25*exp(-a)/(exp(-a)+1);
       u(i,2) = 0.75 + 0.25*exp(-a)/(exp(-a)+1);
    end
  end


function [ierr] = ...
    monitr(npde, t, dt, dtnew, tlast, nlev, xmin, ymin, dxb, dyb, lgrid, ...
                 istruc, lsol, sol, ierr)

  lenxy=int64(2500);

  global iout xsol ysol usol vsol inds;

  if (tlast && iout == 2)
    for level = int64(1:nlev)
      ipsol = lsol(level);

      % Get grid information
      [npts, x, y, ifail] = ...
         d03rz(...
               level, xmin, ymin, dxb, dyb, lgrid, istruc, lenxy);
      % Get exact solution
      [uex] = pdeiv(npts, npde, t, x, y);

      xsol(inds+1:inds+npts) = x(1:npts);
      ysol(inds+1:inds+npts) = y(1:npts);
      usol(inds+1:inds+npts) = sol(ipsol+1:ipsol+npts);
      vsol(inds+1:inds+npts) = sol(ipsol+npts+1:ipsol+npts+npts);
      inds = inds + npts;
    end
  end
d03rb example results


Statistics
Time =   0.2500
Total number of accepted timesteps = 14
Total number of rejected timesteps = 0

             T o t a l  n u m b e r  o f    
          Residual  Jacobian    Newton  Lin sys
             evals     evals     iters    iters
 At level 
     1        196         14         28         14
     2        196         14         28         22
     3        196         14         28         25
     4        196         14         28         31
     5        141         10         21         29

             M a x i m u m  n u m b e r  o f
              Newton iters    Lin sys iters 
 At level 
     1             2             1
     2             2             1
     3             2             1
     4             2             2
     5             3             2

Statistics
Time =   1.0000
Total number of accepted timesteps = 45
Total number of rejected timesteps = 0

             T o t a l  n u m b e r  o f    
          Residual  Jacobian    Newton  Lin sys
             evals     evals     iters    iters
 At level 
     1        630         45         90         45
     2        630         45         90         78
     3        630         45         90         87
     4        630         45         90        124
     5        575         41         83        122

             M a x i m u m  n u m b e r  o f
              Newton iters    Lin sys iters 
 At level 
     1             2             1
     2             2             1
     3             2             1
     4             2             2
     5             3             2
d03rb_fig1.png
d03rb_fig2.png

PDF version (NAG web site, 64-bit version, 64-bit version)
Chapter Contents
Chapter Introduction
NAG Toolbox

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