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

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

NAG Toolbox: nag_pde_1d_parab_remesh_keller (d03pr)

Purpose

nag_pde_1d_parab_remesh_keller (d03pr) integrates a system of linear or nonlinear, first-order, time-dependent 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 the Keller box scheme (see Keller (1970)) 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).

Syntax

[ts, u, x, rsave, isave, ind, ifail] = d03pr(npde, ts, tout, pdedef, bndary, uvinit, u, x, nleft, ncode, odedef, xi, rtol, atol, itol, norm_p, laopt, algopt, remesh, xfix, nrmesh, dxmesh, trmesh, ipminf, monitf, rsave, isave, itask, itrace, ind, 'npts', npts, 'nxi', nxi, 'neqn', neqn, 'nxfix', nxfix, 'xratio', xratio, 'con', con)
[ts, u, x, rsave, isave, ind, ifail] = nag_pde_1d_parab_remesh_keller(npde, ts, tout, pdedef, bndary, uvinit, u, x, nleft, ncode, odedef, xi, rtol, atol, itol, norm_p, laopt, algopt, remesh, xfix, nrmesh, dxmesh, trmesh, ipminf, monitf, rsave, isave, itask, itrace, ind, 'npts', npts, 'nxi', nxi, 'neqn', neqn, 'nxfix', nxfix, 'xratio', xratio, 'con', con)
Note: the interface to this routine has changed since earlier releases of the toolbox:
Mark 22: lrsave, lisave have been removed from the interface
.

Description

nag_pde_1d_parab_remesh_keller (d03pr) integrates the system of first-order PDEs and coupled ODEs given by the master equations:
Gi (x,t,U,Ux,Ut,V,V.) = 0 ,   i = 1,2,,npde ,   axb,tt0,
Gi (x,t,U,Ux,Ut,V,V.) = 0 ,   i=1,2,,npde ,   axb,tt0,
(1)
Ri(t,V,V.,ξ,U*,Ux * ,Ut * ) = 0,  i = 1,2,,ncode.
Ri(t,V,V.,ξ,U*,Ux*,Ut*)=0,  i=1,2,,ncode.
(2)
In the PDE part of the problem given by (1), the functions GiGi must have the general form
npde ncode
Gi = Pi,j(Uj)/(t) + Qi,jV.j + Si = 0,  i = 1,2,,npde,
j = 1 j = 1
Gi=j=1npdePi,j Uj t +j=1ncodeQi,jV.j+Si=0,  i=1,2,,npde,
(3)
where Pi,jPi,j, Qi,jQi,j and SiSi depend on xx, tt, UU, UxUx and VV.
The vector UU is the set of PDE solution values
U(x,t) = [U1(x,t),,Unpde(x,t)]T,
U(x,t)=[U1(x,t),,Unpde(x,t)]T,
and the vector UxUx is the partial derivative with respect to xx. The vector VV is the set of ODE solution values
V(t) = [V1(t),,Vncode(t)]T,
V(t)=[V1(t),,Vncode(t)]T,
and V.V. denotes its derivative with respect to time.
In the ODE part given by (2), ξξ represents a vector of nξ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*U*, Ux * Ux* and Ut * Ut* are the functions UU, UxUx and UtUt evaluated at these coupling points. Each RiRi may only depend linearly on time derivatives. Hence equation (2) may be written more precisely as
R = ABV.CUt * ,
R=A-BV.-CUt*,
(4)
where R = [R1,,Rncode]TR=[R1,,Rncode]T, AA is a vector of length ncode, BB is an ncode by ncode matrix, CC is an ncode by (nξ × npde)(nξ×npde) matrix and the entries in AA, BB and CC may depend on tt, ξξ, U*U*, Ux * Ux* and VV. In practice you only need to supply a vector of information to define the ODEs and not the matrices BB and CC. (See Section [Parameters] for the specification of odedef.)
The integration in time is from t0t0 to touttout, over the space interval axbaxb, where a = x1a=x1 and b = xnptsb=xnpts are the leftmost and rightmost points of a mesh x1,x2,,xnptsx1,x2,,xnpts defined initially by you and (possibly) adapted automatically during the integration according to user-specified criteria.
The PDE system which is defined by the functions GiGi must be specified in pdedef.
The initial (t = t0)(t=t0) values of the functions U(x,t)U(x,t) and V(t)V(t) must be specified in uvinit. Note that uvinit will be called again following any remeshing, and so U(x,t0)U(x,t0) should be specified for all values of xx in the interval axbaxb, and not just the initial mesh points.
For a first-order system of PDEs, only one boundary condition is required for each PDE component UiUi. The npde boundary conditions are separated into nana at the left-hand boundary x = ax=a, and nbnb at the right-hand boundary x = bx=b, such that na + nb = npdena+nb=npde. The position of the boundary condition for each component should be chosen with care; the general rule is that if the characteristic direction of UiUi at the left-hand boundary (say) points into the interior of the solution domain, then the boundary condition for UiUi should be specified at the left-hand boundary. Incorrect positioning of boundary conditions generally results in initialization or integration difficulties in the underlying time integration functions.
The boundary conditions have the master equation form:
GiL (x,t,U,Ut,V,V.) = 0   at ​ x = a ,   i = 1,2,,na ,
GiL (x,t,U,Ut,V,V.) = 0   at ​ x = a ,   i=1,2,,na ,
(5)
at the left-hand boundary, and
GiR (x,t,U,Ut,V,V.) = 0   at ​ x = b ,   i = 1,2,,nb ,
GiR (x,t,U,Ut,V,V.) = 0   at ​ x = b ,   i=1,2,,nb ,
(6)
at the right-hand boundary.
Note that the functions GiLGiL and GiRGiR must not depend on UxUx, since spatial derivatives are not determined explicitly in the Keller box scheme functions. If the problem involves derivative (Neumann) boundary conditions then it is generally possible to restate such boundary conditions in terms of permissible variables. Also note that GiLGiL and GiRGiR must be linear with respect to time derivatives, so that the boundary conditions have the general form:
npde ncode
Ei,jL(Uj)/(t) + Hi,jLV.j + KiL = 0,  i = 1,2,,na,
j = 1 j = 1
j=1 npde E i,j L Uj t + j=1 ncode H i,j L V.j + KiL = 0 ,   i=1,2,,na ,
(7)
at the left-hand boundary, and
npde ncode
Ei,jR(Uj)/(t) + Hi,jRV.j + KiR = 0,  i = 1,2,,nb,
j = 1 j = 1
j=1 npde E i,j R Uj t + j=1 ncode H i,j R V.j + KiR = 0 ,   i=1,2,,nb ,
(8)
at the right-hand boundary, where Ei,jLEi,jL, Ei,jREi,jR, Hi,jLHi,jL, Hi,jRHi,jR, KiLKiL and KiRKiR depend on x,t,Ux,t,U and VV only.
The boundary conditions must be specified in bndary.
The problem is subject to the following restrictions:
(i) Pi,jPi,j, Qi,jQi,j and SiSi must not depend on any time derivatives;
(ii) t0 < toutt0<tout, so that integration is in the forward direction;
(iii) The evaluation of the function GiGi is done approximately at the mid-points of the mesh x(i)xi, for i = 1,2,,nptsi=1,2,,npts, by calling pdedef for each mid-point in turn. Any discontinuities in the function must therefore be at one or more of the fixed mesh points specified by xfix;
(iv) At least one of the functions Pi,jPi,j must be nonzero so that there is a time derivative present in the PDE problem.
The algebraic-differential equation system which is defined by the functions RiRi must be specified in odedef. You must also specify the coupling points ξξ in the array xi.
The first-order equations are approximated by a system of ODEs in time for the values of UiUi at mesh points. In this method of lines approach the Keller box scheme is applied to each PDE in the space variable only, resulting in a system of ODEs in time for the values of UiUi at each mesh point. In total there are npde × npts + ncodenpde×npts+ncode 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 monitf which specifies in an analytic or numeric 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 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 along 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.

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
Keller H B (1970) A new difference scheme for parabolic problems Numerical Solutions of Partial Differential Equations (ed J Bramble) 2 327–350 Academic Press
Pennington S V and Berzins M (1994) New NAG Library software for first-order partial differential equations ACM Trans. Math. Softw. 20 63–99

