NAG CL Interface
d03pjc (dim1_​parab_​dae_​coll)

Settings help

CL Name Style:


1 Purpose

d03pjc integrates a system of linear or nonlinear parabolic partial differential equations (PDEs), in one space variable with scope for coupled ordinary differential equations (ODEs). The spatial discretization is performed using a Chebyshev C0 collocation method, and the method of lines is employed to reduce the PDEs to a system of ODEs. The resulting system is solved using a backward differentiation formula (BDF) method or a Theta method (switching between Newton's method and functional iteration).

2 Specification

#include <nag.h>
void  d03pjc (Integer npde, Integer m, double *ts, double tout,
void (*pdedef)(Integer npde, double t, const double x[], Integer nptl, const double u[], const double ux[], Integer nv, const double v[], const double vdot[], double p[], double q[], double r[], Integer *ires, Nag_Comm *comm),
void (*bndary)(Integer npde, double t, const double u[], const double ux[], Integer nv, const double v[], const double vdot[], Integer ibnd, double beta[], double gamma[], Integer *ires, Nag_Comm *comm),
double u[], Integer nbkpts, const double xbkpts[], Integer npoly, Integer npts, double x[], Integer nv,
void (*odedef)(Integer npde, double t, Integer nv, const double v[], const double vdot[], Integer nxi, const double xi[], const double ucp[], const double ucpx[], const double rcp[], const double ucpt[], const double ucptx[], double f[], Integer *ires, Nag_Comm *comm),
Integer nxi, const double xi[], Integer neqn,
void (*uvinit)(Integer npde, Integer npts, const double x[], double u[], Integer nv, double v[], Nag_Comm *comm),
const double rtol[], const double atol[], Integer itol, Nag_NormType norm, Nag_LinAlgOption laopt, const double algopt[], 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)
The function may be called by the names: d03pjc, nag_pde_dim1_parab_dae_coll or nag_pde_parab_1d_coll_ode.

3 Description

d03pjc integrates the system of parabolic-elliptic equations and coupled ODEs
j=1npdePi,j Uj t +Qi=x-m x (xmRi),  i=1,2,,npde,  axb,tt0, (1)
Fi(t,V,V.,ξ,U*,Ux*,R*,Ut*,Uxt*)=0,  i=1,2,,nv, (2)
where (1) defines the PDE part and (2) generalizes the coupled ODE part of the problem.
In (1), Pi,j and Ri depend on x, t, U, Ux, and V; Qi depends on x, t, U, Ux, V and linearly on V.. The vector U is the set of PDE solution values
U (x,t) = [ U 1 (x,t) ,, U npde (x,t) ] T ,  
and the vector Ux is the partial derivative with respect to x. Note that Pi,j, Qi and Ri must not depend on U t . The vector V is the set of ODE solution values
V(t)=[V1(t),,Vnv(t)]T,  
and V. denotes its derivative with respect to time.
In (2), ξ represents a vector of nξ spatial coupling points at which the ODEs are coupled to the PDEs. These points may or may not be equal to some of the PDE spatial mesh points. U*, Ux*, R*, Ut* and Uxt* are the functions U, Ux, R, Ut and Uxt evaluated at these coupling points. Each Fi may only depend linearly on time derivatives. Hence the equation (2) may be written more precisely as
F=G-AV.-B ( Ut* Uxt* ) , (3)
where F = [F1,,Fnv] T , G is a vector of length nv, A is an nv by nv matrix, B is an nv by (nξ×npde) matrix and the entries in G, A and B may depend on t, ξ, U*, Ux* and V. In practice you need only supply a vector of information to define the ODEs and not the matrices A and B. (See Section 5 for the specification of odedef.)
The integration in time is from t0 to tout, over the space interval axb, where a=x1 and b=xnbkpts are the leftmost and rightmost of a user-defined set of break-points x1,x2,,xnbkpts. The coordinate system in space is defined by the value of m; m=0 for Cartesian coordinates, m=1 for cylindrical polar coordinates and m=2 for spherical polar coordinates.
The PDE system which is defined by the functions Pi,j, Qi and Ri must be specified in pdedef.
The initial values of the functions U(x,t) and V(t) must be given at t=t0. These values are calculated in uvinit.
The functions Ri which may be thought of as fluxes, are also used in the definition of the boundary conditions. The boundary conditions must have the form
βi(x,t)Ri(x,t,U,Ux,V)=γi(x,t,U,Ux,V,V.),  i=1,2,,npde, (4)
where x=a or x=b. The functions γi may only depend linearly on V..
The boundary conditions must be specified in bndary.
The algebraic-differential equation system which is defined by the functions Fi must be specified in odedef. You must also specify the coupling points ξ in the array xi. Thus, the problem is subject to the following restrictions:
  1. (i)in (1), V.j(t), for j=1,2,,nv, may only appear linearly in the functions Qi, for i=1,2,,npde, with a similar restriction for γ;
  2. (ii)Pi,j and the flux Ri must not depend on any time derivatives;
  3. (iii)t0<tout, so that integration is in the forward direction;
  4. (iv)the evaluation of the functions Pi,j, Qi and Ri is done at both the break-points and internally selected points for each element in turn, that is Pi,j, Qi and Ri are evaluated twice at each break-point. Any discontinuities in these functions must, therefore, be at one or more of the mesh points;
  5. (v)at least one of the functions Pi,j must be nonzero so that there is a time derivative present in the PDE problem;
  6. (vi)if m>0 and x1=0.0, which is the left boundary point, then it must be ensured that the PDE solution is bounded at this point. This can be done either by specifying the solution at x=0.0 or by specifying a zero flux there, that is βi=1.0 and γi=0.0.
The parabolic equations are approximated by a system of ODEs in time for the values of Ui at the mesh points. This ODE system is obtained by approximating the PDE solution between each pair of break-points by a Chebyshev polynomial of degree npoly. The interval between each pair of break-points is treated by d03pjc as an element, and on this element, a polynomial and its space and time derivatives are made to satisfy the system of PDEs at npoly-1 spatial points, which are chosen internally by the code and the break-points. The user-defined break-points and the internally selected points together define the mesh. The smallest value that npoly can take is one, in which case, the solution is approximated by piecewise linear polynomials between consecutive break-points and the method is similar to an ordinary finite element method.
In total there are (nbkpts-1)×npoly+1 mesh points in the spatial direction, and npde×((nbkpts-1)×npoly+1)+nv ODEs in the time direction; one ODE at each break-point for each PDE component, npoly-1 ODEs for each PDE component between each pair of break-points, and nv coupled ODEs. The system is then integrated forwards in time using a Backward Differentiation Formula (BDF) method or a Theta method.

4 References

Berzins M (1990) Developments in the NAG Library software for parabolic equations Scientific Software Systems (eds J C Mason and M G Cox) 59–72 Chapman and Hall
Berzins M and Dew P M (1991) Algorithm 690: Chebyshev polynomial software for elliptic-parabolic systems of PDEs ACM Trans. Math. Software 17 178–206
Berzins M, Dew P M and Furzeland R M (1988) Software tools for time-dependent equations in simulation and optimization of large systems Proc. IMA Conf. Simulation and Optimization (ed A J Osiadcz) 35–50 Clarendon Press, Oxford
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
Zaturska N B, Drazin P G and Banks W H H (1988) On the flow of a viscous fluid driven along a channel by a suction at porous walls Fluid Dynamics Research 4

5 Arguments

1: npde Integer Input
On entry: the number of PDEs to be solved.
Constraint: npde1.
2: m Integer Input
On entry: the coordinate system used:
m=0
Indicates Cartesian coordinates.
m=1
Indicates cylindrical polar coordinates.
m=2
Indicates spherical polar coordinates.
Constraint: m=0, 1 or 2.
3: ts double * Input/Output
On entry: the initial value of the independent variable t.
On exit: the value of t corresponding to the solution values in u. Normally ts=tout.
Constraint: ts<tout.
4: tout double Input
On entry: the final value of t to which the integration is to be carried out.
5: pdedef function, supplied by the user External Function
pdedef must compute the functions Pi,j, Qi and Ri which define the system of PDEs. The functions may depend on x, t, U, Ux and V; Qi may depend linearly on V.. The functions must be evaluated at a set of points.
The specification of pdedef is:
void  pdedef (Integer npde, double t, const double x[], Integer nptl, const double u[], const double ux[], Integer nv, const double v[], const double vdot[], double p[], double q[], double r[], Integer *ires, Nag_Comm *comm)
1: npde Integer Input
On entry: the number of PDEs in the system.
2: t double Input
On entry: the current value of the independent variable t.
3: x[dim] const double Input
On entry: contains a set of mesh points at which Pi,j, Qi and Ri are to be evaluated. x[0] and x[nptl-1] contain successive user-supplied break-points and the elements of the array will satisfy x[0]<x[1]<<x[nptl-1].
4: nptl Integer Input
On entry: the number of points at which evaluations are required (the value of npoly+1).
5: u[npde×nptl] const double Input
Note: the (i,j)th element of the matrix U is stored in u[(j-1)×npde+i-1].
On entry: u[(j-1)×npde+i-1] contains the value of the component Ui(x,t) where x=x[j-1], for i=1,2,,npde and j=1,2,,nptl.
6: ux[npde×nptl] const double Input
Note: the (i,j)th element of the matrix is stored in ux[(j-1)×npde+i-1].
On entry: ux[(j-1)×npde+i-1] contains the value of the component Ui(x,t) x where x=x[j-1], for i=1,2,,npde and j=1,2,,nptl.
7: nv Integer Input
On entry: the number of coupled ODEs in the system.
8: v[nv] const double Input
On entry: if nv>0, v[i-1] contains the value of the component Vi(t), for i=1,2,,nv.
9: vdot[nv] const double Input
On entry: if nv>0, vdot[i-1] contains the value of component V.i(t), for i=1,2,,nv.
Note: V.i(t), for i=1,2,,nv, may only appear linearly in Qj, for j=1,2,,npde.
10: p[dim] double Output
Note: the dimension, dim, of the array p is npde×npde×nptl.
where P(i,j,k) appears in this document, it refers to the array element p[(k-1)×npde×npde+(j-1)×npde+i-1].
On exit: p[ npde×npde×(k-1)+ npde×(j-1)+ (i-1)] must be set to the value of Pi,j(x,t,U,Ux,V) where x=x[k-1], for i=1,2,,npde, j=1,2,,npde and k=1,2,,nptl.
11: q[npde×nptl] double Output
Note: the (i,j)th element of the matrix Q is stored in q[(j-1)×npde+i-1].
On exit: q[(j-1)×npde+i-1] must be set to the value of Qi(x,t,U,Ux,V,V.) where x=x[j-1], for i=1,2,,npde and j=1,2,,nptl.
12: r[npde×nptl] double Output
Note: the (i,j)th element of the matrix R is stored in r[(j-1)×npde+i-1].
On exit: r[(j-1)×npde+i-1] must be set to the value of Ri(x,t,U,Ux,V) where x=x[i-1], for i=1,2,,npde and j=1,2,,nptl.
13: ires Integer * Input/Output
On entry: set to −1 or 1.
On exit: should usually remain unchanged. However, you may set ires to force the integration function to take certain actions as described below:
ires=2
Indicates to the integrator that control should be passed back immediately to the calling function with the error indicator set to fail.code= NE_USER_STOP.
ires=3
Indicates to the integrator that the current time step should be abandoned and a smaller time step used instead. You may wish to set ires=3 when a physically meaningless input or output value has been generated. If you consecutively set ires=3, d03pjc returns to the calling function with the error indicator set to fail.code= NE_FAILED_DERIV.
14: comm Nag_Comm *
Pointer to structure of type Nag_Comm; the following members are relevant to pdedef.
userdouble *
iuserInteger *
pPointer 
The type Pointer will be void *. Before calling d03pjc you may allocate memory and initialize these pointers with various quantities for use by pdedef when called from d03pjc (see Section 3.1.1 in the Introduction to the NAG Library CL Interface).
Note: pdedef should not return floating-point NaN (Not a Number) or infinity values, since these are not handled by d03pjc. If your code inadvertently does return any NaNs or infinities, d03pjc is likely to produce unexpected results.
6: bndary function, supplied by the user External Function
bndary must compute the functions βi and γi which define the boundary conditions as in equation (4).
The specification of bndary is:
void  bndary (Integer npde, double t, const double u[], const double ux[], Integer nv, const double v[], const double vdot[], Integer ibnd, double beta[], double gamma[], Integer *ires, Nag_Comm *comm)
1: npde Integer Input
On entry: the number of PDEs in the system.
2: t double Input
On entry: the current value of the independent variable t.
3: u[npde] const double Input
On entry: u[i-1] contains the value of the component Ui(x,t) at the boundary specified by ibnd, for i=1,2,,npde.
4: ux[npde] const double Input
On entry: ux[i-1] contains the value of the component Ui(x,t) x at the boundary specified by ibnd, for i=1,2,,npde.
5: nv Integer Input
On entry: the number of coupled ODEs in the system.
6: v[nv] const double Input
On entry: if nv>0, v[i-1] contains the value of the component Vi(t), for i=1,2,,nv.
7: vdot[nv] const double Input
On entry: if nv>0, vdot[i-1] contains the value of component V.i(t), for i=1,2,,nv.
Note: V.i(t), for i=1,2,,nv, may only appear linearly in Qj, for j=1,2,,npde.
8: ibnd Integer Input
On entry: specifies which boundary conditions are to be evaluated.
ibnd=0
bndary must set up the coefficients of the left-hand boundary, x=a.
ibnd0
bndary must set up the coefficients of the right-hand boundary, x=b.
9: beta[npde] double Output
On exit: beta[i-1] must be set to the value of βi(x,t) at the boundary specified by ibnd, for i=1,2,,npde.
10: gamma[npde] double Output
On exit: gamma[i-1] must be set to the value of γi(x,t,U,Ux,V,V.) at the boundary specified by ibnd, for i=1,2,,npde.
11: ires Integer * Input/Output
On entry: set to −1 or 1.
On exit: should usually remain unchanged. However, you may set ires to force the integration function to take certain actions as described below:
ires=2
Indicates to the integrator that control should be passed back immediately to the calling function with the error indicator set to fail.code= NE_USER_STOP.
ires=3
Indicates to the integrator that the current time step should be abandoned and a smaller time step used instead. You may wish to set ires=3 when a physically meaningless input or output value has been generated. If you consecutively set ires=3, d03pjc returns to the calling function with the error indicator set to fail.code= NE_FAILED_DERIV.
12: comm Nag_Comm *
Pointer to structure of type Nag_Comm; the following members are relevant to bndary.
userdouble *
iuserInteger *
pPointer 
The type Pointer will be void *. Before calling d03pjc you may allocate memory and initialize these pointers with various quantities for use by bndary when called from d03pjc (see Section 3.1.1 in the Introduction to the NAG Library CL Interface).
Note: bndary should not return floating-point NaN (Not a Number) or infinity values, since these are not handled by d03pjc. If your code inadvertently does return any NaNs or infinities, d03pjc is likely to produce unexpected results.
7: u[neqn] double Input/Output
On entry: if ind=1 the value of u must be unchanged from the previous call.
On exit: the computed solution Ui(xj,t), for i=1,2,,npde and j=1,2,,npts, and Vk(t), for k=1,2,,nv, evaluated at t=ts, as follows:
  • u[npde×(j-1)+i-1] contain Ui(xj,t), for i=1,2,,npde and j=1,2,,npts, and
  • u[npts×npde+i-1] contain Vi(t), for i=1,2,,nv.
8: nbkpts Integer Input
On entry: the number of break-points in the interval [a,b].
Constraint: nbkpts2.
9: xbkpts[nbkpts] const double Input
On entry: the values of the break-points in the space direction. xbkpts[0] must specify the left-hand boundary, a, and xbkpts[nbkpts-1] must specify the right-hand boundary, b.
Constraint: xbkpts[0]<xbkpts[1]<<xbkpts[nbkpts-1].
10: npoly Integer Input
On entry: the degree of the Chebyshev polynomial to be used in approximating the PDE solution between each pair of break-points.
Constraint: 1npoly49.
11: npts Integer Input
On entry: the number of mesh points in the interval [a,b].
Constraint: npts=(nbkpts-1)×npoly+1.
12: x[npts] double Output
On exit: the mesh points chosen by d03pjc in the spatial direction. The values of x will satisfy x[0]<x[1]<<x[npts-1].
13: nv Integer Input
On entry: the number of coupled ODE components.
Constraint: nv0.
14: odedef function, supplied by the user External Function
odedef must evaluate the functions F, which define the system of ODEs, as given in (3).
odedef will never be called and the NAG defined null void function pointer, NULLFN, can be supplied in the call to d03pjc.
The specification of odedef is:
void  odedef (Integer npde, double t, Integer nv, const double v[], const double vdot[], Integer nxi, const double xi[], const double ucp[], const double ucpx[], const double rcp[], const double ucpt[], const double ucptx[], double f[], Integer *ires, Nag_Comm *comm)
1: npde Integer Input
On entry: the number of PDEs in the system.
2: t double Input
On entry: the current value of the independent variable t.
3: nv Integer Input
On entry: the number of coupled ODEs in the system.
4: v[nv] const double Input
On entry: if nv>0, v[i-1] contains the value of the component Vi(t), for i=1,2,,nv.
5: vdot[nv] const double Input
On entry: if nv>0, vdot[i-1] contains the value of component V.i(t), for i=1,2,,nv.
6: nxi Integer Input
On entry: the number of ODE/PDE coupling points.
7: xi[nxi] const double Input
On entry: if nxi>0, xi[i-1] contains the ODE/PDE coupling points, ξi, for i=1,2,,nxi.
8: ucp[npde×nxi] const double Input
Note: the (i,j)th element of the matrix is stored in ucp[(j-1)×npde+i-1].
On entry: if nxi>0, ucp[(j-1)×npde+i-1] contains the value of Ui(x,t) at the coupling point x=ξj, for i=1,2,,npde and j=1,2,,nxi.
9: ucpx[npde×nxi] const double Input
Note: the (i,j)th element of the matrix is stored in ucpx[(j-1)×npde+i-1].
On entry: if nxi>0, ucpx[(j-1)×npde+i-1] contains the value of Ui(x,t) x at the coupling point x=ξj, for i=1,2,,npde and j=1,2,,nxi.
10: rcp[npde×nxi] const double Input
Note: the (i,j)th element of the matrix is stored in rcp[(j-1)×npde+i-1].
On entry: rcp[(j-1)×npde+i-1] contains the value of the flux Ri at the coupling point x=ξj, for i=1,2,,npde and j=1,2,,nxi.
11: ucpt[npde×nxi] const double Input
Note: the (i,j)th element of the matrix is stored in ucpt[(j-1)×npde+i-1].
On entry: if nxi>0, ucpt[(j-1)×npde+i-1] contains the value of Ui t at the coupling point x=ξj, for i=1,2,,npde and j=1,2,,nxi.
12: ucptx[npde×nxi] const double Input
Note: the (i,j)th element of the matrix is stored in ucptx[(j-1)×npde+i-1].
On entry: ucptx[(j-1)×npde+i-1] contains the value of 2Ui x t at the coupling point x=ξj, for i=1,2,,npde and j=1,2,,nxi.
13: f[nv] double Output
On exit: f[i-1] must contain the ith component of F, for i=1,2,,nv, where F is defined as
F=G-AV.-B ( Ut* Uxt* ) , (5)
or
F=-AV.-B ( Ut* Uxt* ) . (6)
The definition of F is determined by the input value of ires.
14: ires Integer * Input/Output
On entry: the form of F that must be returned in the array f.
ires=1
Equation (5) must be used.
ires=−1
Equation (6) must be used.
On exit: should usually remain unchanged. However, you may reset ires to force the integration function to take certain actions as described below:
ires=2
Indicates to the integrator that control should be passed back immediately to the calling function with the error indicator set to fail.code= NE_USER_STOP.
ires=3
Indicates to the integrator that the current time step should be abandoned and a smaller time step used instead. You may wish to set ires=3 when a physically meaningless input or output value has been generated. If you consecutively set ires=3, d03pjc returns to the calling function with the error indicator set to fail.code= NE_FAILED_DERIV.
15: comm Nag_Comm *
Pointer to structure of type Nag_Comm; the following members are relevant to odedef.
userdouble *
iuserInteger *
pPointer 
The type Pointer will be void *. Before calling d03pjc you may allocate memory and initialize these pointers with various quantities for use by odedef when called from d03pjc (see Section 3.1.1 in the Introduction to the NAG Library CL Interface).
Note: odedef should not return floating-point NaN (Not a Number) or infinity values, since these are not handled by d03pjc. If your code inadvertently does return any NaNs or infinities, d03pjc is likely to produce unexpected results.
15: nxi Integer Input
On entry: the number of ODE/PDE coupling points.
Constraints:
  • if nv=0, nxi=0;
  • if nv>0, nxi0.
16: xi[nxi] const double Input
On entry: xi[i-1], for i=1,2,,nxi, must be set to the ODE/PDE coupling points.
Constraint: xbkpts[0]xi[0]<xi[1]<<xi[nxi-1]xbkpts[nbkpts-1].
17: neqn Integer Input
On entry: the number of ODEs in the time direction.
Constraint: neqn=npde×npts+nv.
18: uvinit function, supplied by the user External Function
uvinit must compute the initial values of the PDE and the ODE components Ui(xj,t0), for i=1,2,,npde and j=1,2,,npts, and Vk(t0), for k=1,2,,nv.
The specification of uvinit is:
void  uvinit (Integer npde, Integer npts, const double x[], double u[], Integer nv, double v[], Nag_Comm *comm)
1: npde Integer Input
On entry: the number of PDEs in the system.
2: npts Integer Input
On entry: the number of mesh points in the interval [a,b].
3: x[npts] const double Input
On entry: x[i-1], for i=1,2,,npts, contains the current values of the space variable xi.
4: u[npde×npts] double Output
Note: the (i,j)th element of the matrix U is stored in u[(j-1)×npde+i-1].
On exit: u[(j-1)×npde+i-1] contains the value of the component Ui(xj,t0), for i=1,2,,npde and j=1,2,,npts.
5: nv Integer Input
On entry: the number of coupled ODEs in the system.
6: v[nv] double Output
On exit: v[i-1] contains the value of component Vi(t0), for i=1,2,,nv.
7: comm Nag_Comm *
Pointer to structure of type Nag_Comm; the following members are relevant to uvinit.
userdouble *
iuserInteger *
pPointer 
The type Pointer will be void *. Before calling d03pjc you may allocate memory and initialize these pointers with various quantities for use by uvinit when called from d03pjc (see Section 3.1.1 in the Introduction to the NAG Library CL Interface).
Note: uvinit should not return floating-point NaN (Not a Number) or infinity values, since these are not handled by d03pjc. If your code inadvertently does return any NaNs or infinities, d03pjc is likely to produce unexpected results.
19: rtol[dim] const double Input
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.
20: atol[dim] const double Input
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.
21: itol Integer Input
On entry: a value to indicate the form of the local error test. itol indicates to d03pjc whether to interpret either or both of rtol or atol as a vector or scalar. The error test to be satisfied is ei/wi<1.0, where wi is defined as follows:
itol rtol atol wi
1 scalar scalar rtol[0]×|Ui|+atol[0]
2 scalar vector rtol[0]×|Ui|+atol[i-1]
3 vector scalar rtol[i-1]×|Ui|+atol[0]
4 vector vector rtol[i-1]×|Ui|+atol[i-1]
In the above, ei denotes the estimated local error for the ith component of the coupled PDE/ODE system in time, u[i-1], for i=1,2,,neqn.
The choice of norm used is defined by the argument norm.
Constraint: 1itol4.
22: norm Nag_NormType Input
On entry: the type of norm to be used.
norm=Nag_MaxNorm
Maximum norm.
norm=Nag_TwoNorm
Averaged L2 norm.
If unorm denotes the norm of the vector u of length neqn, then for the averaged L2 norm
unorm=1neqni=1neqn(u[i-1]/wi)2,  
while for the maximum norm
u norm = maxi|u[i-1]/wi| .  
See the description of itol for the formulation of the weight vector w.
Constraint: norm=Nag_MaxNorm or Nag_TwoNorm.
23: laopt Nag_LinAlgOption Input
On entry: the type of matrix algebra required.
laopt=Nag_LinAlgFull
Full matrix methods to be used.
laopt=Nag_LinAlgBand
Banded matrix methods to be used.
laopt=Nag_LinAlgSparse
Sparse matrix methods to be used.
Constraint: laopt=Nag_LinAlgFull, Nag_LinAlgBand or Nag_LinAlgSparse.
Note: you are recommended to use the banded option when no coupled ODEs are present (i.e., nv=0).
24: algopt[30] const double Input
On entry: may be set to control various options available in the integrator. If you wish to employ all the default options, algopt[0] should be set to 0.0. Default values will also be used for any other elements of algopt set to zero. The permissible values, default values, and meanings are as follows:
algopt[0]
Selects the ODE integration method to be used. If algopt[0]=1.0, a BDF method is used and if algopt[0]=2.0, a Theta method is used. The default value is algopt[0]=1.0.
If algopt[0]=2.0, algopt[i-1], for i=2,3,4 are not used.
algopt[1]
Specifies the maximum order of the BDF integration formula to be used. algopt[1] may be 1.0, 2.0, 3.0, 4.0 or 5.0. The default value is algopt[1]=5.0.
algopt[2]
Specifies what method is to be used to solve the system of nonlinear equations arising on each step of the BDF method. If algopt[2]=1.0 a modified Newton iteration is used and if algopt[2]=2.0 a functional iteration method is used. If functional iteration is selected and the integrator encounters difficulty, there is an automatic switch to the modified Newton iteration. The default value is algopt[2]=1.0.
algopt[3]
Specifies whether or not the Petzold error test is to be employed. The Petzold error test results in extra overhead but is more suitable when algebraic equations are present, such as Pi,j=0.0, for j=1,2,,npde, for some i or when there is no V.i(t) dependence in the coupled ODE system. If algopt[3]=1.0, the Petzold test is used. If algopt[3]=2.0, the Petzold test is not used. The default value is algopt[3]=1.0.
If algopt[0]=1.0, algopt[i-1], for i=5,6,7, are not used.
algopt[4]
Specifies the value of Theta to be used in the Theta integration method. 0.51algopt[4]0.99. The default value is algopt[4]=0.55.
algopt[5]
Specifies what method is to be used to solve the system of nonlinear equations arising on each step of the Theta method. If algopt[5]=1.0, a modified Newton iteration is used and if algopt[5]=2.0, a functional iteration method is used. The default value is algopt[5]=1.0.
algopt[6]
Specifies whether or not the integrator is allowed to switch automatically between modified Newton and functional iteration methods in order to be more efficient. If algopt[6]=1.0, switching is allowed and if algopt[6]=2.0, switching is not allowed. The default value is algopt[6]=1.0.
algopt[10]
Specifies a point in the time direction, tcrit, beyond which integration must not be attempted. The use of tcrit is described under the argument itask. If algopt[0]0.0, a value of 0.0 for algopt[10], say, should be specified even if itask subsequently specifies that tcrit will not be used.
algopt[11]
Specifies the minimum absolute step size to be allowed in the time integration. If this option is not required, algopt[11] should be set to 0.0.
algopt[12]
Specifies the maximum absolute step size to be allowed in the time integration. If this option is not required, algopt[12] should be set to 0.0.
algopt[13]
Specifies the initial step size to be attempted by the integrator. If algopt[13]=0.0, the initial step size is calculated internally.
algopt[14]
Specifies the maximum number of steps to be attempted by the integrator in any one call. If algopt[14]=0.0, no limit is imposed.
algopt[22]
Specifies what method is to be used to solve the nonlinear equations at the initial point to initialize the values of U, Ut, V and V.. If algopt[22]=1.0, a modified Newton iteration is used and if algopt[22]=2.0, functional iteration is used. The default value is algopt[22]=1.0.
algopt[28] and algopt[29] are used only for the sparse matrix algebra option, laopt=Nag_LinAlgSparse.
algopt[28]
Governs the choice of pivots during the decomposition of the first Jacobian matrix. It should lie in the range 0.0<algopt[28]<1.0, with smaller values biasing the algorithm towards maintaining sparsity at the expense of numerical stability. If algopt[28] lies outside this range then the default value is used. If the functions regard the Jacobian matrix as numerically singular then increasing algopt[28] towards 1.0 may help, but at the cost of increased fill-in. The default value is algopt[28]=0.1.
algopt[29]
Is used as a relative pivot threshold during subsequent Jacobian decompositions (see algopt[28]) below which an internal error is invoked. If algopt[29] is greater than 1.0 no check is made on the pivot size, and this may be a necessary option if the Jacobian is found to be numerically singular (see algopt[28]). The default value is algopt[29]=0.0001.
25: rsave[lrsave] double Communication 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.
26: lrsave Integer Input
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, lrsave(3mlu+1)×neqn+nwkres+lenode.
If laopt=Nag_LinAlgSparse, lrsave4neqn+11neqn/2+1+nwkres+lenode.
Where mlu is the lower or upper half bandwidths such that
for PDE problems only (no coupled ODEs),
mlu=3npde-1;
for coupled PDE/ODE problems,
mlu=neqn-1.
Where nwkres is defined by
if nv>0​ and ​nxi>0,
nwkres=3 (npoly+1)2 + (npoly+1) × [npde2+6npde+nbkpts+1] + 8npde+nxi × (5npde+1) + nv+3 ;
if nv>0​ and ​nxi=0,
nwkres=3 (npoly+1)2 + (npoly+1) × [npde2+6npde+nbkpts+1] + 13npde+nv+4 ;
if nv=0,
nwkres=3 (npoly+1)2 + (npoly+1) × [npde2+6npde+nbkpts+1] + 13npde+5 .
Where lenode is defined by
if the BDF method is used,
lenode= (6+int(algopt[1])) × neqn+50 ;
if the Theta method is used,
lenode= 9neqn+50 .
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.
27: isave[lisave] Integer Communication Array
If ind=0, isave need not be set on entry.
If ind=1, isave must be unchanged from the previous call to the function because it contains required information about the iteration required for subsequent calls. In particular:
isave[0]
Contains the number of steps taken in time.
isave[1]
Contains the number of residual evaluations of the resulting ODE system used. One such evaluation involves computing the PDE functions at all the mesh points, as well as one evaluation of the functions in the boundary conditions.
isave[2]
Contains the number of Jacobian evaluations performed by the time integrator.
isave[3]
Contains the order of the ODE method last used in the time integration.
isave[4]
Contains the number of Newton iterations performed by the time integrator. Each iteration involves residual evaluation of the resulting ODE system followed by a back-substitution using the LU decomposition of the Jacobian matrix.
28: lisave Integer Input
On entry: the dimension of the array isave. Its size depends on the type of matrix algebra selected:
  • if laopt=Nag_LinAlgFull, lisave24;
  • if laopt=Nag_LinAlgBand, lisaveneqn+24;
  • if laopt=Nag_LinAlgSparse, lisave25×neqn+24.
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.
29: itask Integer Input
On entry: specifies the task to be performed by the ODE integrator.
itask=1
Normal computation of output values u at t=tout.
itask=2
One step and return.
itask=3
Stop at first internal integration point at or beyond t=tout.
itask=4
Normal computation of output values u at t=tout but without overshooting t=tcrit where tcrit is described under the argument algopt.
itask=5
Take one step in the time direction and return, without passing tcrit, where tcrit is described under the argument algopt.
Constraint: itask=1, 2, 3, 4 or 5.
30: itrace Integer Input
On entry: the level of trace information required from d03pjc 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, −1 is assumed and similarly if itrace>3, 3 is assumed.
The advisory messages are given in greater detail as itrace increases.
31: outfile const char * Input
On entry: the name of a file to which diagnostic output will be directed. If outfile is NULL the diagnostic output will be directed to standard output.
32: ind Integer * 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 argument tout should be reset between calls to d03pjc.
Constraint: ind=0 or 1.
On exit: ind=1.
33: comm Nag_Comm *
The NAG communication argument (see Section 3.1.1 in the Introduction to the NAG Library CL Interface).
34: saved Nag_D03_Save * Communication Structure
saved must remain unchanged following a previous call to a Chapter D03 function and prior to any subsequent call to a Chapter D03 function.
35: fail NagError * Input/Output
The NAG error argument (see Section 7 in the Introduction to the NAG Library CL Interface).

6 Error Indicators and Warnings

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

7 Accuracy

d03pjc 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 argument atol and rtol.

8 Parallelism and Performance

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

9 Further Comments

The argument specification allows you to include equations with only first-order derivatives in the space direction but there is no guarantee that the method of integration will be satisfactory for such systems. The position and nature of the boundary conditions in particular are critical in defining a stable problem.
The time taken depends on the complexity of the parabolic system and on the accuracy requested.

10 Example

This example provides a simple coupled system of one PDE and one ODE.
( V 1 ) 2 U 1 t -x V 1 V . 1 U 1 x = 2 U 1 x 2 V . 1 = V 1 U 1 + U 1 x +1 +t ,  
for t[10−4,0.1×2i],  i=1,2,,5,x[0,1].
The left boundary condition at x=0 is
U1 x =-V1expt.  
The right boundary condition at x=1 is
U1=-V1V.1.  
The initial conditions at t=10−4 are defined by the exact solution:
V1=t,   and  U1(x,t)=exp{t(1-x)}-1.0,  x[0,1],  
and the coupling point is at ξ1=1.0.

10.1 Program Text

Program Text (d03pjce.c)

10.2 Program Data

None.

10.3 Program Results

Program Results (d03pjce.r)
GnuplotProduced by GNUPLOT 5.0 patchlevel 0 Example Program Parabolic PDE Coupled with ODE using Collocation and BDF U(x,t) gnuplot_plot_1 gnuplot_plot_2 0.1 0.5 1 2 3 Time (logscale) 0 0.2 0.4 0.6 0.8 1 x −5 0 5 10 15 20 25