NAG C Library Function Document

nag_pde_parab_1d_fd_ode_remesh (d03ppc)

1
Purpose

nag_pde_parab_1d_fd_ode_remesh (d03ppc) integrates a system of linear or nonlinear parabolic partial differential equations (PDEs) in one space variable, with scope for coupled ordinary differential equations (ODEs), and automatic adaptive spatial remeshing. The spatial discretization is performed using finite differences, and the method of lines is employed to reduce the PDEs to a system of ODEs. The resulting system is solved using a Backward Differentiation Formula (BDF) method or a Theta method (switching between Newton's method and functional iteration).

2
Specification

#include <nag.h>
#include <nagd03.h>
void  nag_pde_parab_1d_fd_ode_remesh (Integer npde, Integer m, double *ts, double tout,
void (*pdedef)(Integer npde, double t, double x, const double u[], const double ux[], Integer nv, const double v[], const double vdot[], double p[], double q[], double r[], Integer *ires, Nag_Comm *comm),
void (*bndary)(Integer npde, double t, const double u[], const double ux[], Integer nv, const double v[], const double vdot[], Integer ibnd, double beta[], double gamma[], Integer *ires, Nag_Comm *comm),
void (*uvinit)(Integer npde, Integer npts, Integer nxi, const double x[], const double xi[], double u[], Integer nv, double v[], Nag_Comm *comm),
double u[], Integer npts, double x[], Integer nv,
void (*odedef)(Integer npde, double t, Integer nv, const double v[], const double vdot[], Integer nxi, const double xi[], const double ucp[], const double ucpx[], const double rcp[], const double ucpt[], const double ucptx[], double f[], Integer *ires, Nag_Comm *comm),
Integer nxi, const double xi[], Integer neqn, const double rtol[], const double atol[], Integer itol, Nag_NormType norm, Nag_LinAlgOption laopt, const double algopt[], Nag_Boolean remesh, Integer nxfix, const double xfix[], Integer nrmesh, double dxmesh, double trmesh, Integer ipminf, double xratio, double con,
void (*monitf)(double t, Integer npts, Integer npde, const double x[], const double u[], const double r[], double fmon[], Nag_Comm *comm),
double rsave[], Integer lrsave, Integer isave[], Integer lisave, Integer itask, Integer itrace, const char *outfile, Integer *ind, Nag_Comm *comm, Nag_D03_Save *saved, NagError *fail)

3
Description

nag_pde_parab_1d_fd_ode_remesh (d03ppc) integrates the system of parabolic-elliptic equations and coupled ODEs
j=1npdePi,j Uj t +Qi=x-m x xmRi,  i=1,2,,npde ,   axb,tt0, (1)
Fit,V,V.,ξ,U*,Ux*,R*,Ut*,Uxt*=0,  i=1,2,,nv, (2)
where (1) defines the PDE part and (2) generalizes the coupled ODE part of the problem.
In (1), Pi,j and Ri depend on x, t, U, Ux, and V; Qi depends on x, t, U, Ux, V and linearly on V.. The vector U is the set of PDE solution values
Ux,t=U1x,t,,Unpdex,tT,  
and the vector Ux is the partial derivative with respect to x. The vector V is the set of ODE solution values
Vt=V1t,,VnvtT,  
and V. denotes its derivative with respect to time.
In (2), ξ represents a vector of nξ spatial coupling points at which the ODEs are coupled to the PDEs. These points may or may not be equal to some of the PDE spatial mesh points. U*, Ux*, R*, Ut* and Uxt* are the functions U, Ux, R, Ut and Uxt evaluated at these coupling points. Each Fi may only depend linearly on time derivatives. Hence the equation (2) may be written more precisely as
F=G-AV.-B Ut* Uxt* , (3)
where F=F1,,FnvT, G is a vector of length nv, A is an nv by nv matrix, B is an nv by nξ×npde matrix and the entries in G, A and B may depend on t, ξ, U*, Ux* and V. In practice you only need to supply a vector of information to define the ODEs and not the matrices A and B. (See Section 5 for the specification of odedef.)
The integration in time is from t0 to tout, over the space interval axb, where a=x1 and b=xnpts are the leftmost and rightmost points of a mesh x1,x2,,xnpts defined initially by you and (possibly) adapted automatically during the integration according to user-specified criteria. The coordinate system in space is defined by the following values of m; m=0 for Cartesian coordinates, m=1 for cylindrical polar coordinates and m=2 for spherical polar coordinates.
The PDE system which is defined by the functions Pi,j, Qi and Ri must be specified in pdedef.
The initial t=t0 values of the functions Ux,t and Vt must be specified in uvinit. Note that uvinit will be called again following any initial remeshing, and so Ux,t0 should be specified for all values of x in the interval axb, and not just the initial mesh points.
The functions Ri which may be thought of as fluxes, are also used in the definition of the boundary conditions. The boundary conditions must have the form
βix,tRix,t,U,Ux,V=γix,t,U,Ux,V,V.,  i=1,2,,npde, (4)
where x=a or x=b.
The boundary conditions must be specified in bndary. The function γi may depend linearly on V..
The problem is subject to the following restrictions:
(i) In (1), V.jt, for j=1,2,,nv, may only appear linearly in the functions Qi, for i=1,2,,npde, with a similar restriction for γ;
(ii) Pi,j and the flux Ri must not depend on any time derivatives;
(iii) t0<tout, so that integration is in the forward direction;
(iv) The evaluation of the terms Pi,j, Qi and Ri is done approximately at the mid-points of the mesh x[i-1], for i=1,2,,npts, by calling the pdedef for each mid-point in turn. Any discontinuities in these functions must therefore be at one or more of the fixed mesh points specified by xfix;
(v) At least one of the functions Pi,j must be nonzero so that there is a time derivative present in the PDE problem;
(vi) If m>0 and x1=0.0, which is the left boundary point, then it must be ensured that the PDE solution is bounded at this point. This can be done by either specifying the solution at x=0.0 or by specifying a zero flux there, that is βi=1.0 and γi=0.0. See also Section 9.
The algebraic-differential equation system which is defined by the functions Fi must be specified in odedef. You must also specify the coupling points ξ in the array xi.
The parabolic equations are approximated by a system of ODEs in time for the values of Ui at mesh points. For simple problems in Cartesian coordinates, this system is obtained by replacing the space derivatives by the usual central, three-point finite difference formula. However, for polar and spherical problems, or problems with nonlinear coefficients, the space derivatives are replaced by a modified three-point formula which maintains second order accuracy. In total there are npde×npts+nv ODEs in time direction. This system is then integrated forwards in time using a Backward Differentiation Formula (BDF) or a Theta method.
The adaptive space remeshing can be used to generate meshes that automatically follow the changing time-dependent nature of the solution, generally resulting in a more efficient and accurate solution using fewer mesh points than may be necessary with a fixed uniform or non-uniform mesh. Problems with travelling wavefronts or variable-width boundary layers for example will benefit from using a moving adaptive mesh. The discrete time-step method used here (developed by Furzeland (1984)) automatically creates a new mesh based on the current solution profile at certain time-steps, and the solution is then interpolated onto the new mesh and the integration continues.
The method requires you to supply a monitf which specifies in an analytical or numerical form the particular aspect of the solution behaviour you wish to track. This so-called monitor function is used to choose a mesh which equally distributes the integral of the monitor function over the domain. A typical choice of monitor function is the second space derivative of the solution value at each point (or some combination of the second space derivatives if there is more than one solution component), which results in refinement in regions where the solution gradient is changing most rapidly.
You must specify the frequency of mesh updates together with certain other criteria such as adjacent mesh ratios. Remeshing can be expensive and you are encouraged to experiment with the different options in order to achieve an efficient solution which adequately tracks the desired features of the solution.
Note that unless the monitor function for the initial solution values is zero at all user-specified initial mesh points, a new initial mesh is calculated and adopted according to the user-specified remeshing criteria. uvinit will then be called again to determine the initial solution values at the new mesh points (there is no interpolation at this stage) and the integration proceeds.

4
References

Berzins M (1990) Developments in the NAG Library software for parabolic equations Scientific Software Systems (eds J C Mason and M G Cox) 59–72 Chapman and Hall
Berzins M, Dew P M and Furzeland R M (1989) Developing software for time-dependent problems using the method of lines and differential-algebraic integrators Appl. Numer. Math. 5 375–397
Berzins M and Furzeland R M (1992) An adaptive theta method for the solution of stiff and nonstiff differential equations Appl. Numer. Math. 9 1–19
Furzeland R M (1984) The construction of adaptive space meshes TNER.85.022 Thornton Research Centre, Chester
Skeel R D and Berzins M (1990) A method for the spatial discretization of parabolic equations in one space variable SIAM J. Sci. Statist. Comput. 11(1) 1–32

5
Arguments

1:     npde IntegerInput
On entry: the number of PDEs to be solved.
Constraint: npde1.
2:     m IntegerInput
On entry: the coordinate system used:
m=0
Indicates Cartesian coordinates.
m=1
Indicates cylindrical polar coordinates.
m=2
Indicates spherical polar coordinates.
Constraint: m=0, 1 or 2.
3:     ts double *Input/Output
On entry: the initial value of the independent variable t.
On exit: the value of t corresponding to the solution values in u. Normally ts=tout.
Constraint: ts<tout.
4:     tout doubleInput
On entry: the final value of t to which the integration is to be carried out.
5:     pdedef function, supplied by the userExternal Function
pdedef must evaluate the functions Pi,j, Qi and Ri which define the system of PDEs. The functions may depend on x, t, U, Ux and V. Qi may depend linearly on V.. pdedef is called approximately midway between each pair of mesh points in turn by nag_pde_parab_1d_fd_ode_remesh (d03ppc).
The specification of pdedef is:
void  pdedef (Integer npde, double t, double x, const double u[], const double ux[], Integer nv, const double v[], const double vdot[], double p[], double q[], double r[], Integer *ires, Nag_Comm *comm)
1:     npde IntegerInput
On entry: the number of PDEs in the system.
2:     t doubleInput
On entry: the current value of the independent variable t.
3:     x doubleInput
On entry: the current value of the space variable x.
4:     u[npde] const doubleInput
On entry: u[i-1] contains the value of the component Uix,t, for i=1,2,,npde.
5:     ux[npde] const doubleInput
On entry: ux[i-1] contains the value of the component Uix,t x , for i=1,2,,npde.
6:     nv IntegerInput
On entry: the number of coupled ODEs in the system.
7:     v[nv] const doubleInput
On entry: if nv>0, v[i-1] contains the value of the component Vit, for i=1,2,,nv.
8:     vdot[nv] const doubleInput
On entry: if nv>0, vdot[i-1] contains the value of component V.it, for i=1,2,,nv.
Note: V.it, for i=1,2,,nv, may only appear linearly in Qj, for j=1,2,,npde.
9:     p[npde×npde] doubleOutput
Note: the i,jth element of the matrix P is stored in p[j-1×npde+i-1].
On exit: p[j-1×npde+i-1] must be set to the value of Pi,jx,t,U,Ux,V, for i=1,2,,npde and j=1,2,,npde.
10:   q[npde] doubleOutput
On exit: q[i-1] must be set to the value of Qix,t,U,Ux,V,V., for i=1,2,,npde.
11:   r[npde] doubleOutput
On exit: r[i-1] must be set to the value of Rix,t,U,Ux,V, for i=1,2,,npde.
12:   ires Integer *Input/Output
On entry: set to -1 or 1.
On exit: should usually remain unchanged. However, you may set ires to force the integration function to take certain actions as described below:
ires=2
Indicates to the integrator that control should be passed back immediately to the calling function with the error indicator set to fail.code= NE_USER_STOP.
ires=3
Indicates to the integrator that the current time step should be abandoned and a smaller time step used instead. You may wish to set ires=3 when a physically meaningless input or output value has been generated. If you consecutively set ires=3, nag_pde_parab_1d_fd_ode_remesh (d03ppc) returns to the calling function with the error indicator set to fail.code= NE_FAILED_DERIV.
13:   comm Nag_Comm *
Pointer to structure of type Nag_Comm; the following members are relevant to pdedef.
userdouble *
iuserInteger *
pPointer 
The type Pointer will be void *. Before calling nag_pde_parab_1d_fd_ode_remesh (d03ppc) you may allocate memory and initialize these pointers with various quantities for use by pdedef when called from nag_pde_parab_1d_fd_ode_remesh (d03ppc) (see Section 3.3.1.1 in How to Use the NAG Library and its Documentation).
Note: pdedef should not return floating-point NaN (Not a Number) or infinity values, since these are not handled by nag_pde_parab_1d_fd_ode_remesh (d03ppc). If your code inadvertently does return any NaNs or infinities, nag_pde_parab_1d_fd_ode_remesh (d03ppc) is likely to produce unexpected results.
6:     bndary function, supplied by the userExternal Function
bndary must evaluate the functions βi and γi which describe the boundary conditions, as given in (4).
The specification of bndary is:
void  bndary (Integer npde, double t, const double u[], const double ux[], Integer nv, const double v[], const double vdot[], Integer ibnd, double beta[], double gamma[], Integer *ires, Nag_Comm *comm)
1:     npde IntegerInput
On entry: the number of PDEs in the system.
2:     t doubleInput
On entry: the current value of the independent variable t.
3:     u[npde] const doubleInput
On entry: u[i-1] contains the value of the component Uix,t at the boundary specified by ibnd, for i=1,2,,npde.
4:     ux[npde] const doubleInput
On entry: ux[i-1] contains the value of the component Uix,t x  at the boundary specified by ibnd, for i=1,2,,npde.
5:     nv IntegerInput
On entry: the number of coupled ODEs in the system.
6:     v[nv] const doubleInput
On entry: if nv>0, v[i-1] contains the value of the component Vit, for i=1,2,,nv.
7:     vdot[nv] const doubleInput
On entry: vdot[i-1] contains the value of component V.it, for i=1,2,,nv.
Note: V.it, for i=1,2,,nv, may only appear linearly in γj, for j=1,2,,npde.
8:     ibnd IntegerInput
On entry: specifies which boundary conditions are to be evaluated.
ibnd=0
bndary must set up the coefficients of the left-hand boundary, x=a.
ibnd0
bndary must set up the coefficients of the right-hand boundary, x=b.
9:     beta[npde] doubleOutput
On exit: beta[i-1] must be set to the value of βix,t at the boundary specified by ibnd, for i=1,2,,npde.
10:   gamma[npde] doubleOutput
On exit: gamma[i-1] must be set to the value of γix,t,U,Ux,V,V. at the boundary specified by ibnd, for i=1,2,,npde.
11:   ires Integer *Input/Output
On entry: set to -1 or 1.
On exit: should usually remain unchanged. However, you may set ires to force the integration function to take certain actions as described below:
ires=2
Indicates to the integrator that control should be passed back immediately to the calling function with the error indicator set to fail.code= NE_USER_STOP.
ires=3
Indicates to the integrator that the current time step should be abandoned and a smaller time step used instead. You may wish to set ires=3 when a physically meaningless input or output value has been generated. If you consecutively set ires=3, nag_pde_parab_1d_fd_ode_remesh (d03ppc) returns to the calling function with the error indicator set to fail.code= NE_FAILED_DERIV.
12:   comm Nag_Comm *
Pointer to structure of type Nag_Comm; the following members are relevant to bndary.
userdouble *
iuserInteger *
pPointer 
The type Pointer will be void *. Before calling nag_pde_parab_1d_fd_ode_remesh (d03ppc) you may allocate memory and initialize these pointers with various quantities for use by bndary when called from nag_pde_parab_1d_fd_ode_remesh (d03ppc) (see Section 3.3.1.1 in How to Use the NAG Library and its Documentation).
Note: bndary should not return floating-point NaN (Not a Number) or infinity values, since these are not handled by nag_pde_parab_1d_fd_ode_remesh (d03ppc). If your code inadvertently does return any NaNs or infinities, nag_pde_parab_1d_fd_ode_remesh (d03ppc) is likely to produce unexpected results.
7:     uvinit function, supplied by the userExternal Function
uvinit must supply the initial t=t0 values of Ux,t and Vt for all values of x in the interval axb.
The specification of uvinit is:
void  uvinit (Integer npde, Integer npts, Integer nxi, const double x[], const double xi[], double u[], Integer nv, double v[], Nag_Comm *comm)
1:     npde IntegerInput
On entry: the number of PDEs in the system.
2:     npts IntegerInput
On entry: the number of mesh points in the interval a,b.
3:     nxi IntegerInput
On entry: the number of ODE/PDE coupling points.
4:     x[npts] const doubleInput
On entry: the current mesh. x[i-1] contains the value of xi, for i=1,2,,npts.
5:     xi[nxi] const doubleInput
On entry: if nxi>0, xi[i-1] contains the value of the ODE/PDE coupling point, ξi, for i=1,2,,nxi.
6:     u[npde×npts] doubleOutput
Note: the i,jth element of the matrix U is stored in u[j-1×npde+i-1].
On exit: u[j-1×npde+i-1] contains the value of the component Uixj,t0, for i=1,2,,npde and j=1,2,,npts.
7:     nv IntegerInput
On entry: the number of coupled ODEs in the system.
8:     v[nv] doubleOutput
On exit: v[i-1] contains the value of component Vit0, for i=1,2,,nv.
9:     comm Nag_Comm *
Pointer to structure of type Nag_Comm; the following members are relevant to uvinit.
userdouble *
iuserInteger *
pPointer 
The type Pointer will be void *. Before calling nag_pde_parab_1d_fd_ode_remesh (d03ppc) you may allocate memory and initialize these pointers with various quantities for use by uvinit when called from nag_pde_parab_1d_fd_ode_remesh (d03ppc) (see Section 3.3.1.1 in How to Use the NAG Library and its Documentation).
Note: uvinit should not return floating-point NaN (Not a Number) or infinity values, since these are not handled by nag_pde_parab_1d_fd_ode_remesh (d03ppc). If your code inadvertently does return any NaNs or infinities, nag_pde_parab_1d_fd_ode_remesh (d03ppc) is likely to produce unexpected results.
8:     u[neqn] doubleInput/Output
On entry: if ind=1, the value of u must be unchanged from the previous call.
On exit: u[npde×j-1+i-1] contains the computed solution Uixj,t, for i=1,2,,npde and j=1,2,,npts, and u[npts×npde+k-1] contains Vkt, for k=1,2,,nv, evaluated at t=ts.
9:     npts IntegerInput
On entry: the number of mesh points in the interval a,b.
Constraint: npts3.
10:   x[npts] doubleInput/Output
On entry: the initial mesh points in the space direction. x[0] must specify the left-hand boundary, a, and x[npts-1] must specify the right-hand boundary, b.
Constraint: x[0]<x[1]<<x[npts-1].
On exit: the final values of the mesh points.
11:   nv IntegerInput
On entry: the number of coupled ODE in the system.
Constraint: nv0.
12:   odedef function, supplied by the userExternal Function
odedef must evaluate the functions F, which define the system of ODEs, as given in (3).
odedef will never be called and the NAG defined null void function pointer, NULLFN, can be supplied in the call to nag_pde_parab_1d_fd_ode_remesh (d03ppc).
The specification of odedef is:
void  odedef (Integer npde, double t, Integer nv, const double v[], const double vdot[], Integer nxi, const double xi[], const double ucp[], const double ucpx[], const double rcp[], const double ucpt[], const double ucptx[], double f[], Integer *ires, Nag_Comm *comm)
1:     npde IntegerInput
On entry: the number of PDEs in the system.
2:     t doubleInput
On entry: the current value of the independent variable t.
3:     nv IntegerInput
On entry: the number of coupled ODEs in the system.
4:     v[nv] const doubleInput
On entry: if nv>0, v[i-1] contains the value of the component Vit, for i=1,2,,nv.
5:     vdot[nv] const doubleInput
On entry: if nv>0, vdot[i-1] contains the value of component V.it, for i=1,2,,nv.
6:     nxi IntegerInput
On entry: the number of ODE/PDE coupling points.
7:     xi[nxi] const doubleInput
On entry: if nxi>0, xi[i-1] contains the ODE/PDE coupling points, ξi, for i=1,2,,nxi.
8:     ucp[npde×nxi] const doubleInput
Note: the i,jth element of the matrix is stored in ucp[j-1×npde+i-1].
On entry: if nxi>0, ucp[j-1×npde+i-1] contains the value of Uix,t at the coupling point x=ξj, for i=1,2,,npde and j=1,2,,nxi.
9:     ucpx[npde×nxi] const doubleInput
Note: the i,jth element of the matrix is stored in ucpx[j-1×npde+i-1].
On entry: if nxi>0, ucpx[j-1×npde+i-1] contains the value of Uix,t x  at the coupling point x=ξj, for i=1,2,,npde and j=1,2,,nxi.
10:   rcp[npde×nxi] const doubleInput
Note: the i,jth element of the matrix is stored in rcp[j-1×npde+i-1].
On entry: rcp[j-1×npde+i-1] contains the value of the flux Ri at the coupling point x=ξj, for i=1,2,,npde and j=1,2,,nxi.
11:   ucpt[npde×nxi] const doubleInput
Note: the i,jth element of the matrix is stored in ucpt[j-1×npde+i-1].
On entry: if nxi>0, ucpt[j-1×npde+i-1] contains the value of Ui t  at the coupling point x=ξj, for i=1,2,,npde and j=1,2,,nxi.
12:   ucptx[npde×nxi] const doubleInput
Note: the i,jth element of the matrix is stored in ucptx[j-1×npde+i-1].
On entry: ucptx[j-1×npde+i-1] contains the value of 2Ui x t  at the coupling point x=ξj, for i=1,2,,npde and j=1,2,,nxi.
13:   f[nv] doubleOutput
On exit: f[i-1] must contain the ith component of F, for i=1,2,,nv, where F is defined as
F=G-AV.-B Ut* Uxt* , (5)
or
F=-AV.-B Ut* Uxt* . (6)
The definition of F is determined by the input value of ires.
14:   ires Integer *Input/Output
On entry: the form of F that must be returned in the array f.
ires=1
Equation (5) must be used.
ires=-1
Equation (6) must be used.
On exit: should usually remain unchanged. However, you may reset ires to force the integration function to take certain actions as described below:
ires=2
Indicates to the integrator that control should be passed back immediately to the calling function with the error indicator set to fail.code= NE_USER_STOP.
ires=3
Indicates to the integrator that the current time step should be abandoned and a smaller time step used instead. You may wish to set ires=3 when a physically meaningless input or output value has been generated. If you consecutively set ires=3, nag_pde_parab_1d_fd_ode_remesh (d03ppc) returns to the calling function with the error indicator set to fail.code= NE_FAILED_DERIV.
15:   comm Nag_Comm *
Pointer to structure of type Nag_Comm; the following members are relevant to odedef.
userdouble *
iuserInteger *
pPointer 
The type Pointer will be void *. Before calling nag_pde_parab_1d_fd_ode_remesh (d03ppc) you may allocate memory and initialize these pointers with various quantities for use by odedef when called from nag_pde_parab_1d_fd_ode_remesh (d03ppc) (see Section 3.3.1.1 in How to Use the NAG Library and its Documentation).
Note: odedef should not return floating-point NaN (Not a Number) or infinity values, since these are not handled by nag_pde_parab_1d_fd_ode_remesh (d03ppc). If your code inadvertently does return any NaNs or infinities, nag_pde_parab_1d_fd_ode_remesh (d03ppc) is likely to produce unexpected results.
13:   nxi IntegerInput
On entry: the number of ODE/PDE coupling points.
Constraints:
  • if nv=0, nxi=0;
  • if nv>0, nxi0.
14:   xi[nxi] const doubleInput
On entry: if nxi>0, xi[i-1], for i=1,2,,nxi, must be set to the ODE/PDE coupling points.
Constraint: x[0]xi[0]<xi[1]<<xi[nxi-1]x[npts-1].
15:   neqn IntegerInput
On entry: the number of ODEs in the time direction.
Constraint: neqn=npde×npts+nv.
16:   rtol[dim] const doubleInput
Note: the dimension, dim, of the array rtol must be at least
  • 1 when itol=1 or 2;
  • neqn when itol=3 or 4.
On entry: the relative local error tolerance.
Constraint: rtol[i-1]0.0 for all relevant i.
17:   atol[dim] const doubleInput
Note: the dimension, dim, of the array atol must be at least
  • 1 when itol=1 or 3;
  • neqn when itol=2 or 4.
On entry: the absolute local error tolerance.
Constraints:
  • atol[i-1]0.0 for all relevant i;
  • Corresponding elements of atol and rtol cannot both be 0.0.
18:   itol IntegerInput
On entry: a value to indicate the form of the local error test. itol indicates to nag_pde_parab_1d_fd_ode_remesh (d03ppc) whether to interpret either or both of rtol or atol as a vector or scalar. The error test to be satisfied is ei/wi<1.0, where wi is defined as follows:
itolrtolatolwi
1scalarscalarrtol[0]×Ui+atol[0]
2scalarvectorrtol[0]×Ui+atol[i-1]
3vectorscalarrtol[i-1]×Ui+atol[0]
4vectorvectorrtol[i-1]×Ui+atol[i-1]
In the above, ei denotes the estimated local error for the ith component of the coupled PDE/ODE system in time, u[i-1], for i=1,2,,neqn.
The choice of norm used is defined by the argument norm.
Constraint: 1itol4.
19:   norm Nag_NormTypeInput
On entry: the type of norm to be used.
norm=Nag_MaxNorm
Maximum norm.
norm=Nag_TwoNorm
Averaged L2 norm.
If unorm denotes the norm of the vector u of length neqn, then for the averaged L2 norm
unorm=1neqni=1neqnu[i-1]/wi2,  
while for the maximum norm
u norm = maxi u[i-1] / wi .  
See the description of itol for the formulation of the weight vector w.
Constraint: norm=Nag_MaxNorm or Nag_TwoNorm.
20:   laopt Nag_LinAlgOptionInput
On entry: the type of matrix algebra required.
laopt=Nag_LinAlgFull
Full matrix methods to be used.
laopt=Nag_LinAlgBand
Banded matrix methods to be used.
laopt=Nag_LinAlgSparse
Sparse matrix methods to be used.
Constraint: laopt=Nag_LinAlgFull, Nag_LinAlgBand or Nag_LinAlgSparse.
Note: you are recommended to use the banded option when no coupled ODEs are present (i.e., nv=0).
21:   algopt[30] const doubleInput
On entry: may be set to control various options available in the integrator. If you wish to employ all the default options, algopt[0] should be set to 0.0. Default values will also be used for any other elements of algopt set to zero. The permissible values, default values, and meanings are as follows:
algopt[0]
Selects the ODE integration method to be used. If algopt[0]=1.0, a BDF method is used and if algopt[0]=2.0, a Theta method is used. The default value is algopt[0]=1.0.
If algopt[0]=2.0, algopt[i-1], for i=2,3,4 are not used.
algopt[1]
Specifies the maximum order of the BDF integration formula to be used. algopt[1] may be 1.0, 2.0, 3.0, 4.0 or 5.0. The default value is algopt[1]=5.0.
algopt[2]
Specifies what method is to be used to solve the system of nonlinear equations arising on each step of the BDF method. If algopt[2]=1.0 a modified Newton iteration is used and if algopt[2]=2.0 a functional iteration method is used. If functional iteration is selected and the integrator encounters difficulty, there is an automatic switch to the modified Newton iteration. The default value is algopt[2]=1.0.
algopt[3]
Specifies whether or not the Petzold error test is to be employed. The Petzold error test results in extra overhead but is more suitable when algebraic equations are present, such as Pi,j=0.0, for j=1,2,,npde, for some i or when there is no V.it dependence in the coupled ODE system. If algopt[3]=1.0, the Petzold test is used. If algopt[3]=2.0, the Petzold test is not used. The default value is algopt[3]=1.0.
If algopt[0]=1.0, algopt[i-1], for i=5,6,7, are not used.
algopt[4]
Specifies the value of Theta to be used in the Theta integration method. 0.51algopt[4]0.99. The default value is algopt[4]=0.55.
algopt[5]
Specifies what method is to be used to solve the system of nonlinear equations arising on each step of the Theta method. If algopt[5]=1.0, a modified Newton iteration is used and if algopt[5]=2.0, a functional iteration method is used. The default value is algopt[5]=1.0.
algopt[6]
Specifies whether or not the integrator is allowed to switch automatically between modified Newton and functional iteration methods in order to be more efficient. If algopt[6]=1.0, switching is allowed and if algopt[6]=2.0, switching is not allowed. The default value is algopt[6]=1.0.
algopt[10]
Specifies a point in the time direction, tcrit, beyond which integration must not be attempted. The use of tcrit is described under the argument itask. If algopt[0]0.0, a value of 0.0 for algopt[10], say, should be specified even if itask subsequently specifies that tcrit will not be used.
algopt[11]
Specifies the minimum absolute step size to be allowed in the time integration. If this option is not required, algopt[11] should be set to 0.0.
algopt[12]
Specifies the maximum absolute step size to be allowed in the time integration. If this option is not required, algopt[12] should be set to 0.0.
algopt[13]
Specifies the initial step size to be attempted by the integrator. If algopt[13]=0.0, the initial step size is calculated internally.
algopt[14]
Specifies the maximum number of steps to be attempted by the integrator in any one call. If algopt[14]=0.0, no limit is imposed.
algopt[22]
Specifies what method is to be used to solve the nonlinear equations at the initial point to initialize the values of U, Ut, V and V.. If algopt[22]=1.0, a modified Newton iteration is used and if algopt[22]=2.0, functional iteration is used. The default value is algopt[22]=1.0.
algopt[28] and algopt[29] are used only for the sparse matrix algebra option, laopt=Nag_LinAlgSparse.
algopt[28]
Governs the choice of pivots during the decomposition of the first Jacobian matrix. It should lie in the range 0.0<algopt[28]<1.0, with smaller values biasing the algorithm towards maintaining sparsity at the expense of numerical stability. If algopt[28] lies outside this range then the default value is used. If the functions regard the Jacobian matrix as numerically singular then increasing algopt[28] towards 1.0 may help, but at the cost of increased fill-in. The default value is algopt[28]=0.1.
algopt[29]
Is used as a relative pivot threshold during subsequent Jacobian decompositions (see algopt[28]) below which an internal error is invoked. If algopt[29] is greater than 1.0 no check is made on the pivot size, and this may be a necessary option if the Jacobian is found to be numerically singular (see algopt[28]). The default value is algopt[29]=0.0001.
22:   remesh Nag_BooleanInput
On entry: indicates whether or not spatial remeshing should be performed.
remesh=Nag_TRUE
Indicates that spatial remeshing should be performed as specified.
remesh=Nag_FALSE
Indicates that spatial remeshing should be suppressed.
Note: remesh should not be changed between consecutive calls to nag_pde_parab_1d_fd_ode_remesh (d03ppc). Remeshing can be switched off or on at specified times by using appropriate values for the arguments nrmesh and trmesh at each call.
23:   nxfix IntegerInput
On entry: the number of fixed mesh points.
Constraint: 0nxfixnpts-2.
Note: the end points x[0] and x[npts-1] are fixed automatically and hence should not be specified as fixed points.
24:   xfix[nxfix] const doubleInput
On entry: xfix[i-1], for i=1,2,,nxfix, must contain the value of the x coordinate at the ith fixed mesh point.
Constraints:
  • xfix[i-1]<xfix[i], for i=1,2,,nxfix-1;
  • each fixed mesh point must coincide with a user-supplied initial mesh point, that is xfix[i-1]=x[j-1] for some j, 2jnpts-1.
Note: the positions of the fixed mesh points in the array x remain fixed during remeshing, and so the number of mesh points between adjacent fixed points (or between fixed points and end points) does not change. You should take this into account when choosing the initial mesh distribution.
25:   nrmesh IntegerInput
On entry: specifies the spatial remeshing frequency and criteria for the calculation and adoption of a new mesh.
nrmesh<0
Indicates that a new mesh is adopted according to the argument dxmesh. The mesh is tested every nrmesh timesteps.
nrmesh=0
Indicates that remeshing should take place just once at the end of the first time step reached when t>trmesh.
nrmesh>0
Indicates that remeshing will take place every nrmesh time steps, with no testing using dxmesh.
Note: nrmesh may be changed between consecutive calls to nag_pde_parab_1d_fd_ode_remesh (d03ppc) to give greater flexibility over the times of remeshing.
26:   dxmesh doubleInput
On entry: determines whether a new mesh is adopted when nrmesh is set less than zero. A possible new mesh is calculated at the end of every nrmesh time steps, but is adopted only if
xinew>xi old +dxmesh×xi+1 old -xi old  
or
xinew<xi old -dxmesh×xi old -xi- 1 old  
dxmesh thus imposes a lower limit on the difference between one mesh and the next.
Constraint: dxmesh0.0.
27:   trmesh doubleInput
On entry: specifies when remeshing will take place when nrmesh is set to zero. Remeshing will occur just once at the end of the first time step reached when t is greater than trmesh.
Note: trmesh may be changed between consecutive calls to nag_pde_parab_1d_fd_ode_remesh (d03ppc) to force remeshing at several specified times.
28:   ipminf IntegerInput
On entry: the level of trace information regarding the adaptive remeshing.
ipminf=0
No trace information.
ipminf=1
Brief summary of mesh characteristics.
ipminf=2
More detailed information, including old and new mesh points, mesh sizes and monitor function values.
Constraint: ipminf=0, 1 or 2.
29:   xratio doubleInput
On entry: an input bound on the adjacent mesh ratio (greater than 1.0 and typically in the range 1.5 to 3.0). The remeshing functions will attempt to ensure that
xi-xi-1/xratio<xi+1-xi<xratio×xi-xi-1. 
Suggested value: xratio=1.5.
Constraint: xratio>1.0.
30:   con doubleInput
On entry: an input bound on the sub-integral of the monitor function Fmonx over each space step. The remeshing functions will attempt to ensure that
xixi+1Fmonxdxconx1xnptsFmonxdx,  
(see Furzeland (1984)). con gives you more control over the mesh distribution e.g., decreasing con allows more clustering. A typical value is 2/npts-1, but you are encouraged to experiment with different values. Its value is not critical and the mesh should be qualitatively correct for all values in the range given below.
Suggested value: con=2.0/npts-1.
Constraint: 0.1/npts-1con10.0/npts-1.
31:   monitf function, supplied by the userExternal Function
monitf must supply and evaluate a remesh monitor function to indicate the solution behaviour of interest.
If you specify remesh=Nag_FALSE, i.e., no remeshing, monitf will not be called and may be specified as NULLFN.
The specification of monitf is:
void  monitf (double t, Integer npts, Integer npde, const double x[], const double u[], const double r[], double fmon[], Nag_Comm *comm)
1:     t doubleInput
On entry: the current value of the independent variable t.
2:     npts IntegerInput
On entry: the number of mesh points in the interval a,b.
3:     npde IntegerInput
On entry: the number of PDEs in the system.
4:     x[npts] const doubleInput
On entry: the current mesh. x[i-1] contains the value of xi, for i=1,2,,npts.
5:     u[npde×npts] const doubleInput
Note: the i,jth element of the matrix U is stored in u[j-1×npde+i-1].
On entry: u[j-1×npde+i-1] contains the value of Uix,t at x=x[j-1] and time t, for i=1,2,,npde and j=1,2,,npts.
6:     r[npde×npts] const doubleInput
Note: the i,jth element of the matrix R is stored in r[j-1×npde+i-1].
On entry: r[j-1×npde+i-1] contains the value of Rix,t,U,Ux,V at x=x[j-1] and time t, for i=1,2,,npde and j=1,2,,npts.
7:     fmon[npts] doubleOutput
On exit: fmon[i-1] must contain the value of the monitor function Fmonx at mesh point x=x[i-1].
Constraint: fmon[i-1]0.0.
8:     comm Nag_Comm *
Pointer to structure of type Nag_Comm; the following members are relevant to monitf.
userdouble *
iuserInteger *
pPointer 
The type Pointer will be void *. Before calling nag_pde_parab_1d_fd_ode_remesh (d03ppc) you may allocate memory and initialize these pointers with various quantities for use by monitf when called from nag_pde_parab_1d_fd_ode_remesh (d03ppc) (see Section 3.3.1.1 in How to Use the NAG Library and its Documentation).
Note: monitf should not return floating-point NaN (Not a Number) or infinity values, since these are not handled by nag_pde_parab_1d_fd_ode_remesh (d03ppc). If your code inadvertently does return any NaNs or infinities, nag_pde_parab_1d_fd_ode_remesh (d03ppc) is likely to produce unexpected results.
32:   rsave[lrsave] doubleCommunication Array
If ind=0, rsave need not be set on entry.
If ind=1, rsave must be unchanged from the previous call to the function because it contains required information about the iteration.
33:   lrsave IntegerInput
On entry: the dimension of the array rsave. Its size depends on the type of matrix algebra selected.
If laopt=Nag_LinAlgFull, lrsaveneqn×neqn+neqn+nwkres+lenode.
If laopt=Nag_LinAlgBand, lrsave3×mlu+1×neqn+nwkres+lenode.
If laopt=Nag_LinAlgSparse, lrsave4×neqn+11×neqn/2+1+nwkres+lenode.
Where
mlu is the lower or upper half bandwidths such that
mlu=2×npde-1, for PDE problems only; or
mlu=neqn-1, for coupled PDE/ODE problems.
nwkres= npde×3×npde+6×nxi+npts+15+nxi+nv+7×npts+nxfix+1, when ​nv>0​ and ​nxi>0; npde×3×npde+npts+21+nv+7×npts+nxfix+2, when ​nv>0​ and ​nxi=0; or npde×3×npde+npts+21+7×npts+nxfix+3, when ​nv=0.  
lenode= 6+intalgopt[1]×neqn+50, when the BDF method is used; 9×neqn+50, when the Theta method is used.  
Note: when using the sparse option, the value of lrsave may be too small when supplied to the integrator. An estimate of the minimum size of lrsave is printed on the current error message unit if itrace>0 and the function returns with fail.code= NE_INT_2.
34:   isave[lisave] IntegerCommunication Array
If ind=0, isave need not be set on entry.
If ind=1, isave must be unchanged from the previous call to the function because it contains required information about the iteration required for subsequent calls. In particular:
isave[0]
Contains the number of steps taken in time.
isave[1]
Contains the number of residual evaluations of the resulting ODE system used. One such evaluation involves computing the PDE functions at all the mesh points, as well as one evaluation of the functions in the boundary conditions.
isave[2]
Contains the number of Jacobian evaluations performed by the time integrator.
isave[3]
Contains the order of the ODE method last used in the time integration.
isave[4]
Contains the number of Newton iterations performed by the time integrator. Each iteration involves residual evaluation of the resulting ODE system followed by a back-substitution using the LU decomposition of the Jacobian matrix.
The rest of the array is used as workspace.
35:   lisave IntegerInput
On entry: the dimension of the array isave.
Its size depends on the type of matrix algebra selected:
  • if laopt=Nag_LinAlgBand, lisaveneqn+25+nxfix;
  • if laopt=Nag_LinAlgFull, lisave25+nxfix;
  • if laopt=Nag_LinAlgSparse, lisave25×neqn+25+nxfix.
Note: when using the sparse option, the value of lisave may be too small when supplied to the integrator. An estimate of the minimum size of lisave is printed if itrace>0 and the function returns with fail.code= NE_INT_2.
36:   itask IntegerInput
On entry: specifies the task to be performed by the ODE integrator.
itask=1
Normal computation of output values u at t=tout.
itask=2
One step and return.
itask=3
Stop at first internal integration point at or beyond t=tout.
itask=4
Normal computation of output values u at t=tout but without overshooting t=tcrit where tcrit is described under the argument algopt.
itask=5
Take one step in the time direction and return, without passing tcrit, where tcrit is described under the argument algopt.
Constraint: itask=1, 2, 3, 4 or 5.
37:   itrace IntegerInput
On entry: the level of trace information required from nag_pde_parab_1d_fd_ode_remesh (d03ppc) and the underlying ODE solver:
itrace-1
No output is generated.
itrace=0
Only warning messages from the PDE solver are printed.
itrace=1
Output from the underlying ODE solver is printed. This output contains details of Jacobian entries, the nonlinear iteration and the time integration during the computation of the ODE system.
itrace=2
Output from the underlying ODE solver is similar to that produced when itrace=1, except that the advisory messages are given in greater detail.
itrace3
Output from the underlying ODE solver is similar to that produced when itrace=2, except that the advisory messages are given in greater detail.
38:   outfile const char *Input
On entry: the name of a file to which diagnostic output will be directed. If outfile is NULL the diagnostic output will be directed to standard output.
39:   ind Integer *Input/Output
On entry: must be set to 0 or 1.
ind=0
Starts or restarts the integration in time.
ind=1
Continues the integration after an earlier exit from the function. In this case, only the arguments tout and fail and the remeshing arguments nrmesh, dxmesh, trmesh, xratio and con may be reset between calls to nag_pde_parab_1d_fd_ode_remesh (d03ppc).
Constraint: 0ind1.
On exit: ind=1.
40:   comm Nag_Comm *
The NAG communication argument (see Section 3.3.1.1 in How to Use the NAG Library and its Documentation).
41:   saved Nag_D03_Save *Communication Structure
saved must remain unchanged following a previous call to a Chapter d03 function and prior to any subsequent call to a Chapter d03 function.
42:   fail NagError *Input/Output
The NAG error argument (see Section 3.7 in How to Use the NAG Library and its Documentation).

6
Error Indicators and Warnings

NE_ACC_IN_DOUBT
Integration completed, but small changes in atol or rtol are unlikely to result in a changed solution.
The required task has been completed, but it is estimated that a small change in atol and rtol is unlikely to produce any change in the computed solution. (Only applies when you are not operating in one step mode, that is when itask2 or 5.)
NE_ALLOC_FAIL
Dynamic memory allocation failed.
See Section 2.3.1.2 in How to Use the NAG Library and its Documentation for further information.
NE_BAD_PARAM
On entry, argument value had an illegal value.
NE_FAILED_DERIV
In setting up the ODE system an internal auxiliary was unable to initialize the derivative. This could be due to your setting ires=3 in pdedef or bndary.
NE_FAILED_START
atol and rtol were too small to start integration.
Underlying ODE solver cannot make further progress from the point ts with the supplied values of atol and rtol. ts=value.
NE_FAILED_STEP
Error during Jacobian formulation for ODE system. Increase itrace for further details.
Repeated errors in an attempted step of underlying ODE solver. Integration was successful as far as ts: ts=value.
In the underlying ODE solver, there were repeated error test failures on an attempted step, before completing the requested task, but the integration was successful as far as t=ts. The problem may have a singularity, or the error requirement may be inappropriate.
NE_INCOMPAT_PARAM
On entry, con=value, npts=value.
Constraint: con10.0/npts-1.
On entry, con=value, npts=value.
Constraint: con0.1/npts-1.
On entry, m=value and x[0]=value.
Constraint: m0 or x[0]0.0 
On entry, the point xfix[i-1] does not coincide with any x[j-1]: i=value and xfix[i-1]=value.
NE_INT
ires set to an invalid value in call to pdedef, bndary, or odedef.
On entry, ind=value.
Constraint: ind=0 or 1.
On entry, ipminf=value.
Constraint: ipminf=0, 1 or 2.
On entry, itask=value.
Constraint: itask=1, 2, 3, 4 or 5.
On entry, itol=value.
Constraint: itol=1, 2, 3 or 4.
On entry, m=value.
Constraint: m=0, 1 or 2.
On entry, npde=value.
Constraint: npde1.
On entry, npts=value.
Constraint: npts3.
On entry, nv=value.
Constraint: nv0.
On entry, nxfix=value.
Constraint: nxfix0.
On entry, on initial entry ind=1.
Constraint: on initial entry ind=0.
NE_INT_2
On entry, i=value and j=value.
Constraint: corresponding elements atol[i-1] and rtol[j-1] cannot both be 0.0.
On entry, lisave=value.
Constraint: lisavevalue.
On entry, lrsave=value.
Constraint: lrsavevalue.
On entry, nv=value and nxi=value.
Constraint: nxi=0 when nv=0.
On entry, nv=value and nxi=value.
Constraint: nxi0 when nv>0.
On entry, nxfix=value, npts=value.
Constraint: nxfixnpts-2.
When using the sparse option lisave or lrsave is too small: lisave=value, lrsave=value.
NE_INT_4
On entry, neqn=value, npde=value, npts=value and nv=value.
Constraint: neqn=npde×npts+nv.
NE_INTERNAL_ERROR
An internal error has occurred in this function. Check the function call and any array sizes. If the call is correct then please contact NAG for assistance.
See Section 2.7.6 in How to Use the NAG Library and its Documentation for further information.
Serious error in internal call to an auxiliary. Increase itrace for further details.
NE_ITER_FAIL
In solving ODE system, the maximum number of steps algopt[14] has been exceeded. algopt[14]=value.
NE_NO_LICENCE
Your licence key may have expired or may not have been installed correctly.
See Section 2.7.5 in How to Use the NAG Library and its Documentation for further information.
NE_NOT_CLOSE_FILE
Cannot close file value.
NE_NOT_STRICTLY_INCREASING
On entry, i=value, xfix[i]=value and xfix[i-1]=value.
Constraint: xfix[i]>xfix[i-1].
On entry, i=value, x[i-1]=value, j=value and x[j-1]=value.
Constraint: x[0]<x[1]<<x[npts-1].
On entry, i=value, xi[i]=value and xi[i-1]=value.
Constraint: xi[i]>xi[i-1].
NE_NOT_WRITE_FILE
Cannot open file value for writing.
NE_REAL
On entry, dxmesh=value.
Constraint: dxmesh0.0.
On entry, xratio=value.
Constraint: xratio>1.0.
NE_REAL_2
On entry, at least one point in xi lies outside x[0],x[npts-1]: x[0]=value and x[npts-1]=value.
On entry, tout=value and ts=value.
Constraint: tout>ts.
On entry, tout-ts is too small: tout=value and ts=value.
NE_REAL_ARRAY
On entry, i=value and atol[i-1]=value.
Constraint: atol[i-1]0.0.
On entry, i=value and rtol[i-1]=value.
Constraint: rtol[i-1]0.0.
NE_REMESH_CHANGED
remesh has been changed between calls to nag_pde_parab_1d_fd_ode_remesh (d03ppc).
NE_SING_JAC
Singular Jacobian of ODE system. Check problem formulation.
NE_TIME_DERIV_DEP
Flux function appears to depend on time derivatives.
NE_USER_STOP
In evaluating residual of ODE system, ires=2 has been set in pdedef, bndary, or odedef. Integration is successful as far as ts: ts=value.
NE_ZERO_WTS
Zero error weights encountered during time integration.
Some error weights wi became zero during the time integration (see the description of itol). Pure relative error control (atol[i-1]=0.0) was requested on a variable (the ith) which has become zero. The integration was successful as far as t=ts.

7
Accuracy

nag_pde_parab_1d_fd_ode_remesh (d03ppc) controls the accuracy of the integration in the time direction but not the accuracy of the approximation in space. The spatial accuracy depends on both the number of mesh points and on their distribution in space. In the time integration only the local error over a single step is controlled and so the accuracy over a number of steps cannot be guaranteed. You should therefore test the effect of varying the accuracy arguments, atol and rtol.

8
Parallelism and Performance

nag_pde_parab_1d_fd_ode_remesh (d03ppc) is threaded by NAG for parallel execution in multithreaded implementations of the NAG Library.
nag_pde_parab_1d_fd_ode_remesh (d03ppc) makes calls to BLAS and/or LAPACK routines, which may be threaded within the vendor library used by this implementation. Consult the documentation for the vendor library for further information.
Please consult the x06 Chapter Introduction for information on how to control and interrogate the OpenMP environment used within this function. Please also consult the Users' Note for your implementation for any additional implementation-specific information.

9
Further Comments

The argument specification allows you to include equations with only first-order derivatives in the space direction but there is no guarantee that the method of integration will be satisfactory for such systems. The position and nature of the boundary conditions in particular are critical in defining a stable problem. It may be advisable in such cases to reduce the whole system to first-order and to use the Keller box scheme function nag_pde_parab_1d_keller_ode_remesh (d03prc).
The time taken depends on the complexity of the parabolic system, the accuracy requested, and the frequency of the mesh updates. For a given system with fixed accuracy and mesh-update frequency it is approximately proportional to neqn.

10
Example

This example uses Burgers Equation, a common test problem for remeshing algorithms, given by
U t =-U U x +E 2U x2 ,  
for x0,1 and t0,1, where E is a small constant.
The initial and boundary conditions are given by the exact solution
Ux,t=0.1exp-A+0.5exp-B+exp-C exp-A+exp-B+exp-C ,  
where
A = 50Ex- 0.5+ 4.95t, B = 250Ex- 0.5+ 0.75t, C = 500Ex- 0.375.  

10.1
Program Text

Program Text (d03ppce.c)

10.2
Program Data

None.

10.3
Program Results

Program Results (d03ppce.r)

GnuplotProduced by GNUPLOT 5.0 patchlevel 0 Example Program Solution of Burgers Equation using Moving Mesh U(x,t) gnuplot_plot_1 gnuplot_plot_2 0 0.1 0.2 0.3 0.4 0.5 0.6 0.7 0.8 0.9 1 Time 0 0.2 0.4 0.6 0.8 1 x 0.1 0.5 1