Parameters

Compulsory Input Parameters

1:     npde – int64int32nag_int scalar
The number of PDEs to be solved.
Constraint: npde1npde1.
2:     ts – double scalar
The initial value of the independent variable tt.
Constraint: ts < toutts<tout.
3:     tout – double scalar
The final value of tt to which the integration is to be carried out.
4:     pdedef – function handle or string containing name of m-file
pdedef must evaluate the functions GiGi which define the system of PDEs. pdedef is called approximately midway between each pair of mesh points in turn by nag_pde_1d_parab_remesh_keller (d03pr).
[res, ires] = pdedef(npde, t, x, u, udot, ux, ncode, v, vdot, ires)

Input Parameters

1:     npde – int64int32nag_int scalar
The number of PDEs in the system.
2:     t – double scalar
The current value of the independent variable tt.
3:     x – double scalar
The current value of the space variable xx.
4:     u(npde) – double array
u(i)ui contains the value of the component Ui(x,t)Ui(x,t), for i = 1,2,,npdei=1,2,,npde.
5:     udot(npde) – double array
udot(i)udoti contains the value of the component (Ui(x,t))/(t) Ui(x,t) t , for i = 1,2,,npdei=1,2,,npde.
6:     ux(npde) – double array
ux(i)uxi contains the value of the component (Ui(x,t))/(x) Ui(x,t) x , for i = 1,2,,npdei=1,2,,npde.
7:     ncode – int64int32nag_int scalar
The number of coupled ODEs in the system.
8:     v(ncode) – double array
If ncode > 0ncode>0, v(i)vi contains the value of the component Vi(t)Vi(t), for i = 1,2,,ncodei=1,2,,ncode.
9:     vdot(ncode) – double array
If ncode > 0ncode>0, vdot(i)vdoti contains the value of component V.i(t)V.i(t), for i = 1,2,,ncodei=1,2,,ncode.
10:   ires – int64int32nag_int scalar
The form of GiGi that must be returned in the array res.
ires = 1ires=-1
Equation (9) must be used.
ires = 1ires=1
Equation (10) must be used.

Output Parameters

1:     res(npde) – double array
res(i)resi must contain the iith component of GG, for i = 1,2,,npdei=1,2,,npde, where GG is defined as
npde ncode
Gi = Pi,j(Uj)/(t) + Qi,jV.j,
j = 1 j = 1
Gi=j=1npdePi,j Uj t +j=1ncodeQi,jV.j,
(9)
i.e., only terms depending explicitly on time derivatives, or
npde ncode
Gi = Pi,j(Uj)/(t) + Qi,jV.j + Si,
j = 1 j = 1
Gi=j=1npdePi,j Uj t +j=1ncodeQi,jV.j+Si,
(10)
i.e., all terms in equation (3).
The definition of GG is determined by the input value of ires.
2:     ires – int64int32nag_int scalar
Should usually remain unchanged. However, you may set ires to force the integration function to take certain actions, as described below:
ires = 2ires=2
Indicates to the integrator that control should be passed back immediately to the calling (sub)routine with the error indicator set to ifail = 6ifail=6.
ires = 3ires=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 = 3ires=3 when a physically meaningless input or output value has been generated. If you consecutively set ires = 3ires=3, then nag_pde_1d_parab_remesh_keller (d03pr) returns to the calling function with the error indicator set to ifail = 4ifail=4.
5:     bndary – function handle or string containing name of m-file
bndary must evaluate the functions GiLGiL and GiRGiR which describe the boundary conditions, as given in (5) and (6).
[res, ires] = bndary(npde, t, ibnd, nobc, u, udot, ncode, v, vdot, ires)

Input Parameters

1:     npde – int64int32nag_int scalar
The number of PDEs in the system.
2:     t – double scalar
The current value of the independent variable tt.
3:     ibnd – int64int32nag_int scalar
Specifies which boundary conditions are to be evaluated.
ibnd = 0ibnd=0
bndary must compute the left-hand boundary condition at x = ax=a.
ibnd0ibnd0
bndary must compute of the right-hand boundary condition at x = bx=b.
4:     nobc – int64int32nag_int scalar
Specifies the number nana of boundary conditions at the boundary specified by ibnd.
5:     u(npde) – double array
u(i)ui contains the value of the component Ui(x,t)Ui(x,t) at the boundary specified by ibnd, for i = 1,2,,npdei=1,2,,npde.
6:     udot(npde) – double array
udot(i)udoti contains the value of the component (Ui(x,t))/(t) Ui(x,t) t , for i = 1,2,,npdei=1,2,,npde.
7:     ncode – int64int32nag_int scalar
The number of coupled ODEs in the system.
8:     v(ncode) – double array
If ncode > 0ncode>0, v(i)vi contains the value of the component Vi(t)Vi(t), for i = 1,2,,ncodei=1,2,,ncode.
9:     vdot(ncode) – double array
If ncode > 0ncode>0, vdot(i)vdoti contains the value of component V.i(t)V.i(t), for i = 1,2,,ncodei=1,2,,ncode.
Note: vdot(i)vdoti, for i = 1,2,,ncodei=1,2,,ncode, may only appear linearly as in (11) and (12).
10:   ires – int64int32nag_int scalar
The form of GiLGiL (or GiRGiR) that must be returned in the array res.
ires = 1ires=-1
Equation (11) must be used.
ires = 1ires=1
Equation (12) must be used.

Output Parameters

1:     res(nobc) – double array
res(i)resi must contain the iith component of GLGL or GRGR, depending on the value of ibnd, for i = 1,2,,nobci=1,2,,nobc, where GLGL is defined as
npde ncode
GiL = Ei,jL(Uj)/(t) + Hi,jLV.j,
j = 1 j = 1
GiL=j=1npdeEi,jL Uj t +j=1ncodeHi,jLV.j,
(11)
i.e., only terms depending explicitly on time derivatives, or
npde ncode
GiL = Ei,jL(Uj)/(t) + Hi,jLV.j + KiL,
j = 1 j = 1
GiL=j=1npdeEi,jL Uj t +j=1ncodeHi,jLV.j+KiL,
(12)
i.e., all terms in equation (7), and similarly for GiRGiR.
The definitions of GLGL and GRGR are determined by the input value of ires.
2:     ires – int64int32nag_int scalar
Should usually remain unchanged. However, you may set ires to force the integration function to take certain actions as described below:
ires = 2ires=2
Indicates to the integrator that control should be passed back immediately to the calling (sub)routine with the error indicator set to ifail = 6ifail=6.
ires = 3ires=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 = 3ires=3 when a physically meaningless input or output value has been generated. If you consecutively set ires = 3ires=3, then nag_pde_1d_parab_remesh_keller (d03pr) returns to the calling function with the error indicator set to ifail = 4ifail=4.
6:     uvinit – function handle or string containing name of m-file
uvinit must supply the initial (t = t0)(t=t0) values of U(x,t)U(x,t) and V(t)V(t) for all values of xx in the interval [a,b][a,b].
[u, v] = uvinit(npde, npts, nxi, x, xi, ncode)

