nag_pde_parab_1d_cd_ode_remesh (d03psc) (PDF version)
d03 Chapter Contents
d03 Chapter Introduction
NAG Library Manual

NAG Library Function Document

nag_pde_parab_1d_cd_ode_remesh (d03psc)

+ Contents

    1  Purpose
    7  Accuracy

1  Purpose

nag_pde_parab_1d_cd_ode_remesh (d03psc) integrates a system of linear or nonlinear convection-diffusion equations in one space dimension, with optional source terms and scope for coupled ordinary differential equations (ODEs). The system must be posed in conservative form. This function also includes the option of automatic adaptive spatial remeshing. Convection terms are discretized using a sophisticated upwind scheme involving a user-supplied numerical flux function based on the solution of a Riemann problem at each mesh point. The method of lines is employed to reduce the partial differential equations (PDEs) to a system of ODEs, and the resulting system is solved using a backward differentiation formula (BDF) method or a Theta method.

2  Specification

#include <nag.h>
#include <nagd03.h>
void  nag_pde_parab_1d_cd_ode_remesh (Integer npde, double *ts, double tout,
void (*pdedef)(Integer npde, double t, double x, const double u[], const double ux[], Integer ncode, const double v[], const double vdot[], double p[], double c[], double d[], double s[], Integer *ires, Nag_Comm *comm),
void (*numflx)(Integer npde, double t, double x, Integer ncode, const double v[], const double uleft[], const double uright[], double flux[], Integer *ires, Nag_Comm *comm, Nag_D03_Save *saved),
void (*bndary)(Integer npde, Integer npts, double t, const double x[], const double u[], Integer ncode, const double v[], const double vdot[], Integer ibnd, double g[], Integer *ires, Nag_Comm *comm),
void (*uvinit)(Integer npde, Integer npts, Integer nxi, const double x[], const double xi[], double u[], Integer ncode, double v[], Nag_Comm *comm),
double u[], Integer npts, double x[], Integer ncode,
void (*odedef)(Integer npde, double t, Integer ncode, const double v[], const double vdot[], Integer nxi, const double xi[], const double ucp[], const double ucpx[], const double ucpt[], double r[], 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[], 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_cd_ode_remesh (d03psc) integrates the system of convection-diffusion equations in conservative form:
j=1npdePi,j Uj t + Fi x =Ci Di x +Si, (1)
or the hyperbolic convection-only system:
Ui t + Fi x =0, (2)
for i=1,2,,npde, axb, tt0, where the vector U is the set of PDE solution values
Ux,t=U1x,t,,Unpdex,tT.
The optional coupled ODEs are of the general form
Rit,V,V.,ξ,U*,Ux*,Ut*=0,  i=1,2,,ncode, (3)
where the vector V is the set of ODE solution values
Vt=V1t,,VncodetT,
V. denotes its derivative with respect to time, and Ux is the spatial derivative of U.
In (2), Pi,j, Fi and Ci depend on x, t, U and V; Di depends on x, t, U, Ux and V; and Si depends on x, t, U, V and linearly on V.. Note that Pi,j, Fi, Ci and Si must not depend on any space derivatives, and Pi,j, Fi, Ci and Di must not depend on any time derivatives. In terms of conservation laws, Fi, CiDi x  and Si are the convective flux, diffusion and source terms respectively.
In (3), ξ 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 PDE spatial mesh points. U*, Ux* and Ut* are the functions U, Ux and Ut evaluated at these coupling points. Each Ri may depend only linearly on time derivatives. Hence (3) may be written more precisely as
R=L-MV.-NUt*, (4)
where R=R1,,RncodeT, L is a vector of length ncode, M is an ncode by ncode matrix, N is an ncode by nξ×npde matrix and the entries in L, M and N 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 L, M and N. (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 user-defined mesh x1,x2,,xnpts defined initially by you and (possibly) adapted automatically during the integration according to user-specified criteria.
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 PDEs are approximated by a system of ODEs in time for the values of Ui at mesh points using a spatial discretization method similar to the central-difference scheme used in nag_pde_parab_1d_fd (d03pcc)nag_pde_parab_1d_fd_ode (d03phc) and nag_pde_parab_1d_fd_ode_remesh (d03ppc), but with the flux Fi replaced by a numerical flux, which is a representation of the flux taking into account the direction of the flow of information at that point (i.e., the direction of the characteristics). Simple central differencing of the numerical flux then becomes a sophisticated upwind scheme in which the correct direction of upwinding is automatically achieved.
The numerical flux, F^i say, must be calculated by you in terms of the left and right values of the solution vector U (denoted by UL and UR respectively), at each mid-point of the mesh xj-12=xj-1+xj/2, for j=2,3,,npts. The left and right values are calculated by nag_pde_parab_1d_cd_ode_remesh (d03psc) from two adjacent mesh points using a standard upwind technique combined with a Van Leer slope-limiter (see LeVeque (1990)). The physically correct value for F^i is derived from the solution of the Riemann problem given by
Ui t + Fi y =0, (5)
where y=x-xj-12, i.e., y=0 corresponds to x=xj-12, with discontinuous initial values U=UL for y<0 and U=UR for y>0, using an approximate Riemann solver. This applies for either of the systems (1) or (2); the numerical flux is independent of the functions Pi,j, Ci, Di and Si. A description of several approximate Riemann solvers can be found in LeVeque (1990) and Berzins et al. (1989). Roe's scheme (see Roe (1981)) is perhaps the easiest to understand and use, and a brief summary follows. Consider the system of PDEs Ut+Fx=0 or equivalently Ut+AUx=0. Provided the system is linear in U, i.e., the Jacobian matrix A does not depend on U, the numerical flux F^ is given by
F^=12 FL+FR-12k=1npdeαkλkek, (6)
where FL (FR) is the flux F calculated at the left (right) value of U, denoted by UL (UR); the λk are the eigenvalues of A; the ek are the right eigenvectors of A; and the αk are defined by
UR-UL=k=1npdeαkek. (7)
Examples are given in the documents for nag_pde_parab_1d_cd (d03pfc) and nag_pde_parab_1d_cd_ode (d03plc).
If the system is nonlinear, Roe's scheme requires that a linearized Jacobian is found (see Roe (1981)).
The functions Pi,j, Ci, Di and Si (but not Fi) must be specified in pdedef. The numerical flux F^i must be supplied in numflx. For problems in the form (2), NULL may be used for pdedef. In this case, a default function sets the matrix with entries Pi,j to the identity matrix, and the functions Ci, Di and Si to zero.
For second-order problems, i.e., diffusion terms are present, a boundary condition is required for each PDE at both boundaries for the problem to be well-posed. If there are no diffusion terms present, then the continuous PDE problem generally requires exactly one boundary condition for each PDE, that is npde boundary conditions in total. However, in common with most discretization schemes for first-order problems, a numerical boundary condition is required at the other boundary for each PDE. In order to be consistent with the characteristic directions of the PDE system, the numerical boundary conditions must be derived from the solution inside the domain in some manner (see below). You must supply both types of boundary conditions, i.e., a total of npde conditions at each boundary point.
The position of each boundary condition should be chosen with care. In simple terms, if information is flowing into the domain then a physical boundary condition is required at that boundary, and a numerical boundary condition is required at the other boundary. In many cases the boundary conditions are simple, e.g., for the linear advection equation. In general you should calculate the characteristics of the PDE system and specify a physical boundary condition for each of the characteristic variables associated with incoming characteristics, and a numerical boundary condition for each outgoing characteristic.
A common way of providing numerical boundary conditions is to extrapolate the characteristic variables from the inside of the domain (note that when using banded matrix algebra the fixed bandwidth means that only linear extrapolation is allowed, i.e., using information at just two interior points adjacent to the boundary). For problems in which the solution is known to be uniform (in space) towards a boundary during the period of integration then extrapolation is unnecessary; the numerical boundary condition can be supplied as the known solution at the boundary. Another method of supplying numerical boundary conditions involves the solution of the characteristic equations associated with the outgoing characteristics. Examples of both methods can be found in the documents for nag_pde_parab_1d_cd (d03pfc) and nag_pde_parab_1d_cd_ode (d03plc).
The boundary conditions must be specified in bndary in the form
GiLx,t,U,V,V.=0  at ​x=a,  i=1,2,,npde, (8)
at the left-hand boundary, and
GiRx,t,U,V,V.=0  at ​x=b,  i=1,2,,npde, (9)
at the right-hand boundary.
Note that spatial derivatives at the boundary are not passed explicitly to bndary, but they can be calculated using values of U at and adjacent to the boundaries if required. However, it should be noted that instabilities may occur if such one-sided differencing opposes the characteristic direction at the boundary.
The algebraic-differential equation system which is defined by the functions Ri must be specified in odedef. You must also specify the coupling points ξ (if any) in the array xi.
In total there are npde×npts+ncode ODEs in the time direction. This system is then integrated forwards in time using a BDF or Theta method, optionally switching between Newton's method and functional iteration (see Berzins et al. (1989) and the references therein).
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 by the function 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.
The problem is subject to the following restrictions:
(i) In (1), V.jt, for j=1,2,,ncode, may only appear linearly in the functions Si, for i=1,2,,npde, with a similar restriction for GiL and GiR;
(ii) Pi,j, Fi, Ci and Si must not depend on any space derivatives; and Pi,j, Ci, Di and Fi 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, Ci, Di and Si is done by calling the pdedef at a point approximately midway between each pair of mesh points 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.
For further details of the scheme, see Pennington and Berzins (1994) and the references therein.

4  References

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
Furzeland R M (1984) The construction of adaptive space meshes TNER.85.022 Thornton Research Centre, Chester
Hirsch C (1990) Numerical Computation of Internal and External Flows, Volume 2: Computational Methods for Inviscid and Viscous Flows John Wiley
LeVeque R J (1990) Numerical Methods for Conservation Laws Birkhäuser Verlag
Pennington S V and Berzins M (1994) New NAG Library software for first-order partial differential equations ACM Trans. Math. Softw. 20 63–99
Roe P L (1981) Approximate Riemann solvers, parameter vectors, and difference schemes J. Comput. Phys. 43 357–372

5  Arguments

1:     npdeIntegerInput
On entry: the number of PDEs to be solved.
Constraint: npde1.
2:     tsdouble *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.
3:     toutdoubleInput
On entry: the final value of t to which the integration is to be carried out.
4:     pdedeffunction, supplied by the userExternal Function
pdedef must evaluate the functions Pi,j, Ci, Di and Si which partially define the system of PDEs. Pi,j and Ci may depend on x, t, U and V; Di may depend on x, t, U, Ux and V; and Si may depend on x, t, U, V and linearly on V.. pdedef is called approximately midway between each pair of mesh points in turn by nag_pde_parab_1d_cd_ode_remesh (d03psc). The argument may be specified as NULL for problems in the form (2).
The specification of pdedef is:
void  pdedef (Integer npde, double t, double x, const double u[], const double ux[], Integer ncode, const double v[], const double vdot[], double p[], double c[], double d[], double s[], Integer *ires, Nag_Comm *comm)
1:     npdeIntegerInput
On entry: the number of PDEs in the system.
2:     tdoubleInput
On entry: the current value of the independent variable t.
3:     xdoubleInput
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:     ncodeIntegerInput
On entry: the number of coupled ODEs in the system.
7:     v[ncode]const doubleInput
On entry: if ncode>0, v[i-1] contains the value of the component Vit, for i=1,2,,ncode.
8:     vdot[ncode]const doubleInput
On entry: if ncode>0, vdot[i-1] contains the value of component V.it, for i=1,2,,ncode.
Note: V.it, for i=1,2,,ncode, may only appear linearly in Sj, for j=1,2,,npde.
9:     p[npde×npde]doubleOutput
On exit: p[npde×j-1+i-1] must be set to the value of Pi,jx,t,U,V, for i=1,2,,npde and j=1,2,,npde.
10:   c[npde]doubleOutput
On exit: c[i-1] must be set to the value of Cix,t,U,V, for i=1,2,,npde.
11:   d[npde]doubleOutput
On exit: d[i-1] must be set to the value of Dix,t,U,Ux,V, for i=1,2,,npde.
12:   s[npde]doubleOutput
On exit: s[i-1] must be set to the value of Six,t,U,V,V., for i=1,2,,npde.
13:   iresInteger *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, then nag_pde_parab_1d_cd_ode_remesh (d03psc) returns to the calling function with the error indicator set to fail.code= NE_FAILED_DERIV.
14:   commNag_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_cd_ode_remesh (d03psc) you may allocate memory and initialize these pointers with various quantities for use by pdedef when called from nag_pde_parab_1d_cd_ode_remesh (d03psc) (see Section 3.2.1.1 in the Essential Introduction).
5:     numflxfunction, supplied by the userExternal Function
numflx must supply the numerical flux for each PDE given the left and right values of the solution vector u. numflx is called approximately midway between each pair of mesh points in turn by nag_pde_parab_1d_cd_ode_remesh (d03psc).
The specification of numflx is:
void  numflx (Integer npde, double t, double x, Integer ncode, const double v[], const double uleft[], const double uright[], double flux[], Integer *ires, Nag_Comm *comm, Nag_D03_Save *saved)
1:     npdeIntegerInput
On entry: the number of PDEs in the system.
2:     tdoubleInput
On entry: the current value of the independent variable t.
3:     xdoubleInput
On entry: the current value of the space variable x.
4:     ncodeIntegerInput
On entry: the number of coupled ODEs in the system.
5:     v[ncode]const doubleInput
On entry: if ncode>0, v[i-1] contains the value of the component Vit, for i=1,2,,ncode.
6:     uleft[npde]const doubleInput
On entry: uleft[i-1] contains the left value of the component Uix, for i=1,2,,npde.
7:     uright[npde]const doubleInput
On entry: uright[i-1] contains the right value of the component Uix, for i=1,2,,npde.
8:     flux[npde]doubleOutput
On exit: flux[i-1] must be set to the numerical flux F^i, for i=1,2,,npde.
9:     iresInteger *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, then nag_pde_parab_1d_cd_ode_remesh (d03psc) returns to the calling function with the error indicator set to fail.code= NE_FAILED_DERIV.
10:   commNag_Comm *
Pointer to structure of type Nag_Comm; the following members are relevant to numflx.
userdouble *
iuserInteger *
pPointer 
The type Pointer will be void *. Before calling nag_pde_parab_1d_cd_ode_remesh (d03psc) you may allocate memory and initialize these pointers with various quantities for use by numflx when called from nag_pde_parab_1d_cd_ode_remesh (d03psc) (see Section 3.2.1.1 in the Essential Introduction).
11:   savedNag_D03_Save *Communication Structure
If numflx calls one of the approximate Riemann solvers nag_pde_parab_1d_euler_roe (d03puc)nag_pde_parab_1d_euler_osher (d03pvc)nag_pde_parab_1d_euler_hll (d03pwc) or nag_pde_parab_1d_euler_exact (d03pxc) then saved is used to pass data concerning the computation to the solver. You should not change the components of saved.
6:     bndaryfunction, supplied by the userExternal Function
bndary must evaluate the functions GiL and GiR which describe the physical and numerical boundary conditions, as given by (8) and (9).
The specification of bndary is:
void  bndary (Integer npde, Integer npts, double t, const double x[], const double u[], Integer ncode, const double v[], const double vdot[], Integer ibnd, double g[], Integer *ires, Nag_Comm *comm)
1:     npdeIntegerInput
On entry: the number of PDEs in the system.
2:     nptsIntegerInput
On entry: the number of mesh points in the interval a,b.
3:     tdoubleInput
On entry: the current value of the independent variable t.
4:     x[npts]const doubleInput
On entry: the mesh points in the spatial direction. x[0] corresponds to the left-hand boundary, a, and x[npts-1] corresponds to the right-hand boundary, b.
5:     u[npde×npts]const doubleInput
On entry: u[npde×j-1+i-1] contains the value of the component Uix,t at x=x[j-1], for i=1,2,,npde and j=1,2,,npts.
Note:  if banded matrix algebra is to be used then the functions GiL and GiR may depend on the value of Uix,t at the boundary point and the two adjacent points only.
6:     ncodeIntegerInput
On entry: the number of coupled ODEs in the system.
7:     v[ncode]const doubleInput
On entry: if ncode>0, v[i-1] contains the value of the component Vit, for i=1,2,,ncode.
8:     vdot[ncode]const doubleInput
On entry: if ncode>0, vdot[i-1] contains the value of component V.it, for i=1,2,,ncode.
Note: V.it, for i=1,2,,ncode, may only appear linearly in GjL and GjR, for j=1,2,,npde.
9:     ibndIntegerInput
On entry: specifies which boundary conditions are to be evaluated.
ibnd=0
bndary must evaluate the left-hand boundary condition at x=a.
ibnd0
bndary must evaluate the right-hand boundary condition at x=b.
10:   g[npde]doubleOutput
On exit: g[i-1] must contain the ith component of either GiL or GiR in (8) and (9), depending on the value of ibnd, for i=1,2,,npde.
11:   iresInteger *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, then nag_pde_parab_1d_cd_ode_remesh (d03psc) returns to the calling function with the error indicator set to fail.code= NE_FAILED_DERIV.
12:   commNag_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_cd_ode_remesh (d03psc) you may allocate memory and initialize these pointers with various quantities for use by bndary when called from nag_pde_parab_1d_cd_ode_remesh (d03psc) (see Section 3.2.1.1 in the Essential Introduction).
7:     uvinitfunction, 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 ncode, double v[], Nag_Comm *comm)
1:     npdeIntegerInput
On entry: the number of PDEs in the system.
2:     nptsIntegerInput
On entry: the number of mesh points in the interval [a,b].
3:     nxiIntegerInput
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 ODE/PDE coupling point, ξi, for i=1,2,,nxi.
6:     u[npde×npts]doubleOutput
On exit: if nxi>0, u[npde×j-1+i-1] contains the value of the component Uixj,t0, for i=1,2,,npde and j=1,2,,npts.
7:     ncodeIntegerInput
On entry: the number of coupled ODEs in the system.
8:     v[ncode]doubleOutput
On exit: if ncode>0, v[i-1] must contain the value of component Vit0, for i=1,2,,ncode.
9:     commNag_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_cd_ode_remesh (d03psc) you may allocate memory and initialize these pointers with various quantities for use by uvinit when called from nag_pde_parab_1d_cd_ode_remesh (d03psc) (see Section 3.2.1.1 in the Essential Introduction).
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,,ncode, all evaluated at t=ts.
9:     nptsIntegerInput
On entry: the number of mesh points in the interval a,b.
Constraint: npts3.
10:   x[npts]doubleInput/Output
On entry: the 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:   ncodeIntegerInput
On entry: the number of coupled ODE components.
Constraint: ncode0.
12:   odedeffunction, supplied by the userExternal Function
odedef must evaluate the functions R, which define the system of ODEs, as given in (4).
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_cd_ode_remesh (d03psc).
The specification of odedef is:
void  odedef (Integer npde, double t, Integer ncode, const double v[], const double vdot[], Integer nxi, const double xi[], const double ucp[], const double ucpx[], const double ucpt[], double r[], Integer *ires, Nag_Comm *comm)
1:     npdeIntegerInput
On entry: the number of PDEs in the system.
2:     tdoubleInput
On entry: the current value of the independent variable t.
3:     ncodeIntegerInput
On entry: the number of coupled ODEs in the system.
4:     v[ncode]const doubleInput
On entry: if ncode>0, v[i-1] contains the value of the component Vit, for i=1,2,,ncode.
5:     vdot[ncode]const doubleInput
On entry: if ncode>0, vdot[i-1] contains the value of component V.it, for i=1,2,,ncode.
6:     nxiIntegerInput
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 point, ξi, for i=1,2,,nxi.
8:     ucp[npde×nxi]const doubleInput
On entry: if nxi>0, ucp[npde×j-1+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
On entry: if nxi>0, ucpx[npde×j-1+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:   ucpt[npde×nxi]const doubleInput
On entry: if nxi>0, ucpt[npde×j-1+i-1] contains the value of Ui t  at the coupling point x=ξj, for i=1,2,,npde and j=1,2,,nxi.
11:   r[ncode]doubleOutput
On exit: r[i-1] must contain the ith component of R, for i=1,2,,ncode, where R is defined as
R=L-MV.-NUt*, (10)
or
R=-MV.-NUt*. (11)
The definition of R is determined by the input value of ires.
12:   iresInteger *Input/Output
On entry: the form of R that must be returned in the array r.
ires=1
Equation (10) must be used.
ires=-1
Equation (11) 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, then nag_pde_parab_1d_cd_ode_remesh (d03psc) returns to the calling function with the error indicator set to fail.code= NE_FAILED_DERIV.
13:   commNag_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_cd_ode_remesh (d03psc) you may allocate memory and initialize these pointers with various quantities for use by odedef when called from nag_pde_parab_1d_cd_ode_remesh (d03psc) (see Section 3.2.1.1 in the Essential Introduction).
13:   nxiIntegerInput
On entry: the number of ODE/PDE coupling points.
Constraints:
  • if ncode=0, nxi=0;
  • if ncode>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:   neqnIntegerInput
On entry: the number of ODEs in the time direction.
Constraint: neqn=npde×npts+ncode.
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.
Constraint: atol[i-1]0.0 for all relevant i.
Note: corresponding elements of rtol and atol cannot both be 0.0.
18:   itolIntegerInput
On entry: a value to indicate the form of the local error test. If ei is the estimated local error for u[i-1], for i=1,2,,neqn, and , denotes the norm, then the error test to be satisfied is ei<1.0. itol indicates to nag_pde_parab_1d_cd_ode_remesh (d03psc) whether to interpret either or both of rtol and atol as a vector or scalar in the formation of the weights wi used in the calculation of the norm (see the description of norm):
itol rtol atol wi
1 scalar scalar rtol[0]×u[i-1]+atol[0]
2 scalar vector rtol[0]×u[i-1]+atol[i-1]
3 vector scalar rtol[i-1]×u[i-1]+atol[0]
4 vector vector rtol[i-1]×u[i-1]+atol[i-1]
Constraint: itol=1, 2, 3 or 4.
19:   normNag_NormTypeInput
On entry: the type of norm to be used.
norm=Nag_OneNorm
Averaged L1 norm.
norm=Nag_TwoNorm
Averaged L2 norm.
If Unorm denotes the norm of the vector u of length neqn, then for the averaged L1 norm
Unorm=1neqni=1neqnu[i-1]/wi,
and for the averaged L2 norm
Unorm=1neqni=1neqnu[i-1]/wi2,
See the description of itol for the formulation of the weight vector w.
Constraint: norm=Nag_OneNorm or Nag_TwoNorm.
20:   laoptNag_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 (ncode=0). Also, the banded option should not be used if the boundary conditions involve solution components at points other than the boundary and the immediately adjacent two points.
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, then 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 is algopt[0]=1.0.
If algopt[0]=2.0, then 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, then 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, then the Petzold test is used. If algopt[3]=2.0, then the Petzold test is not used. The default value is algopt[3]=1.0.
If algopt[0]=1.0, then 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, then switching is allowed and if algopt[6]=2.0, then 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, then 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, then 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, i.e., 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 the 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]
Used as the relative pivot threshold during subsequent Jacobian decompositions (see algopt[28]) below which an internal error is invoked. algopt[29] must be greater than zero, otherwise the default value is used. 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 matrix is found to be numerically singular (see algopt[28]). The default value is algopt[29]=0.0001.
22:   remeshNag_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_cd_ode_remesh (d03psc). Remeshing can be switched off or on at specified times by using appropriate values for the arguments nrmesh and trmesh at each call.
23:   nxfixIntegerInput
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[dim]const doubleInput
Note: the dimension, dim, of the array xfix must be at least max1,nxfix.
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[npts-1] 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:   nrmeshIntegerInput
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_cd_ode_remesh (d03psc) to give greater flexibility over the times of remeshing.
26:   dxmeshdoubleInput
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>xiold+dxmesh×xi+1old-xiold
or
xinew<xiold-dxmesh×xiold-xi- 1old
dxmesh thus imposes a lower limit on the difference between one mesh and the next.
Constraint: dxmesh0.0.
27:   trmeshdoubleInput
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_cd_ode_remesh (d03psc) to force remeshing at several specified times.
28:   ipminfIntegerInput
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:   xratiodoubleInput
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:   condoubleInput
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.0/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:   monitffunction, supplied by the userExternal Function
monitf must supply and evaluate a remesh monitor function to indicate the solution behaviour of interest.
monitf will never be called and the NAG defined null void function pointer, NULLFN, can be supplied in the call to nag_pde_parab_1d_cd_ode_remesh (d03psc).
The specification of monitf is:
void  monitf (double t, Integer npts, Integer npde, const double x[], const double u[], double fmon[], Nag_Comm *comm)
1:     tdoubleInput
On entry: the current value of the independent variable t.
2:     nptsIntegerInput
On entry: the number of mesh points in the interval a,b.
3:     npdeIntegerInput
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
On entry: u[npde×j-1+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:     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.
7:     commNag_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_cd_ode_remesh (d03psc) you may allocate memory and initialize these pointers with various quantities for use by monitf when called from nag_pde_parab_1d_cd_ode_remesh (d03psc) (see Section 3.2.1.1 in the Essential Introduction).
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:   lrsaveIntegerInput
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=3×npde-1, for PDE problems only (no coupled ODEs); or
mlu=neqn-1, for coupled PDE/ODE problems.
nwkres= npde×2×npts+6×nxi+3×npde+26+nxi+ncode+7×npts+nxfix+1, when ​ncode>0​ and ​nxi>0; or npde×2×npts+3×npde+32+ncode+7×npts+nxfix+2, when ​ncode>0​ and ​nxi=0; or ​ npde×2×npts+3×npde+32+7×npts+nxfix+3, when ​ncode=0.
lenode= 6+intalgopt[1]×neqn+50, when the BDF method is used; or 9×neqn+50, when the Theta method is used.  
Note: when laopt=Nag_LinAlgSparse, 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.
If ind=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[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 evaluating 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 BDF method last used in the time integration, if applicable. When the Theta method is used, isave[3] contains no useful information.
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.
35:   lisaveIntegerInput
On entry: the dimension of the array isave. Its size depends on the type of matrix algebra selected:
  • if laopt=Nag_LinAlgFull, lisave25;
  • if laopt=Nag_LinAlgBand, lisaveneqn+nxfix+25;
  • if laopt=Nag_LinAlgSparse, lisave25×neqn+nxfix+25.
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:   itaskIntegerInput
On entry: the task to be performed by the ODE integrator.
itask=1
Normal computation of output values u at t=tout (by overshooting and interpolating).
itask=2
Take one step in the time direction 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:   itraceIntegerInput
On entry: the level of trace information required from nag_pde_parab_1d_cd_ode_remesh (d03psc) and the underlying ODE solver. itrace may take the value -1, 0, 1, 2 or 3.
itrace=-1
No output is generated.
itrace=0
Only warning messages from the PDE solver are printed.
itrace>0
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.
If itrace<-1, then -1 is assumed and similarly if itrace>3, then 3 is assumed.
The advisory messages are given in greater detail as itrace increases.
38:   outfileconst 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:   indInteger *Input/Output
On entry: indicates whether this is a continuation call or a new integration.
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, fail, nrmesh and trmesh may be reset between calls to nag_pde_parab_1d_cd_ode_remesh (d03psc).
Constraint: ind=0 or 1.
On exit: ind=1.
40:   commNag_Comm *Communication Structure
The NAG communication argument (see Section 3.2.1.1 in the Essential Introduction).
41:   savedNag_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:   failNagError *Input/Output
The NAG error argument (see Section 3.6 in the Essential Introduction).

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.
NE_ALLOC_FAIL
Dynamic memory allocation failed.
NE_BAD_MONIT
fmon is negative at one or more mesh points, or zero mesh spacing has been obtained due to a poor monitor function.
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, numflx, bndary or odedef.
NE_FAILED_START
atol and rtol were too small to start integration.
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.
Underlying ODE solver cannot make further progress from the point ts with the supplied values of atol and rtol. ts=value.
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, 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 a call to user-supplied functions pdedef, numflx, 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, ncode=value.
Constraint: ncode0.
On entry, npde=value.
Constraint: npde1.
On entry, npts=value.
Constraint: npts3.
On entry, nxfix=value.
Constraint: nxfix0.
NE_INT_2
On entry, corresponding elements atol[I-1] and rtol[J-1] are both zero: I=value and J=value.
On entry, lisave is too small: lisave=value. Minimum possible dimension: value.
On entry, lrsave is too small: lrsave=value. Minimum possible dimension: value.
On entry, ncode=value and nxi=value.
Constraint: nxi=0 when ncode=0.
On entry, ncode=value and nxi=value.
Constraint: nxi0 when ncode>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 ncode=value.
Constraint: neqn=npde×npts+ncode.
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.
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_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, xi[I]=value and xi[I-1]=value.
Constraint: xi[I]>xi[I-1].
On entry, mesh points x appear to be badly ordered: I=value, x[I-1]=value, J=value and x[J-1]=value.
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_cd_ode_remesh (d03psc).
NE_SING_JAC
Singular Jacobian of ODE system. Check problem formulation.
NE_TIME_DERIV_DEP
The functions P, D, or C appear to depend on time derivatives.
NE_USER_STOP
In evaluating residual of ODE system, ires=2 has been set in user-supplied functions pdedef, numflx, bndary or odedef. Integration is successful as far as ts: ts=value.
NE_ZERO_WTS
Zero error weights encountered during time integration.

7  Accuracy

nag_pde_parab_1d_cd_ode_remesh (d03psc) 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_cd_ode_remesh (d03psc) is threaded by NAG for parallel execution in multithreaded implementations of the NAG Library.
nag_pde_parab_1d_cd_ode_remesh (d03psc) 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 Users' Note for your implementation for any additional implementation-specific information.

9  Further Comments

nag_pde_parab_1d_cd_ode_remesh (d03psc) is designed to solve systems of PDEs in conservative form, with optional source terms which are independent of space derivatives, and optional second-order diffusion terms. The use of the function to solve systems which are not naturally in this form is discouraged, and you are advised to use one of the central-difference scheme functions for such problems.
You should be aware of the stability limitations for hyperbolic PDEs. For most problems with small error tolerances the ODE integrator does not attempt unstable time steps, but in some cases a maximum time step should be imposed using algopt[12]. It is worth experimenting with this argument, particularly if the integration appears to progress unrealistically fast (with large time steps). Setting the maximum time step to the minimum mesh size is a safe measure, although in some cases this may be too restrictive.
Problems with source terms should be treated with caution, as it is known that for large source terms stable and reasonable looking solutions can be obtained which are in fact incorrect, exhibiting non-physical speeds of propagation of discontinuities (typically one spatial mesh point per time step). It is essential to employ a very fine mesh for problems with source terms and discontinuities, and to check for non-physical propagation speeds by comparing results for different mesh sizes. Further details and an example can be found in Pennington and Berzins (1994).
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.

10  Example

For this function two examples are presented, with a main program and two example problems given in Example 1 (ex1) and Example 2 (ex2).
Example 1 (ex1)
This example is a simple model of the advection and diffusion of a cloud of material:
U t + W U x = C 2 U x2 ,
for x0,1 and t00.3. In this example the constant wind speed W=1 and the diffusion coefficient C=0.002.
The cloud does not reach the boundaries during the time of integration, and so the two (physical) boundary conditions are simply U0,t=U1,t=0.0, and the initial condition is
U x,0 = sin πx-a b-a ,  axb,
and Ux,0=0 elsewhere, where a=0.2 and b=0.4.
The numerical flux is simply F^=WUL.
The monitor function for remeshing is taken to be the absolute value of the second derivative of U.
Example 2 (ex2)
This example is a linear advection equation with a nonlinear source term and discontinuous initial profile:
u t + u x =-puu-1u-12,
for 0x1 and t0. The discontinuity is modelled by a ramp function of width 0.01 and gradient 100, so that the exact solution at any time t0 is
ux,t=1.0+maxminδ,0,-1,
where δ=1000.1-x+t. The initial profile is given by the exact solution. The characteristic points into the domain at x=0 and out of the domain at x=1, and so a physical boundary condition u0,t=1 is imposed at x=0, with a numerical boundary condition at x=1 which can be specified as u1,t=0 since the discontinuity does not reach x=1 during the time of integration.
The numerical flux is simply F^=UL at all times.
The remeshing monitor function (described below) is chosen to create an increasingly fine mesh towards the discontinuity in order to ensure good resolution of the discontinuity, but without loss of efficiency in the surrounding regions. However, refinement must be limited so that the time step required for stability does not become unrealistically small. The region of refinement must also keep up with the discontinuity as it moves across the domain, and hence it cannot be so small that the discontinuity moves out of the refined region between remeshing.
The above requirements mean that the use of the first or second spatial derivative of U for the monitor function is inappropriate; the large relative size of either derivative in the region of the discontinuity leads to extremely small mesh-spacing in a very limited region, and the solution is then far more expensive than for a very fine fixed mesh.
An alternative monitor function based on a cosine function proves very successful. It is only semi-automatic as it requires some knowledge of the solution (for problems without an exact solution an initial approximate solution can be obtained using a coarse fixed mesh). On each call to monitf the discontinuity is located by finding the maximum spatial derivative of the solution. On the first call the desired width of the region of nonzero monitor function is set (this can be changed at a later time if desired). Then on each call the monitor function is assigned using a cosine function so that it has a value of one at the discontinuity down to zero at the edges of the predetermined region of refinement, and zero outside the region. Thus the monitor function and the subsequent refinement are limited, and the region is large enough to ensure that there is always sufficient refinement at the discontinuity.

10.1  Program Text

Program Text (d03psce.c)

10.2  Program Data

None.

10.3  Program Results

Program Results (d03psce.r)

Produced by GNUPLOT 4.4 patchlevel 0 Example Program 1 Advection and Diffusion of a Cloud of Material u(x,t) 0 0.05 0.1 0.15 0.2 0.25 0.3 Time 0 0.2 0.4 0.6 0.8 1 x 0 0.1 0.2 0.3 0.4 0.5 0.6 0.7 0.8 0.9 1
Produced by GNUPLOT 4.4 patchlevel 0 Example Program 2 Linear Advection Equation with Non-linear Source Term u(x,t) 0 0.05 0.1 0.15 0.2 0.25 0.3 0.35 0.4 Time 0 0.2 0.4 0.6 0.8 1 x -0.2 0 0.2 0.4 0.6 0.8 1 1.2

nag_pde_parab_1d_cd_ode_remesh (d03psc) (PDF version)
d03 Chapter Contents
d03 Chapter Introduction
NAG Library Manual

© The Numerical Algorithms Group Ltd, Oxford, UK. 2014