NAG FL Interface
d03rbf (dim2_gen_order2_rectilinear)
1
Purpose
d03rbf 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.
d03rbf originates from the VLUGR2 package (see
Blom and Verwer (1993) and
Blom et al. (1996)).
2
Specification
Fortran Interface
Subroutine d03rbf ( |
npde, ts, tout, dt, tols, tolt, inidom, pdedef, bndary, pdeiv, monitr, opti, optr, rwk, lenrwk, iwk, leniwk, lwk, lenlwk, itrace, ind, ifail) |
Integer, Intent (In) |
:: |
npde, opti(4), lenrwk, leniwk, lenlwk, itrace |
Integer, Intent (Inout) |
:: |
iwk(leniwk), ind, ifail |
Real (Kind=nag_wp), Intent (In) |
:: |
tout, tols, tolt, optr(3,npde) |
Real (Kind=nag_wp), Intent (Inout) |
:: |
ts, dt(3), rwk(lenrwk) |
Logical, Intent (Out) |
:: |
lwk(lenlwk) |
External |
:: |
inidom, pdedef, bndary, pdeiv, monitr |
|
C Header Interface
#include <nag.h>
void |
d03rbf_ (const Integer *npde, double *ts, const double *tout, double dt[], const double *tols, const double *tolt, void (NAG_CALL *inidom)(const Integer *maxpts, double *xmin, double *xmax, double *ymin, double *ymax, Integer *nx, Integer *ny, Integer *npts, Integer *nrows, Integer *nbnds, Integer *nbpts, Integer lrow[], Integer irow[], Integer icol[], Integer llbnd[], Integer ilbnd[], Integer lbnd[], Integer *ierr), 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 *nbnds, const Integer *nbpts, const Integer llbnd[], const Integer ilbnd[], 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 double *xmin, const double *ymin, const double *dxb, const double *dyb, const Integer lgrid[], const Integer istruc[], 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) |
|
C++ Header Interface
#include <nag.h> extern "C" {
void |
d03rbf_ (const Integer &npde, double &ts, const double &tout, double dt[], const double &tols, const double &tolt, void (NAG_CALL *inidom)(const Integer &maxpts, double &xmin, double &xmax, double &ymin, double &ymax, Integer &nx, Integer &ny, Integer &npts, Integer &nrows, Integer &nbnds, Integer &nbpts, Integer lrow[], Integer irow[], Integer icol[], Integer llbnd[], Integer ilbnd[], Integer lbnd[], Integer &ierr), 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 &nbnds, const Integer &nbpts, const Integer llbnd[], const Integer ilbnd[], 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 double &xmin, const double &ymin, const double &dxb, const double &dyb, const Integer lgrid[], const Integer istruc[], 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 d03rbf or nagf_pde_dim2_gen_order2_rectilinear.
3
Description
d03rbf integrates the system of PDEs:
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
d03raf is used.
The vector
is the set of solution values
and
denotes partial differentiation with respect to
, and similarly for
, etc.
The functions
must be supplied by you in
pdedef. Similarly the initial values of the functions
for
must be specified at
in
pdeiv.
Note that whilst complete generality is offered by the master equations
(1),
d03rbf is not appropriate for all PDEs. In particular, hyperbolic systems should not be solved using this routine. Also, at least one component of
must appear in the system of PDEs.
The boundary conditions must be supplied by you in
bndary in the form
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
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
. In this way a solution is obtained at
.
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
. 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 space is given by the coordinates of its boundaries. The number of points and in the and directions must also be given, corresponding to the number of columns and rows respectively. This is sufficient to determine precisely the coordinates of all virtual grid points. Each virtual grid point is then referred to by integer coordinates , where corresponds to the lower-left corner and corresponds to the upper-right corner. and 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 , 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 contains grid point .
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).
Hence the base grid point with index say is in base row , virtual column 4, and virtual row , i.e., virtual grid integer coordinates ; and the base grid point with index say is in base row , virtual column 2, and virtual row , i.e., virtual grid integer coordinates .
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:
- lrow contains the base grid indices of the starting points of the base grid rows;
- irow contains the virtual row numbers of the base grid rows;
- icol contains the virtual column numbers of the base grid points;
- lbnd contains the grid indices of the boundary edges (without corners) and corner points;
- llbnd contains the starting elements of the boundaries and corners in lbnd.
Finally,
ilbnd contains the types of the boundaries and corners, as follows:
Boundaries:
- 1 – lower boundary
- 2 – left boundary
- 3 – upper boundary
- 4 – right boundary
External corners (
):
- 12 – lower-left corner
- 23 – upper-left corner
- 34 – upper-right corner
- 41 – lower-right corner
Internal corners (
):
- 21 – lower-left corner
- 32 – upper-left corner
- 43 – upper-right corner
- 14 – lower-right corner
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).
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.
For this example we have
nrows = 11
npts = 105
nbnds = 28
nbpts = 72
lrow = (1,4,15,26,37,46,57,68,79,88,97)
irow = (0,1,2,3,4,5,6,7,8,9,10)
icol = (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)
lbnd = (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)
llbnd = (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)
ilbnd = (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)
This particular domain is used in the example in
Section 10, and data statements are used to define the above arrays in that example program. For less complicated domains it is simpler to assign the values of the arrays in do-loops. This also allows flexibility in the number of base grid points.
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.,
and
contain the
coordinates of grid point
, and
contains the
th solution component
at grid point
.
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 routine
d03rzf should be called from
monitr to obtain grid information at a particular level.
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
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
5
Arguments
-
1:
– Integer
Input
-
On entry: the number of PDEs in the system.
Constraint:
.
-
2:
– Real (Kind=nag_wp)
Input/Output
-
On entry: the initial value of the independent variable .
On exit: the value of which has been reached. Normally .
Constraint:
.
-
3:
– Real (Kind=nag_wp)
Input
-
On entry: the final value of to which the integration is to be carried out.
-
4:
– Real (Kind=nag_wp) array
Input/Output
-
On entry: the initial, minimum and maximum time step sizes respectively.
- Specifies the initial time step size to be used on the first entry, i.e., when . If then the default value is used. On subsequent entries (), the value of is not referenced.
- Specifies the minimum time step size to be attempted by the integrator. If the default value is used.
- Specifies the maximum time step size to be attempted by the integrator. If the default value is used.
On exit: contains the time step size for the next time step. and are unchanged or set to their default values if zero on entry.
Constraints:
- if , ;
- if and , and , where the values of and will have been reset to their default values if zero on entry;
- .
-
5:
– Real (Kind=nag_wp)
Input
-
On entry: the space tolerance used in the grid refinement strategy (
in equation
(4)). See
Section 9.2.
Constraint:
.
-
6:
– 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:
.
-
7:
– Subroutine, supplied by the user.
External Procedure
-
inidom must specify the base grid in terms of the data structure described in
Section 3.
inidom is not referenced if, on entry,
.
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 grid points including the corners. If Neumann boundary conditions are to be applied the minimum is .
The specification of
inidom is:
Fortran Interface
Subroutine inidom ( |
maxpts, xmin, xmax, ymin, ymax, nx, ny, npts, nrows, nbnds, nbpts, lrow, irow, icol, llbnd, ilbnd, lbnd, ierr) |
Integer, Intent (In) |
:: |
maxpts |
Integer, Intent (Inout) |
:: |
lrow(*), irow(*), icol(*), llbnd(*), ilbnd(*), lbnd(*), ierr |
Integer, Intent (Out) |
:: |
nx, ny, npts, nrows, nbnds, nbpts |
Real (Kind=nag_wp), Intent (Out) |
:: |
xmin, xmax, ymin, ymax |
|
C Header Interface
void |
inidom_ (const Integer *maxpts, double *xmin, double *xmax, double *ymin, double *ymax, Integer *nx, Integer *ny, Integer *npts, Integer *nrows, Integer *nbnds, Integer *nbpts, Integer lrow[], Integer irow[], Integer icol[], Integer llbnd[], Integer ilbnd[], Integer lbnd[], Integer *ierr) |
|
C++ Header Interface
#include <nag.h> extern "C" {
void |
inidom_ (const Integer &maxpts, double &xmin, double &xmax, double &ymin, double &ymax, Integer &nx, Integer &ny, Integer &npts, Integer &nrows, Integer &nbnds, Integer &nbpts, Integer lrow[], Integer irow[], Integer icol[], Integer llbnd[], Integer ilbnd[], Integer lbnd[], Integer &ierr) |
}
|
-
1:
– Integer
Input
-
On entry: the maximum number of base grid points allowed by the available workspace.
-
2:
– Real (Kind=nag_wp)
Output
-
3:
– Real (Kind=nag_wp)
Output
-
On exit: the extents of the virtual grid in the -direction, i.e., the coordinates of the left and right boundaries respectively.
Constraint:
and
xmax must be sufficiently distinguishable from
xmin for the precision of the machine being used.
-
4:
– Real (Kind=nag_wp)
Output
-
5:
– Real (Kind=nag_wp)
Output
-
On exit: the extents of the virtual grid in the -direction, i.e., the coordinates of the left and right boundaries respectively.
Constraint:
and
ymax must be sufficiently distinguishable from
ymin for the precision of the machine being used.
-
6:
– Integer
Output
-
7:
– Integer
Output
-
On exit: the number of virtual grid points in the - and -direction respectively (including the boundary points).
-
8:
– Integer
Output
-
On exit: 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
to avoid overwriting memory.
Constraint:
and if on exit, .
-
9:
– Integer
Output
-
On exit: the total number of rows of the virtual grid that contain base grid points. This is the maximum base row index.
Constraint:
.
-
10:
– Integer
Output
-
On exit: the total number of physical boundaries and corners in the base grid.
Constraint:
.
-
11:
– Integer
Output
-
On exit: the total number of boundary points in the base grid.
Constraint:
.
-
12:
– Integer array
Output
-
On exit: , for , must contain the base grid index of the first grid point in base grid row .
Constraints:
- , for ;
- , for .
-
13:
– Integer array
Output
-
On exit: , for , must contain the virtual row number that corresponds to base grid row .
Constraints:
- , for ;
- , for .
-
14:
– Integer array
Output
-
On exit: , for , must contain the virtual column number that contains base grid point .
Constraint:
, for .
-
15:
– Integer array
Output
-
On exit:
, for
, must contain the element of
lbnd corresponding to the start of the
th 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.
Constraints:
- , for ;
- , for .
-
16:
– Integer array
Output
-
On exit:
, for
, must contain the type of the
th boundary (or corner), as given in
Section 3.
Constraint:
must be equal to one of the following: , , , , , , , , , , or , for .
-
17:
– Integer array
Output
-
On exit:
, for
, must contain the grid index of the
th 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.
Constraint:
, for .
-
18:
– Integer
Input/Output
-
On entry: will be initialized by
d03rbf to some value prior to internal calls to
inidom.
On exit: if the required number of grid points is larger than
maxpts,
ierr must be set to
to force a termination of the integration and an immediate return to the calling program with
. Otherwise,
ierr should remain unchanged.
inidom must either be a module subprogram USEd by, or declared as EXTERNAL in, the (sub)program from which
d03rbf is called. Arguments denoted as
Input must
not be changed by this procedure.
Note: inidom should not return floating-point NaN (Not a Number) or infinity values, since these are not handled by
d03rbf. If your code inadvertently
does return any NaNs or infinities,
d03rbf is likely to produce unexpected results.
-
8:
– Subroutine, supplied by the user.
External Procedure
-
pdedef must evaluate the functions
, for
, 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[]) |
|
C++ Header Interface
#include <nag.h> extern "C" {
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:
– Integer
Input
-
On entry: the number of grid points in the current grid.
-
2:
– Integer
Input
-
On entry: the number of PDEs in the system.
-
3:
– Real (Kind=nag_wp)
Input
-
On entry: the current value of the independent variable .
-
4:
– Real (Kind=nag_wp) array
Input
-
On entry: contains the coordinate of the th grid point, for .
-
5:
– Real (Kind=nag_wp) array
Input
-
On entry: contains the coordinate of the th grid point, for .
-
6:
– Real (Kind=nag_wp) array
Input
-
On entry: contains the value of the th PDE component at the th grid point, for and .
-
7:
– Real (Kind=nag_wp) array
Input
-
On entry: contains the value of for the th PDE component at the th grid point, for and .
-
8:
– Real (Kind=nag_wp) array
Input
-
On entry: contains the value of for the th PDE component at the th grid point, for and .
-
9:
– Real (Kind=nag_wp) array
Input
-
On entry: contains the value of for the th PDE component at the th grid point, for and .
-
10:
– Real (Kind=nag_wp) array
Input
-
On entry: contains the value of for the th PDE component at the th grid point, for and .
-
11:
– Real (Kind=nag_wp) array
Input
-
On entry: contains the value of for the th PDE component at the th grid point, for and .
-
12:
– Real (Kind=nag_wp) array
Input
-
On entry: contains the value of for the th PDE component at the th grid point, for and .
-
13:
– Real (Kind=nag_wp) array
Output
-
On exit: must contain the value of , for , at the
th grid point, for , 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
d03rbf 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
d03rbf. If your code inadvertently
does return any NaNs or infinities,
d03rbf is likely to produce unexpected results.
-
9:
– Subroutine, supplied by the user.
External Procedure
-
bndary must evaluate the functions
, for
, 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, nbnds, nbpts, llbnd, ilbnd, lbnd, res) |
Integer, Intent (In) |
:: |
npts, npde, nbnds, nbpts, llbnd(nbnds), ilbnd(nbnds), 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 *nbnds, const Integer *nbpts, const Integer llbnd[], const Integer ilbnd[], const Integer lbnd[], double res[]) |
|
C++ Header Interface
#include <nag.h> extern "C" {
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 &nbnds, const Integer &nbpts, const Integer llbnd[], const Integer ilbnd[], const Integer lbnd[], double res[]) |
}
|
-
1:
– Integer
Input
-
On entry: the number of grid points in the current grid.
-
2:
– Integer
Input
-
On entry: the number of PDEs in the system.
-
3:
– Real (Kind=nag_wp)
Input
-
On entry: the current value of the independent variable .
-
4:
– Real (Kind=nag_wp) array
Input
-
On entry: contains the coordinate of the th grid point, for .
-
5:
– Real (Kind=nag_wp) array
Input
-
On entry: contains the coordinate of the th grid point, for .
-
6:
– Real (Kind=nag_wp) array
Input
-
On entry: contains the value of the th PDE component at the th grid point, for and .
-
7:
– Real (Kind=nag_wp) array
Input
-
On entry: contains the value of for the th PDE component at the th grid point, for and .
-
8:
– Real (Kind=nag_wp) array
Input
-
On entry: contains the value of for the th PDE component at the th grid point, for and .
-
9:
– Real (Kind=nag_wp) array
Input
-
On entry: contains the value of for the th PDE component at the th grid point, for and .
-
10:
– Integer
Input
-
On entry: the total number of physical boundaries and corners in the grid.
-
11:
– Integer
Input
-
On entry: the total number of boundary points in the grid.
-
12:
– Integer array
Input
-
On entry:
, for
, contains the element of
lbnd corresponding to the start of the
th boundary (or corner).
-
13:
– Integer array
Input
-
On entry:
, for
, contains the type of the
th boundary, as given in
Section 3.
-
14:
– Integer array
Input
-
On entry:
, contains the grid index of the
th boundary point, where the order of the boundaries is as specified in
llbnd. Hence the
th boundary point has coordinates
and
, and the corresponding solution values are
, for
and
.
-
15:
– Real (Kind=nag_wp) array
Input/Output
-
On entry: contains function values returned by
pdedef.
On exit:
must contain the value of
, for
, at the
th boundary point, for
.
Note: elements of
res corresponding to interior points, i.e., points not included in
lbnd, must
not be altered.
bndary must either be a module subprogram USEd by, or declared as EXTERNAL in, the (sub)program from which
d03rbf 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
d03rbf. If your code inadvertently
does return any NaNs or infinities,
d03rbf is likely to produce unexpected results.
-
10:
– Subroutine, supplied by the user.
External Procedure
-
pdeiv must specify the initial values of the PDE components
at all points in the base grid.
pdeiv is not referenced if, on entry,
.
The specification of
pdeiv is:
Fortran Interface
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[]) |
|
C++ Header Interface
#include <nag.h> extern "C" {
void |
pdeiv_ (const Integer &npts, const Integer &npde, const double &t, const double x[], const double y[], double u[]) |
}
|
-
1:
– Integer
Input
-
On entry: the number of grid points in the base grid.
-
2:
– Integer
Input
-
On entry: the number of PDEs in the system.
-
3:
– Real (Kind=nag_wp)
Input
-
On entry: the (initial) value of the independent variable .
-
4:
– Real (Kind=nag_wp) array
Input
-
On entry: contains the coordinate of the th grid point, for .
-
5:
– Real (Kind=nag_wp) array
Input
-
On entry: contains the coordinate of the th grid point, for .
-
6:
– Real (Kind=nag_wp) array
Output
-
On exit: must contain the value of the th PDE component at the th grid point, for and .
pdeiv must either be a module subprogram USEd by, or declared as EXTERNAL in, the (sub)program from which
d03rbf 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
d03rbf. If your code inadvertently
does return any NaNs or infinities,
d03rbf is likely to produce unexpected results.
-
11:
– Subroutine, supplied by the user.
External Procedure
-
monitr is called by
d03rbf 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.
d03rzf should be called from
monitr in order to extract the number of points and their
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.
The specification of
monitr is:
Fortran Interface
Subroutine monitr ( |
npde, t, dt, dtnew, tlast, nlev, xmin, ymin, dxb, dyb, lgrid, istruc, lsol, sol, ierr) |
Integer, Intent (In) |
:: |
npde, nlev, lgrid(*), istruc(*), lsol(nlev) |
Integer, Intent (Inout) |
:: |
ierr |
Real (Kind=nag_wp), Intent (In) |
:: |
t, dt, dtnew, xmin, ymin, dxb, dyb, 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 double *xmin, const double *ymin, const double *dxb, const double *dyb, const Integer lgrid[], const Integer istruc[], const Integer lsol[], const double sol[], Integer *ierr) |
|
C++ Header Interface
#include <nag.h> extern "C" {
void |
monitr_ (const Integer &npde, const double &t, const double &dt, const double &dtnew, const logical &tlast, const Integer &nlev, const double &xmin, const double &ymin, const double &dxb, const double &dyb, const Integer lgrid[], const Integer istruc[], const Integer lsol[], const double sol[], Integer &ierr) |
}
|
-
1:
– Integer
Input
-
On entry: the number of PDEs in the system.
-
2:
– Real (Kind=nag_wp)
Input
-
On entry: the current value of the independent variable , i.e., the time at the end of the integration step just completed.
-
3:
– Real (Kind=nag_wp)
Input
-
On entry: the current time step size , i.e., the time step size used for the integration step just completed.
-
4:
– Real (Kind=nag_wp)
Input
-
On entry: the time step size that will be used for the next time step.
-
5:
– Logical
Input
-
On entry: indicates if intermediate or final time step.
for an intermediate step,
for the last call to
monitr before returning to your program.
-
6:
– Integer
Input
-
On entry: the number of grid levels used at time
t.
-
7:
– Real (Kind=nag_wp)
Input
-
8:
– Real (Kind=nag_wp)
Input
-
On entry: the coordinates of the lower-left corner of the virtual grid.
-
9:
– Real (Kind=nag_wp)
Input
-
10:
– Real (Kind=nag_wp)
Input
-
On entry: the sizes of the base grid spacing in the - and -direction respectively.
-
11:
– Integer array
Input
-
On entry: contains pointers to the start of the grid structures in
istruc, and must be passed unchanged to
d03rzf in order to extract the grid information.
-
12:
– Integer array
Input
-
On entry: contains the grid structures for each grid level and must be passed unchanged to
d03rzf in order to extract the grid information.
-
13:
– Integer array
Input
-
On entry:
contains the pointer to the solution in
sol at grid level
and time
t. (
actually contains the array index immediately preceding the start of the solution in
sol.)
-
14:
– Real (Kind=nag_wp) array
Input
-
On entry: contains the solution
at time
t for each grid level
in turn, positioned according to
lsol.
More precisely
represents the
th component of the solution at the
th grid point in the
th level, for
,
and
, where
is the number of grid points at level
(obtainable by a call to
d03rzf).
-
15:
– Integer
Input/Output
-
On entry: will be initialized by
d03rbf to some value prior to internal calls to
ierr.
On exit: should be set to
to force a termination of the integration and an immediate return to the calling program with
.
ierr should remain unchanged otherwise.
monitr must either be a module subprogram USEd by, or declared as EXTERNAL in, the (sub)program from which
d03rbf is called. Arguments denoted as
Input must
not be changed by this procedure.
-
12:
– Integer array
Input
-
On entry: may be set to control various options available in the integrator.
- All the default options are employed.
- The default value of
, for , can be obtained by setting .
- Specifies the maximum number of grid levels allowed (including the base grid). . The default value is .
- Specifies the maximum number of Jacobian evaluations allowed during each nonlinear equations solution. . The default value is .
- Specifies the maximum number of Newton iterations in each nonlinear equations solution. . The default value is .
- Specifies the maximum number of iterations in each linear equations solution. . The default value is .
Constraint:
and if , , for .
-
13:
– Real (Kind=nag_wp) array
Input
-
On entry: may be used to specify the optional vectors
,
and
in the space and time monitors (see
Section 9).
If an optional vector is not required then all its components should be set to .
, for
, specifies
, the approximate maximum absolute value of the
th component of
, as used in
(4) and
(7).
, for
.
, for
, specifies
, the weighting factors used in the space monitor (see
(4)) to indicate the relative importance of the
th component of
on the space monitor.
, for
.
, for
, specifies
, the weighting factors used in the time monitor (see
(6)) to indicate the relative importance of the
th component of
on the time monitor.
, for
.
Constraints:
- , for ;
- , for and .
-
14:
– Real (Kind=nag_wp) array
Communication Array
-
15:
– Integer
Input
-
On entry: the dimension of the array
rwk as declared in the (sub)program from which
d03rbf is called.
The required value of
lenrwk cannot be determined exactly in advance, but a suggested value is
where
if
and
otherwise, and
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
and an estimated required size is printed on the current error message unit (see
x04aaf).
Note: the size of
lenrwk cannot be checked upon initial entry to
d03rbf since the number of grid points on the base grid is not known.
-
16:
– Integer array
Communication Array
-
On entry: if
,
iwk need not be set. Otherwise
iwk must remain unchanged from a previous call to
d03rbf.
On exit: the following components of the array
iwk concern the efficiency of the integration. Here,
is the maximum number of grid levels allowed (
if
and
otherwise), and
is a grid level taking the values
, where
is the number of levels used.
- Contains the number of steps taken in time.
- Contains the number of rejected time steps.
- Contains the total number of residual evaluations performed (i.e., the number of times pdedef was called) at grid level .
- Contains the total number of Jacobian evaluations performed at grid level .
- Contains the total number of Newton iterations performed at grid level .
- Contains the total number of linear solver iterations performed at grid level .
- Contains the maximum number of Newton iterations performed at any one time step at grid level .
- Contains the maximum number of linear solver iterations performed at any one time step at grid level .
Note: the total and maximum numbers are cumulative over all calls to d03rbf. 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.
-
17:
– Integer
Input
-
On entry: the dimension of the array
iwk as declared in the (sub)program from which
d03rbf is called.
The required value of
leniwk cannot be determined exactly in advance, but a suggested value is
where
is the expected maximum number of grid points at any one level and
if
and
otherwise. If during the execution the supplied value is found to be too small then the routine returns with
and an estimated required size is printed on the current error message unit (see
x04aaf).
Note: the size of
leniwk cannot be checked upon initial entry to
d03rbf since the number of grid points on the base grid is not known.
-
18:
– Logical array
Workspace
-
19:
– Integer
Input
-
On entry: the dimension of the array
lwk as declared in the (sub)program from which
d03rbf is called.
The required value of
lenlwk cannot be determined exactly in advance, but a suggested value is
where
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
and an estimated required size is printed on the current error message unit (see
x04aaf).
Note: the size of
lenlwk cannot be checked upon initial entry to
d03rbf since the number of grid points on the base grid is not known.
-
20:
– Integer
Input
-
On entry: the level of trace information required from
d03rbf.
itrace may take the value
,
,
,
or
.
- No output is generated.
- Only warning messages are printed.
- 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 , is assumed and similarly if , is assumed.
The advisory messages are given in greater detail as
itrace increases. Setting
allows you to monitor the progress of the integration without possibly excessive information.
-
21:
– Integer
Input/Output
-
On entry: must be set to
or
, alternatively
or
.
- Starts the integration in time. pdedef is assumed to be serial.
- Continues the integration after an earlier exit from the routine. In this case, only the following parameters may be reset between calls to d03rbf: tout, dt, tols, tolt, opti, optr, itrace and ifail. pdedef is assumed to be serial.
-
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 .
-
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 .
Constraint:
or .
On exit:
, if
ind on input was
or
, or
, if
ind on input was
or
.
Note: for users of serial versions of the NAG Library, it is recommended that you only use
or
. See
Section 8 for more information on the use of
ind.
-
22:
– Integer
Input/Output
-
On entry:
ifail must be set to
,
. 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
is recommended. If the output of error messages is undesirable, then the value
is recommended. Otherwise, if you are not familiar with this argument, the recommended value is
.
When the value is used it is essential to test the value of ifail on exit.
On exit:
unless the routine detects an error or a warning has been flagged (see
Section 6).
6
Error Indicators and Warnings
If on entry
or
, explanatory error messages are output on the current error message unit (as defined by
x04aaf).
Errors or warnings detected by the routine:
-
On entry, .
Constraint: if , .
On entry, and . Note that was reset to default if zero on entry.
Constraint: if , .
On entry, and . Note that was reset to default if zero on entry.
Constraint: if , .
On entry, .
Constraint: .
On entry, and .
Constraint: .
On entry, .
Constraint: .
On entry, , and .
Constraint: .
On entry, and too large:
and .
On entry, and too small:
and .
On entry,
ind is not equal to
or
:
.
On entry, and .
Constraint: .
On entry, .
Constraint: .
On entry, , and .
Constraint: if , .
On entry, .
Constraint: .
On entry, .
Constraint: .
On entry, .
Constraint: .
On entry, and .
Constraint: .
On entry, too small:
.
-
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.
-
One or more of the workspace arrays are too small. Try increasing
itrace for more information.
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
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.
-
set to
in
monitr. Integration completed as far as
ts:
.
-
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
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
or decrease the value of
tols.
-
Invalid output argument from
inidom.
One or more of the output arguments of
inidom was incorrectly specified, i.e.,
| , |
or | xmax too close to xmin, |
or | , |
or | ymax too close to ymin, |
or | nx or , |
or | , |
or | , |
or | , |
or | , |
or | , |
or | , |
or | or , for some , |
or | , for some , |
or | or , for some , |
or | , for some , |
or | or , for some , |
or | or , for some , |
or | , for some , |
or | , , , , , , , , , , or , for some , |
or | or , for some . |
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.
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.
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
Section 9, 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
d03rbf requires a user-supplied routine
pdedef to evaluate the functions
, for
. The parallelism within
d03rbf 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
d03rbf):
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
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 d03rbf 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
d03rbf 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
or
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
or
. If
ind is set to
or
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
d03rbf will be included in a future NAG Library release.
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 .
-
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 with coarse base as current grid.
For each grid point
a space monitor
is determined by
where
and
are the grid widths in the
and
directions; and
,
are the
coordinates at grid point
. The argument
is obtained from
where
is the user-supplied space tolerance;
is a weighting factor for the relative importance of the
th PDE component on the space monitor; and
is the approximate maximum absolute value of the
th component. A value for
must be supplied by you. Values for
and
must also be supplied but may be set to the values
if little information about the solution is known.
A new level of refinement is created if
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
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.
The time integration is controlled using a time monitor calculated at each level
up to the maximum level used, given by
where
is the total number of points on grid level
;
;
is the current time step;
is the time derivative of
which is approximated by first-order finite differences;
is the time equivalent of the space weighting factor
; and
is given by
where
is as before, and
is the user-specified time tolerance.
An integration step is rejected and retried at all levels if
10
Example
This example is taken from
Blom and Verwer (1993) and is the two-dimensional Burgers' system
with
on the domain given in
Figure 3. Dirichlet boundary conditions are used on all boundaries using the exact solution
The solution contains a wave front at which propagates in a direction perpendicular to the front with speed .
10.1
Program Text
10.2
Program Data
10.3
Program Results