Input Parameters

1:     npde – int64int32nag_int scalar
The number of PDEs in the system.
2:     npts – int64int32nag_int scalar
The number of mesh points in the interval [a,b][a,b].
3:     nxi – int64int32nag_int scalar
The number of ODE/PDE coupling points.
4:     x(npts) – double array
The current mesh. x(i)xi contains the value of xixi, for i = 1,2,,nptsi=1,2,,npts.
5:     xi(nxi) – double array
If nxi > 0nxi>0, xi(i)xii contains the ODE/PDE coupling point, ξiξi, for i = 1,2,,nxii=1,2,,nxi.
6:     ncode – int64int32nag_int scalar
The number of coupled ODEs in the system.

Output Parameters

1:     u(npde,npts) – double array
If nxi > 0nxi>0, u(i,j)uij contains the value of the component Ui(xj,t0)Ui(xj,t0), for i = 1,2,,npdei=1,2,,npde and j = 1,2,,nptsj=1,2,,npts.
2:     v(ncode) – double array
If ncode > 0ncode>0, v(i)vi must contain the value of component Vi(t0)Vi(t0), for i = 1,2,,ncodei=1,2,,ncode.
7:     u(neqn) – double array
neqn, the dimension of the array, must satisfy the constraint neqn = npde × npts + ncodeneqn=npde×npts+ncode.
If ind = 1ind=1, the value of u must be unchanged from the previous call.
8:     x(npts) – double array
npts, the dimension of the array, must satisfy the constraint npts3npts3.
The initial mesh points in the space direction. x(1)x1 must specify the left-hand boundary, aa, and x(npts)xnpts must specify the right-hand boundary, bb.
Constraint: x(1) < x(2) < < x(npts)x1<x2<<xnpts.
9:     nleft – int64int32nag_int scalar
The number nana of boundary conditions at the left-hand mesh point x(1)x1.
Constraint: 0nleftnpde0nleftnpde.
10:   ncode – int64int32nag_int scalar
The number of coupled ODE components.
Constraint: ncode0ncode0.
11:   odedef – function handle or string containing name of m-file
odedef must evaluate the functions RR, which define the system of ODEs, as given in (4).
If you wish to compute the solution of a system of PDEs only (i.e., ncode = 0ncode=0), odedef must be the string 'd03pek'. (nag_pde_1d_parab_dae_keller_remesh_fd_dummy_odedef (d03pek) is included in the NAG Toolbox.)
[r, ires] = odedef(npde, t, ncode, v, vdot, nxi, xi, ucp, ucpx, ucpt, ires)

Input Parameters

1:     npde – int64int32nag_int scalar
The number of PDEs in the system.
2:     t – double scalar
The current value of the independent variable tt.
3:     ncode – int64int32nag_int scalar
The number of coupled ODEs in the system.
4:     v(ncode) – double array
If ncode > 0ncode>0, v(i)vi contains the value of the component Vi(t)Vi(t), for i = 1,2,,ncodei=1,2,,ncode.
5:     vdot(ncode) – double array
If ncode > 0ncode>0, vdot(i)vdoti contains the value of component V.i(t)V.i(t), for i = 1,2,,ncodei=1,2,,ncode.
6:     nxi – int64int32nag_int scalar
The number of ODE/PDE coupling points.
7:     xi(nxi) – double array
If nxi > 0nxi>0, xi(i)xii contains the ODE/PDE coupling point, ξiξi, for i = 1,2,,nxii=1,2,,nxi.
8:     ucp(npde, : :) – double array
The second dimension of the array must be at least max (1,nxi)max(1,nxi)
If nxi > 0nxi>0, ucp(i,j)ucpij contains the value of Ui(x,t)Ui(x,t) at the coupling point x = ξjx=ξj, for i = 1,2,,npdei=1,2,,npde and j = 1,2,,nxij=1,2,,nxi.
9:     ucpx(npde, : :) – double array
The second dimension of the array must be at least max (1,nxi)max(1,nxi)
If nxi > 0nxi>0, ucpx(i,j)ucpxij contains the value of (Ui(x,t))/(x) Ui(x,t) x  at the coupling point x = ξjx=ξj, for i = 1,2,,npdei=1,2,,npde and j = 1,2,,nxij=1,2,,nxi.
10:   ucpt(npde, : :) – double array
The second dimension of the array must be at least max (1,nxi)max(1,nxi)
If nxi > 0nxi>0, ucpt(i,j)ucptij contains the value of (Ui)/(t) Ui t  at the coupling point x = ξjx=ξj, for i = 1,2,,npdei=1,2,,npde and j = 1,2,,nxij=1,2,,nxi.
11:   ires – int64int32nag_int scalar
The form of RR that must be returned in the array r.
ires = 1ires=-1
Equation (13) must be used.
ires = 1ires=1
Equation (14) must be used.

Output Parameters

1:     r(ncode) – double array
If ncode > 0ncode>0, r(i)ri must contain the iith component of RR, for i = 1,2,,ncodei=1,2,,ncode, where RR is defined as
R = BV.CUt * ,
R=-BV.-CUt*,
(13)
i.e., only terms depending explicitly on time derivatives, or
R = ABV.CUt * ,
R=A-BV.-CUt*,
(14)
i.e., all terms in equation (4). The definition of RR is determined by the input value of ires.
2:     ires – int64int32nag_int scalar
Should usually remain unchanged. However, you may reset ires to force the integration function to take certain actions, as described below:
ires = 2ires=2
Indicates to the integrator that control should be passed back immediately to the calling (sub)routine with the error indicator set to ifail = 6ifail=6.
ires = 3ires=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 = 3ires=3 when a physically meaningless input or output value has been generated. If you consecutively set ires = 3ires=3, then nag_pde_1d_parab_remesh_keller (d03pr) returns to the calling function with the error indicator set to ifail = 4ifail=4.
12:   xi( : :) – double array
Note: the dimension of the array xi must be at least max (1,nxi)max(1,nxi).
xi(i)xii, for i = 1,2,,nxii=1,2,,nxi, must be set to the ODE/PDE coupling points, ξiξi.
Constraint: x(1)xi(1) < xi(2) < < xi(nxi)x(npts)x1xi1<xi2<<xinxixnpts.
13:   rtol( : :) – double array
Note: the dimension of the array rtol must be at least 11 if itol = 1itol=1 or 22 and at least neqnneqn if itol = 3itol=3 or 44.
The relative local error tolerance.
Constraint: rtol(i)0.0rtoli0.0 for all relevant ii.
14:   atol( : :) – double array
Note: the dimension of the array atol must be at least 11 if itol = 1itol=1 or 33 and at least neqnneqn if itol = 2itol=2 or 44.
The absolute local error tolerance.
Constraint: atol(i)0.0atoli0.0 for all relevant ii.
Note: corresponding elements of rtol and atol cannot both be 0.00.0.
15:   itol – int64int32nag_int scalar
A value to indicate the form of the local error test. itol indicates to nag_pde_1d_parab_remesh_keller (d03pr) 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.0ei/wi<1.0, where wiwi is defined as follows:
itol rtol atol wiwi
1 scalar scalar rtol(1) × |u(i)| + atol(1)rtol1×|ui|+atol1
2 scalar vector rtol(1) × |u(i)| + atol(i)rtol1×|ui|+atoli
3 vector scalar rtol(i) × |u(i)| + atol(1)rtoli×|ui|+atol1
4 vector vector rtol(i) × |u(i)| + atol(i)rtoli×|ui|+atoli
In the above, eiei denotes the estimated local error for the iith component of the coupled PDE/ODE system in time, u(i)ui, for i = 1,2,,neqni=1,2,,neqn.
The choice of norm used is defined by the parameter norm_p.
Constraint: itol = 1itol=1, 22, 33 or 44.
16:   norm_p – string (length ≥ 1)
The type of norm to be used.
norm = 'M'norm='M'
Maximum norm.
norm = 'A'norm='A'
Averaged L2L2 norm.
If UnormUnorm denotes the norm of the vector u of length neqn, then for the averaged L2L2 norm
Unorm = sqrt( 1/(neqn) i = 1neqn (U(i) / wi)2 ),
Unorm= 1neqn i=1neqn (U(i)/wi) 2 ,
while for the maximum norm
Unorm = max |u(i) / wi|.
i
Unorm=maxi|ui/wi|.
See the description of itol for the formulation of the weight vector ww.
Constraint: norm = 'M'norm='M' or 'A''A'.
17:   laopt – string (length ≥ 1)
The type of matrix algebra required.
laopt = 'F'laopt='F'
Full matrix methods to be used.
laopt = 'B'laopt='B'
Banded matrix methods to be used.
laopt = 'S'laopt='S'
Sparse matrix methods to be used.
Constraint: laopt = 'F'laopt='F', 'B''B' or 'S''S'.
Note: you are recommended to use the banded option when no coupled ODEs are present (i.e., ncode = 0ncode=0).
18:   algopt(3030) – double array
May be set to control various options available in the integrator. If you wish to employ all the default options, then algopt(1)algopt1 should be set to 0.00.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(1)algopt1
Selects the ODE integration method to be used. If algopt(1) = 1.0algopt1=1.0, a BDF method is used and if algopt(1) = 2.0algopt1=2.0, a Theta method is used. The default value is algopt(1) = 1.0algopt1=1.0.
If algopt(1) = 2.0algopt1=2.0, then algopt(i)algopti, for i = 2,3,4i=2,3,4, are not used.
algopt(2)algopt2
Specifies the maximum order of the BDF integration formula to be used. algopt(2)algopt2 may be 1.01.0, 2.02.0, 3.03.0, 4.04.0 or 5.05.0. The default value is algopt(2) = 5.0algopt2=5.0.
algopt(3)algopt3
Specifies what method is to be used to solve the system of nonlinear equations arising on each step of the BDF method. If algopt(3) = 1.0algopt3=1.0 a modified Newton iteration is used and if algopt(3) = 2.0algopt3=2.0 a functional iteration method is used. If functional iteration is selected and the integrator encounters difficulty, then there is an automatic switch to the modified Newton iteration. The default value is algopt(3) = 1.0algopt3=1.0.
algopt(4)algopt4
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.0Pi,j=0.0, for j = 1,2,,npdej=1,2,,npde, for some ii or when there is no V.i(t)V.i(t) dependence in the coupled ODE system. If algopt(4) = 1.0algopt4=1.0, then the Petzold test is used. If algopt(4) = 2.0algopt4=2.0, then the Petzold test is not used. The default value is algopt(4) = 1.0algopt4=1.0.
If algopt(1) = 1.0algopt1=1.0, then algopt(i)algopti, for i = 5,6,7i=5,6,7, are not used.
algopt(5)algopt5
Specifies the value of Theta to be used in the Theta integration method. 0.51algopt(5)0.990.51algopt50.99. The default value is algopt(5) = 0.55algopt5=0.55.
algopt(6)algopt6
Specifies what method is to be used to solve the system of nonlinear equations arising on each step of the Theta method. If algopt(6) = 1.0algopt6=1.0, a modified Newton iteration is used and if algopt(6) = 2.0algopt6=2.0, a functional iteration method is used. The default value is algopt(6) = 1.0algopt6=1.0.
algopt(7)algopt7
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(7) = 1.0algopt7=1.0, then switching is allowed and if algopt(7) = 2.0algopt7=2.0, then switching is not allowed. The default value is algopt(7) = 1.0algopt7=1.0.
algopt(11)algopt11
Specifies a point in the time direction, tcrittcrit, beyond which integration must not be attempted. The use of tcrittcrit is described under the parameter itask. If algopt(1)0.0algopt10.0, a value of 0.00.0, for algopt(11)algopt11, say, should be specified even if itask subsequently specifies that tcrittcrit will not be used.
algopt(12)algopt12
Specifies the minimum absolute step size to be allowed in the time integration. If this option is not required, algopt(12)algopt12 should be set to 0.00.0.
algopt(13)algopt13
Specifies the maximum absolute step size to be allowed in the time integration. If this option is not required, algopt(13)algopt13 should be set to 0.00.0.
algopt(14)algopt14
Specifies the initial step size to be attempted by the integrator. If algopt(14) = 0.0algopt14=0.0, then the initial step size is calculated internally.
algopt(15)algopt15
Specifies the maximum number of steps to be attempted by the integrator in any one call. If algopt(15) = 0.0algopt15=0.0, then no limit is imposed.
algopt(23)algopt23
Specifies what method is to be used to solve the nonlinear equations at the initial point to initialize the values of UU, UtUt, VV and V.V.. If algopt(23) = 1.0algopt23=1.0, a modified Newton iteration is used and if algopt(23) = 2.0algopt23=2.0, functional iteration is used. The default value is algopt(23) = 1.0algopt23=1.0.
algopt(29)algopt29 and algopt(30)algopt30 are used only for the sparse matrix algebra option, i.e., laopt = 'S'laopt='S'.
algopt(29)algopt29
Governs the choice of pivots during the decomposition of the first Jacobian matrix. It should lie in the range 0.0 < algopt(29) < 1.00.0<algopt29<1.0, with smaller values biasing the algorithm towards maintaining sparsity at the expense of numerical stability. If algopt(29)algopt29 lies outside this range then the default value is used. If the functions regard the Jacobian matrix as numerically singular then increasing algopt(29)algopt29 towards 1.01.0 may help, but at the cost of increased fill-in. The default value is algopt(29) = 0.1algopt29=0.1.
algopt(30)algopt30
Used as a relative pivot threshold during subsequent Jacobian decompositions (see algopt(29)algopt29) below which an internal error is invoked. algopt(30)algopt30 must be greater than zero, otherwise the default value is used. If algopt(30)algopt30 is greater than 1.01.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(29)algopt29). The default value is algopt(30) = 0.0001algopt30=0.0001.
19:   remesh – logical scalar
Indicates whether or not spatial remeshing should be performed.
remesh = trueremesh=true
Indicates that spatial remeshing should be performed as specified.
remesh = falseremesh=false
Indicates that spatial remeshing should be suppressed.
Note:  remesh should not be changed between consecutive calls to nag_pde_1d_parab_remesh_keller (d03pr). Remeshing can be switched off or on at specified times by using appropriate values for the parameters nrmesh and trmesh at each call.
20:   xfix( : :) – double array
Note: the dimension of the array xfix must be at least max (1,nxfix)max(1,nxfix).
xfix(i)xfixi, for i = 1,2,,nxfixi=1,2,,nxfix, must contain the value of the xx coordinate at the iith fixed mesh point.
Constraint: xfix(i) < xfix(i + 1)xfixi<xfixi+1, for i = 1,2,,nxfix1i=1,2,,nxfix-1, and each fixed mesh point must coincide with a user-supplied initial mesh point, that is xfix(i) = x(j)xfixi=xj for some jj, 2jnpts12jnpts-1.
Note: the positions of the fixed mesh points in the array xx 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.
21:   nrmesh – int64int32nag_int scalar
Indicates the form of meshing to be performed.
nrmesh < 0nrmesh<0
Indicates that a new mesh is adopted according to the parameter dxmesh. The mesh is tested every |nrmesh||nrmesh| timesteps.
nrmesh = 0nrmesh=0
Indicates that remeshing should take place just once at the end of the first time step reached when t > trmesht>trmesh.
nrmesh > 0nrmesh>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_1d_parab_remesh_keller (d03pr) to give greater flexibility over the times of remeshing.
22:   dxmesh – double scalar
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||nrmesh| time steps, but is adopted only if
xinew > xiold + dxmesh × (xi + 1oldxiold),
xinew>xiold+dxmesh×(xi+1old-xiold),
or
xinew < xiolddxmesh × (xioldxi 1old).
xinew<xiold-dxmesh×(xiold-xi- 1old).
dxmesh thus imposes a lower limit on the difference between one mesh and the next.
Constraint: dxmesh0.0dxmesh0.0.
23:   trmesh – double scalar
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 tt is greater than trmesh.
Note: trmesh may be changed between consecutive calls to nag_pde_1d_parab_remesh_keller (d03pr) to force remeshing at several specified times.
24:   ipminf – int64int32nag_int scalar
The level of trace information regarding the adaptive remeshing. Details are directed to the current advisory message unit (see nag_file_set_unit_advisory (x04ab)).
ipminf = 0ipminf=0
No trace information.
ipminf = 1ipminf=1
Brief summary of mesh characteristics.
ipminf = 2ipminf=2
More detailed information, including old and new mesh points, mesh sizes and monitor function values.
Constraint: ipminf = 0ipminf=0, 11 or 22.
25:   monitf – function handle or string containing name of m-file
monitf must supply and evaluate a remesh monitor function to indicate the solution behaviour of interest.
If you specify remesh = falseremesh=false, i.e., no remeshing, then monitf will not be called and the string 'd03pel' may be used for monitf. (nag_pde_1d_parab_dae_keller_remesh_fd_dummy_monitf (d03pel) is included in the NAG Toolbox.)
[fmon] = monitf(t, npts, npde, x, u)

Input Parameters

1:     t – double scalar
The current value of the independent variable tt.
2:     npts – int64int32nag_int scalar
The number of mesh points in the interval [a,b][a,b].
3:     npde – int64int32nag_int scalar
The number of PDEs in the system.
4:     x(npts) – double array
The current mesh. x(i)xi contains the value of xixi, for i = 1,2,,nptsi=1,2,,npts.
5:     u(npde,npts) – double array
The second dimension of the array must be at least npde × nptsnpde×npts
u(i,j)uij contains the value of Ui(x,t)Ui(x,t) at x = x(j)x=xj and time tt, for i = 1,2,,npdei=1,2,,npde and j = 1,2,,nptsj=1,2,,npts.

Output Parameters

1:     fmon(npts) – double array
fmon(i)fmoni must contain the value of the monitor function Fmon(x)Fmon(x) at mesh point x = x(i)x=xi.
26:   rsave(lrsave) – double array
If ind = 0ind=0, rsave need not be set on entry.
If ind = 1ind=1, rsave must be unchanged from the previous call to the function because it contains required information about the iteration.
27:   isave(lisave) – int64int32nag_int array
If ind = 0ind=0, isave need not be set.
If ind = 1ind=1, isave must be unchanged from the previous call to the function because it contains required information about the iteration. In particular the following components of the array isave concern the efficiency of the integration:
isave(1)isave1
Contains the number of steps taken in time.
isave(2)isave2
Contains the number of residual evaluations of the resulting ODE system used. One such evaluation involves evaluating the PDE functions at all the mesh points, as well as one evaluation of the functions in the boundary conditions.
isave(3)isave3
Contains the number of Jacobian evaluations performed by the time integrator.
isave(4)isave4
Contains the order of the ODE method last used in the time integration.
isave(5)isave5
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 LULU decomposition of the Jacobian matrix.
The rest of the array is used as workspace.
28:   itask – int64int32nag_int scalar
The task to be performed by the ODE integrator.
itask = 1itask=1
Normal computation of output values u at t = toutt=tout (by overshooting and interpolating).
itask = 2itask=2
Take one step in the time direction and return.
itask = 3itask=3
Stop at first internal integration point at or beyond t = toutt=tout.
itask = 4itask=4
Normal computation of output values u at t = toutt=tout but without overshooting t = tcritt=tcrit where tcrittcrit is described under the parameter algopt.
itask = 5itask=5
Take one step in the time direction and return, without passing tcrittcrit, where tcrittcrit is described under the parameter algopt.
Constraint: itask = 1itask=1, 22, 33, 44 or 55.
29:   itrace – int64int32nag_int scalar
The level of trace information required from nag_pde_1d_parab_remesh_keller (d03pr) and the underlying ODE solver as follows:
itrace1itrace-1
No output is generated.
itrace = 0itrace=0
Only warning messages from the PDE solver are printed on the current error message unit (see nag_file_set_unit_error (x04aa)).
itrace = 1itrace=1
Output from the underlying ODE solver is printed on the current advisory message unit (see nag_file_set_unit_advisory (x04ab)). This output contains details of Jacobian entries, the nonlinear iteration and the time integration during the computation of the ODE system.
itrace = 2itrace=2
Output from the underlying ODE solver is similar to that produced when itrace = 1itrace=1, except that the advisory messages are given in greater detail.
itrace3itrace3
The output from the underlying ODE solver is similar to that produced when itrace = 2itrace=2, except that the advisory messages are given in greater detail.
You are advised to set itrace = 0itrace=0, unless you are experienced with sub-chapter D02M–N.
30:   ind – int64int32nag_int scalar
Indicates whether this is a continuation call or a new integration.
ind = 0ind=0
Starts or restarts the integration in time.
ind = 1ind=1
Continues the integration after an earlier exit from the function. In this case, only the parameters tout and ifail and the remeshing parameters nrmesh, dxmesh, trmesh, xratio and con may be reset between calls to nag_pde_1d_parab_remesh_keller (d03pr).
Constraint: ind = 0ind=0 or 11.

Optional Input Parameters

1:     npts – int64int32nag_int scalar
Default: The dimension of the array x.
The number of mesh points in the interval [a,ba,b].
Constraint: npts3npts3.
2:     nxi – int64int32nag_int scalar
Default: The dimension of the array xi.
The number of ODE/PDE coupling points.
Constraints:
3:     neqn – int64int32nag_int scalar
Default: The dimension of the array u.
The number of ODEs in the time direction.
Constraint: neqn = npde × npts + ncodeneqn=npde×npts+ncode.
4:     nxfix – int64int32nag_int scalar
Default: The dimension of the array xfix.
The number of fixed mesh points.
Constraint: 0nxfixnpts20nxfixnpts-2.
Note: the end points x(1)x1 and x(npts)xnpts are fixed automatically and hence should not be specified as fixed points.
5:     xratio – double scalar
Input bound on adjacent mesh ratio (greater than 1.01.0 and typically in the range 1.51.5 to 3.03.0). The remeshing functions will attempt to ensure that
(xixi1) / xratio < xi + 1xi < xratio × (xixi1).
(xi-xi-1)/xratio<xi+1-xi<xratio×(xi-xi-1).
Default: 1.51.5
Constraint: xratio > 1.0xratio>1.0.
6:     con – double scalar
An input bound on the sub-integral of the monitor function Fmon(x)Fmon(x) over each space step. The remeshing functions will attempt to ensure that
xi + 1 xnpts
Fmon(x)dxconFmon(x)dx,
x1 x1
x1xi+1Fmon(x)dxconx1xnptsFmon(x)dx,
(see Furzeland (1984)). con gives you more control over the mesh distribution e.g., decreasing con allows more clustering. A typical value is 2 / (npts1)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.
Default: 2.0 / (npts1)2.0/(npts-1)
Constraint: 0.1 / (npts1)con10.0 / (npts1)0.1/(npts-1)con10.0/(npts-1).

Input Parameters Omitted from the MATLAB Interface

lrsave lisave

Output Parameters

1:     ts – double scalar
The value of tt corresponding to the solution values in u. Normally ts = toutts=tout.
2:     u(neqn) – double array
u(npde × (j1) + i)unpde×(j-1)+i contains the computed solution Ui(xj,t)Ui(xj,t), for i = 1,2,,npdei=1,2,,npde and j = 1,2,,nptsj=1,2,,npts, evaluated at t = tst=ts.
3:     x(npts) – double array
The final values of the mesh points.
4:     rsave(lrsave) – double array
If ind = 1ind=1, rsave must be unchanged from the previous call to the function because it contains required information about the iteration.
5:     isave(lisave) – int64int32nag_int array
If ind = 1ind=1, isave must be unchanged from the previous call to the function because it contains required information about the iteration. In particular the following components of the array isave concern the efficiency of the integration:
isave(1)isave1
Contains the number of steps taken in time.
isave(2)isave2
Contains the number of residual evaluations of the resulting ODE system used. One such evaluation involves evaluating the PDE functions at all the mesh points, as well as one evaluation of the functions in the boundary conditions.
isave(3)isave3
Contains the number of Jacobian evaluations performed by the time integrator.
isave(4)isave4
Contains the order of the ODE method last used in the time integration.
isave(5)isave5
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 LULU decomposition of the Jacobian matrix.
The rest of the array is used as workspace.
6:     ind – int64int32nag_int scalar
ind = 1ind=1.
7:     ifail – int64int32nag_int scalar
ifail = 0ifail=0 unless the function detects an error (see [Error Indicators and Warnings]).

Error Indicators and Warnings

Errors or warnings detected by the function:

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

  ifail = 1ifail=1
On entry,(toutts)(tout-ts) is too small,
oritask1itask1, 22, 33, 44 or 55,
orat least one of the coupling points defined in array xi is outside the interval [x(1),x(npts)x1,xnpts],
ornpts < 3npts<3,
ornpde < 1npde<1,
ornleft not in the range 00 to npde,
ornorm'A'norm'A' or 'M''M',
orlaopt'F'laopt'F', 'B''B' or 'S''S',
oritol1itol1, 22, 33 or 44,
orind0ind0 or 11,
ormesh points x(i)xi are badly ordered,
orlrsave is too small,
orlisave is too small,
orncode and nxi are incorrectly defined,
orind = 1ind=1 on initial entry to nag_pde_1d_parab_remesh_keller (d03pr),
oran element of rtol or atol < 0.0atol<0.0,
orcorresponding elements of rtol and atol are both 0.00.0,
orneqnnpde × npts + ncodeneqnnpde×npts+ncode,
ornxfix not in the range 00 to npts2npts-2,
orfixed mesh point(s) do not coincide with any of the user-supplied mesh points,
ordxmesh < 0.0dxmesh<0.0,
oripminf0ipminf0, 11 or 22,
orxratio1.0xratio1.0,
orcon not in the range 0.1 / (npts1)0.1/(npts-1) to 10 / (npts1)10/(npts-1).
W ifail = 2ifail=2
The underlying ODE solver cannot make any further progress, with the values of atol and rtol, across the integration range from the current point t = tst=ts. The components of u contain the computed values at the current point t = tst=ts.
W ifail = 3ifail=3
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 = tst=ts. The problem may have a singularity, or the error requirement may be inappropriate. Incorrect positioning of boundary conditions may also result in this error.
  ifail = 4ifail=4
In setting up the ODE system, the internal initialization function was unable to initialize the derivative of the ODE system. This could be due to the fact that ires was repeatedly set to 33 in one of pdedef, bndary or odedef, when the residual in the underlying ODE solver was being evaluated. Incorrect positioning of boundary conditions may also result in this error.
  ifail = 5ifail=5
In solving the ODE system, a singular Jacobian has been encountered. You should check their problem formulation.
W ifail = 6ifail=6
When evaluating the residual in solving the ODE system, ires was set to 22 in one of pdedef, bndary or odedef. Integration was successful as far as t = tst=ts.
  ifail = 7ifail=7
The values of atol and rtol are so small that the function is unable to start the integration in time.
  ifail = 8ifail=8
In either, pdedef, bndary or odedef, ires was set to an invalid value.
  ifail = 9ifail=9 (nag_ode_ivp_stiff_imp_revcom (d02nn))
A serious error has occurred in an internal call to the specified function. Check the problem specification an all parameters and array dimensions. Setting itrace = 1itrace=1 may provide more information. If the problem persists, contact NAG.
W ifail = 10ifail=10
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 itask2itask2 or 55.)
  ifail = 11ifail=11
An error occurred during Jacobian formulation of the ODE system (a more detailed error description may be directed to the current advisory message unit). If using the sparse matrix algebra option, the values of algopt(29)algopt29 and algopt(30)algopt30 may be inappropriate.
  ifail = 12ifail=12
In solving the ODE system, the maximum number of steps specified in algopt(15)algopt15 have been taken.
W ifail = 13ifail=13
Some error weights wiwi became zero during the time integration (see the description of itol). Pure relative error control (atol(i) = 0.0)(atoli=0.0) was requested on a variable (the iith) which has become zero. The integration was successful as far as t = tst=ts.
  ifail = 14ifail=14
Not applicable.
  ifail = 15ifail=15
When using the sparse option, the value of lisave or lrsave was insufficient (more detailed information may be directed to the current error message unit).
  ifail = 16ifail=16
remesh has been changed between calls to nag_pde_1d_parab_remesh_keller (d03pr).

Accuracy

nag_pde_1d_parab_remesh_keller (d03pr) 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 parameters, atol and rtol.

Further Comments

The Keller box scheme can be used to solve higher-order problems which have been reduced to first-order by the introduction of new variables (see the example in Section [Example]). In general, a second-order problem can be solved with slightly greater accuracy using the Keller box scheme instead of a finite difference scheme (nag_pde_1d_parab_remesh_fd (d03pp) for example), but at the expense of increased CPU time due to the larger number of function evaluations required.
It should be noted that the Keller box scheme, in common with other central-difference schemes, may be unsuitable for some hyperbolic first-order problems such as the apparently simple linear advection equation Ut + aUx = 0Ut+aUx=0, where aa is a constant, resulting in spurious oscillations due to the lack of dissipation. This type of problem requires a discretization scheme with upwind weighting (nag_pde_1d_parab_convdiff_remesh (d03ps) for example), or the addition of a second-order artificial dissipation term.
The time taken depends on the complexity of the 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.

Example

function nag_pde_1d_parab_remesh_keller_example
npde = int64(2);
ts = 0;
tout = 0.05;
u = zeros(122, 1);
x = zeros(61,1);
for i=2:61
 x(i) = (i-1)/60;
end
nleft = int64(1);
ncode = int64(0);
xi = [];
rtol = [5e-05];
atol = [5e-05];
itol = int64(1);
normtype = 'A';
laopt = 'F';
algopt = zeros(30, 1);
remesh = true;
xfix = [];
nrmesh = int64(3);
dxmesh = 0;
trmesh = 0;
ipminf = int64(0);
rsave = zeros(17004, 1);
isave = zeros(25, 1, 'int64');
itask = int64(1);
itrace = int64(0);
ind = int64(0);
[tsOut, uOut, xOut, rsaveOut, isaveOut, indOut, ifail] = ...
    nag_pde_1d_parab_remesh_keller(npde, ts, tout, @pdedef, @bndary, @uvinit, ...
    u, x, nleft, ncode, @odedef, xi, rtol, atol, itol, normtype, ...
    laopt, algopt, remesh, xfix, nrmesh, dxmesh, trmesh, ipminf, ...
    @monitf, rsave, isave, itask, itrace, ind, 'xratio', 1.2, 'con', 5/60);
 tsOut
 uOut
 xOut
 indOut
 ifail

function [res, ires] = pdedef(npde, t, x, u, udot, ux, ncode, v, vdot, ires)
  res = zeros(npde, 1);

  if (ires == -1)
    res(1) = udot(1);
    res(2) = udot(2);
  else
    res(1) = udot(1) + ux(1) + ux(2);
    res(2) = udot(2) + 4.0d0*ux(1) + ux(2);
  end

function [res, ires] = bndary(npde, t, ibnd, nobc, u, udot, ncode, v, vdot, ires)
  res = zeros(nobc, 1);

  pp = 2*pi;
  if (ibnd == 0)
    if (ires ~= -1)
      res(1) = u(1) - ...
               0.5*(exp(t)+exp(-3*t)) - 0.25*(sin(pp*9*t^2)-sin(pp*t^2)) - 2*t^2;
    end
  else
    if (ires ~= -1)
      res(1) = u(2) - ...
               (exp(1-3*t)-exp(1+t) +0.5*(sin(pp*(1-3*t)^2)+sin(pp*(1+t) ^2))+1+5*t^2-2*t);
    end
  end

function [u, v] = uvinit(npde, npts, nxi, x, xi, ncode)
  u = zeros(npde, npts);
  v = zeros(ncode, 1);

  for i = 1:double(npts)
    u(1,i) = exp(x(i));
    u(2,i) = x(i)^2 + sin(2*pi*x(i)^2);
  end

function [f, ires] = odedef(npde, t, ncode, v, vdot, nxi, xi, ucp, ucpx, ucpt, ires)
  f = zeros(ncode, 1);
function [fmon] = monitf(t, npts, npde, x, u)
  fmon = zeros(npts, 1);

  for i = 2:double(npts) - 1
    h1 = x(i) - x(i-1);
    h2 = x(i+1) - x(i);
    h3 = 0.5d0*(x(i+1)-x(i-1));
    % second derivatives ..
    d2x1 = abs(((u(1,i+1)-u(1,i))/h2-(u(1,i)-u(1,i-1))/h1)/h3);
    d2x2 = abs(((u(2,i+1)-u(2,i))/h2-(u(2,i)-u(2,i-1))/h1)/h3);
    fmon(i) = max(d2x1,d2x2);
  end
  fmon(1) = fmon(2);
  fmon(npts) = fmon(double(npts)-1);
 

tsOut =

    0.0500


uOut =

    0.9923
   -0.0997
    1.0006
   -0.1223
    1.0100
   -0.1288
    1.0206
   -0.1191
    1.0327
   -0.0933
    1.0460
   -0.0532
    1.0606
   -0.0017
    1.0767
    0.0586
    1.0948
    0.1255
    1.1150
    0.1964
    1.1379
    0.2691
    1.1636
    0.3414
    1.1924
    0.4112
    1.2245
    0.4769
    1.2599
    0.5371
    1.2987
    0.5906
    1.3405
    0.6365
    1.3852
    0.6741
    1.4325
    0.7031
    1.4819
    0.7235
    1.5330
    0.7353
    1.5853
    0.7389
    1.6384
    0.7347
    1.6920
    0.7233
    1.7467
    0.7049
    1.8031
    0.6793
    1.8620
    0.6456
    1.9239
    0.6030
    1.9888
    0.5505
    2.0549
    0.4885
    2.1187
    0.4204
    2.1771
    0.3498
    2.2295
    0.2790
    2.2754
    0.2096
    2.3147
    0.1429
    2.3474
    0.0802
    2.3737
    0.0223
    2.3939
   -0.0301
    2.4084
   -0.0767
    2.4176
   -0.1172
    2.4220
   -0.1516
    2.4222
   -0.1799
    2.4187
   -0.2024
    2.4120
   -0.2193
    2.4027
   -0.2310
    2.3911
   -0.2380
    2.3776
   -0.2405
    2.3623
   -0.2388
    2.3454
   -0.2329
    2.3268
   -0.2228
    2.3070
   -0.2085
    2.2857
   -0.1896
    2.2634
   -0.1660
    2.2406
   -0.1379
    2.2169
   -0.1041
    2.1932
   -0.0649
    2.1700
   -0.0195
    2.1481
    0.0323
    2.1295
    0.0882
    2.1157
    0.1450
    2.1071
    0.2021


xOut =

         0
    0.0337
    0.0674
    0.1013
    0.1351
    0.1680
    0.1995
    0.2297
    0.2586
    0.2862
    0.3127
    0.3380
    0.3621
    0.3853
    0.4074
    0.4286
    0.4489
    0.4684
    0.4870
    0.5048
    0.5219
    0.5382
    0.5538
    0.5689
    0.5837
    0.5986
    0.6137
    0.6294
    0.6457
    0.6627
    0.6795
    0.6956
    0.7111
    0.7258
    0.7400
    0.7536
    0.7665
    0.7789
    0.7908
    0.8022
    0.8131
    0.8235
    0.8335
    0.8430
    0.8522
    0.8610
    0.8695
    0.8780
    0.8864
    0.8948
    0.9032
    0.9118
    0.9206
    0.9294
    0.9387
    0.9483
    0.9584
    0.9689
    0.9796
    0.9899
    1.0000


indOut =

                    1


ifail =

                    0


function d03pr_example
npde = int64(2);
ts = 0;
tout = 0.05;
u = zeros(122, 1);
x = zeros(61,1);
for i=2:61
 x(i) = (i-1)/60;
end
nleft = int64(1);
ncode = int64(0);
xi = [];
rtol = [5e-05];
atol = [5e-05];
itol = int64(1);
normtype = 'A';
laopt = 'F';
algopt = zeros(30, 1);
remesh = true;
xfix = [];
nrmesh = int64(3);
dxmesh = 0;
trmesh = 0;
ipminf = int64(0);
rsave = zeros(17004, 1);
isave = zeros(25, 1, 'int64');
itask = int64(1);
itrace = int64(0);
ind = int64(0);
[tsOut, uOut, xOut, rsaveOut, isaveOut, indOut, ifail] = ...
    d03pr(npde, ts, tout, @pdedef, @bndary, @uvinit, ...
    u, x, nleft, ncode, @odedef, xi, rtol, atol, itol, normtype, ...
    laopt, algopt, remesh, xfix, nrmesh, dxmesh, trmesh, ipminf, ...
    @monitf, rsave, isave, itask, itrace, ind, 'xratio', 1.2, 'con', 5/60);
 tsOut
 uOut
 xOut
 indOut
 ifail

function [res, ires] = pdedef(npde, t, x, u, udot, ux, ncode, v, vdot, ires)
  res = zeros(npde, 1);

  if (ires == -1)
    res(1) = udot(1);
    res(2) = udot(2);
  else
    res(1) = udot(1) + ux(1) + ux(2);
    res(2) = udot(2) + 4.0d0*ux(1) + ux(2);
  end

function [res, ires] = bndary(npde, t, ibnd, nobc, u, udot, ncode, v, vdot, ires)
  res = zeros(nobc, 1);

  pp = 2*pi;
  if (ibnd == 0)
    if (ires ~= -1)
      res(1) = u(1) - ...
               0.5*(exp(t)+exp(-3*t)) - 0.25*(sin(pp*9*t^2)-sin(pp*t^2)) - 2*t^2;
    end
  else
    if (ires ~= -1)
      res(1) = u(2) - ...
               (exp(1-3*t)-exp(1+t) +0.5*(sin(pp*(1-3*t)^2)+sin(pp*(1+t) ^2))+1+5*t^2-2*t);
    end
  end

function [u, v] = uvinit(npde, npts, nxi, x, xi, ncode)
  u = zeros(npde, npts);
  v = zeros(ncode, 1);

  for i = 1:double(npts)
    u(1,i) = exp(x(i));
    u(2,i) = x(i)^2 + sin(2*pi*x(i)^2);
  end

function [f, ires] = odedef(npde, t, ncode, v, vdot, nxi, xi, ucp, ucpx, ucpt, ires)
  f = zeros(ncode, 1);
function [fmon] = monitf(t, npts, npde, x, u)
  fmon = zeros(npts, 1);

  for i = 2:double(npts) - 1
    h1 = x(i) - x(i-1);
    h2 = x(i+1) - x(i);
    h3 = 0.5d0*(x(i+1)-x(i-1));
    % second derivatives ..
    d2x1 = abs(((u(1,i+1)-u(1,i))/h2-(u(1,i)-u(1,i-1))/h1)/h3);
    d2x2 = abs(((u(2,i+1)-u(2,i))/h2-(u(2,i)-u(2,i-1))/h1)/h3);
    fmon(i) = max(d2x1,d2x2);
  end
  fmon(1) = fmon(2);
  fmon(npts) = fmon(double(npts)-1);
 

tsOut =

    0.0500


uOut =

    0.9923
   -0.0997
    1.0006
   -0.1223
    1.0100
   -0.1288
    1.0206
   -0.1191
    1.0327
   -0.0933
    1.0460
   -0.0532
    1.0606
   -0.0017
    1.0767
    0.0586
    1.0948
    0.1255
    1.1150
    0.1964
    1.1379
    0.2691
    1.1636
    0.3414
    1.1924
    0.4112
    1.2245
    0.4769
    1.2599
    0.5371
    1.2987
    0.5906
    1.3405
    0.6365
    1.3852
    0.6741
    1.4325
    0.7031
    1.4819
    0.7235
    1.5330
    0.7353
    1.5853
    0.7389
    1.6384
    0.7347
    1.6920
    0.7233
    1.7467
    0.7049
    1.8031
    0.6793
    1.8620
    0.6456
    1.9239
    0.6030
    1.9888
    0.5505
    2.0549
    0.4885
    2.1187
    0.4204
    2.1771
    0.3498
    2.2295
    0.2790
    2.2754
    0.2096
    2.3147
    0.1429
    2.3474
    0.0802
    2.3737
    0.0223
    2.3939
   -0.0301
    2.4084
   -0.0767
    2.4176
   -0.1172
    2.4220
   -0.1516
    2.4222
   -0.1799
    2.4187
   -0.2024
    2.4120
   -0.2193
    2.4027
   -0.2310
    2.3911
   -0.2380
    2.3776
   -0.2405
    2.3623
   -0.2388
    2.3454
   -0.2329
    2.3268
   -0.2228
    2.3070
   -0.2085
    2.2857
   -0.1896
    2.2634
   -0.1660
    2.2406
   -0.1379
    2.2169
   -0.1041
    2.1932
   -0.0649
    2.1700
   -0.0195
    2.1481
    0.0323
    2.1295
    0.0882
    2.1157
    0.1450
    2.1071
    0.2021


xOut =

         0
    0.0337
    0.0674
    0.1013
    0.1351
    0.1680
    0.1995
    0.2297
    0.2586
    0.2862
    0.3127
    0.3380
    0.3621
    0.3853
    0.4074
    0.4286
    0.4489
    0.4684
    0.4870
    0.5048
    0.5219
    0.5382
    0.5538
    0.5689
    0.5837
    0.5986
    0.6137
    0.6294
    0.6457
    0.6627
    0.6795
    0.6956
    0.7111
    0.7258
    0.7400
    0.7536
    0.7665
    0.7789
    0.7908
    0.8022
    0.8131
    0.8235
    0.8335
    0.8430
    0.8522
    0.8610
    0.8695
    0.8780
    0.8864
    0.8948
    0.9032
    0.9118
    0.9206
    0.9294
    0.9387
    0.9483
    0.9584
    0.9689
    0.9796
    0.9899
    1.0000


indOut =

                    1


ifail =

                    0



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

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