NAG CL Interface
d03plc (dim1_​parab_​convdiff_​dae)

1 Purpose

d03plc 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. 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>
void  d03plc (Integer npde, double *ts, double tout,
void (*pdedef)(Integer npde, double t, double x, const double u[], const double ux[], Integer nv, const double v[], const double vdot[], double p[], double c[], double d[], double s[], Integer *ires, Nag_Comm *comm),
void (*numflx)(Integer npde, double t, double x, Integer nv, 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 nv, const double v[], const double vdot[], Integer ibnd, double g[], Integer *ires, Nag_Comm *comm),
double u[], Integer npts, const 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 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[], 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: d03plc, nag_pde_dim1_parab_convdiff_dae or nag_pde_parab_1d_cd_ode.

3 Description

d03plc 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,,nv, (3)
where the vector V is the set of ODE solution values
Vt=V1t,,VnvtT,  
V. denotes its derivative with respect to time, and Ux is the spatial derivative of U.
In (1), 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,,RnvT, L is a vector of length nv, M is an nv by nv matrix, N is an nv 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. The initial values of the functions Ux,t and Vt must be given at t=t0.
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 d03pcc, d03phc and 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 vector, 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 d03plc 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=1 npde αk ek . (7)
An example is given in Section 10 and in the d03pfc documentation.
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 a separate numflx. For problems in the form (2)) the NAG defined null void function pointer, NULLFN, can be supplied in the call to d03plc.
The boundary condition specification has sufficient flexibility to allow for different types of problems. For second-order problems, i.e., Di depending on Ux, a boundary condition is required for each PDE at both boundaries for the problem to be well-posed. If there are no second-order 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 condition, 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 Section 10 and in the d03pfc documentation.
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.
The problem is subject to the following restrictions:
  1. (i)In (1), V.jt, for j=1,2,,nv, may only appear linearly in the functions Si, for i=1,2,,npde, with a similar restriction for GiL and GiR;
  2. (ii)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;
  3. (iii)t0<tout, so that integration is in the forward direction;
  4. (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 mesh points x1,x2,,xnpts;
  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.
In total there are npde×npts+nv 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)).
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
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
Sod G A (1978) A survey of several finite difference methods for systems of nonlinear hyperbolic conservation laws J. Comput. Phys. 27 1–31

5 Arguments

1: npde Integer Input
On entry: the number of PDEs to be solved.
Constraint: npde1.
2: 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.
3: tout double Input
On entry: the final value of t to which the integration is to be carried out.
4: pdedef function, supplied by the user External 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 d03plc. For problems in the form (2)) the NAG defined null void function pointer, NULLFN, can be supplied in the call to d03plc.
The specification of pdedef is:
void  pdedef (Integer npde, double t, double x, const double u[], const double ux[], Integer nv, const double v[], const double vdot[], double p[], double c[], double d[], double s[], 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 double Input
On entry: the current value of the space variable x.
4: u[npde] const double Input
On entry: u[i-1] contains the value of the component Uix,t, for i=1,2,,npde.
5: ux[npde] const double Input
On entry: ux[i-1] contains the value of the component Uix,t x , for i=1,2,,npde.
6: nv Integer Input
On entry: the number of coupled ODEs in the system.
7: v[nv] const double Input
On entry: if nv>0, v[i-1] contains the value of the component Vit, for i=1,2,,nv.
8: vdot[nv] const double Input
On entry: if nv>0, vdot[i-1] contains the value of component V.it, for i=1,2,,nv.
Note: V.it, for i=1,2,,nv, may only appear linearly in Sj, for j=1,2,,npde.
9: p[npde×npde] double Output
Note: the i,jth element of the matrix P is stored in p[j-1×npde+i-1].
On exit: p[j-1×npde+i-1] must be set to the value of Pi,jx,t,U,V, for i=1,2,,npde and j=1,2,,npde.
10: c[npde] double Output
On exit: c[i-1] must be set to the value of Cix,t,U,V, for i=1,2,,npde.
11: d[npde] double Output
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] double Output
On exit: s[i-1] must be set to the value of Six,t,U,V,V., for i=1,2,,npde.
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, d03plc 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 d03plc you may allocate memory and initialize these pointers with various quantities for use by pdedef when called from d03plc (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 d03plc. If your code inadvertently does return any NaNs or infinities, d03plc is likely to produce unexpected results.
5: numflx function, supplied by the user External 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 d03plc.
The specification of numflx is:
void  numflx (Integer npde, double t, double x, Integer nv, const double v[], const double uleft[], const double uright[], double flux[], Integer *ires, Nag_Comm *comm, Nag_D03_Save *saved)
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 double Input
On entry: the current value of the space variable x.
4: nv Integer Input
On entry: the number of coupled ODEs in the system.
5: v[nv] const double Input
On entry: if nv>0, v[i-1] contains the value of the component Vit, for i=1,2,,nv.
6: uleft[npde] const double Input
On entry: uleft[i-1] contains the left value of the component Uix, for i=1,2,,npde.
7: uright[npde] const double Input
On entry: uright[i-1] contains the right value of the component Uix, for i=1,2,,npde.
8: flux[npde] double Output
On exit: flux[i-1] must be set to the numerical flux F^i, for i=1,2,,npde.
9: 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, d03plc returns to the calling function with the error indicator set to fail.code= NE_FAILED_DERIV.
10: comm Nag_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 d03plc you may allocate memory and initialize these pointers with various quantities for use by numflx when called from d03plc (see Section 3.1.1 in the Introduction to the NAG Library CL Interface).
11: saved Nag_D03_Save * Communication Structure
If numflx calls one of the approximate Riemann solvers d03puc, d03pvc, d03pwc or d03pxc then saved is used to pass data concerning the computation to the solver. You should not change the components of saved.
Note: numflx should not return floating-point NaN (Not a Number) or infinity values, since these are not handled by d03plc. If your code inadvertently does return any NaNs or infinities, d03plc is likely to produce unexpected results.
6: bndary function, supplied by the user External 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 nv, const double v[], const double vdot[], Integer ibnd, double g[], Integer *ires, 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: t double Input
On entry: the current value of the independent variable t.
4: x[npts] const double Input
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 double Input
Note: the i,jth element of the matrix U is stored in u[j-1×npde+i-1].
On entry: u[j-1×npde+i-1] contains the value of 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: nv Integer Input
On entry: the number of coupled ODEs in the system.
7: v[nv] const double Input
On entry: if nv>0, v[i-1] contains the value of the component Vit, for i=1,2,,nv.
8: vdot[nv] const double Input
On entry: if nv>0, vdot[i-1] contains the value of component V.it, for i=1,2,,nv.
Note: V.it, for i=1,2,,nv, may only appear linearly in GjL and GjR, for j=1,2,,npde.
9: ibnd Integer Input
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] double Output
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: 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, d03plc 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 d03plc you may allocate memory and initialize these pointers with various quantities for use by bndary when called from d03plc (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 d03plc. If your code inadvertently does return any NaNs or infinities, d03plc is likely to produce unexpected results.
7: u[neqn] double Input/Output
On entry: the initial values of the dependent variables defined as follows:
  • u[npde×j-1+i-1] contain Uixj,t0, for i=1,2,,npde and j=1,2,,npts, and
  • u[npts×npde+i-1] contain Vit0, for i=1,2,,nv.
On exit: the computed solution Uixj,t, for i=1,2,,npde and j=1,2,,npts, and Vkt, for k=1,2,,nv, evaluated at t=ts, as follows:
  • u[npde×j-1+i-1] contain Uixj,t, for i=1,2,,npde and j=1,2,,npts, and
  • u[npts×npde+i-1] contain Vit, for i=1,2,,nv.
8: npts Integer Input
On entry: the number of mesh points in the interval a,b.
Constraint: npts3.
9: x[npts] const double Input
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].
10: nv Integer Input
On entry: the number of coupled ODE components.
Constraint: nv0.
11: odedef function, supplied by the user External Function
odedef must evaluate the functions R, which define the system of ODEs, as given in (4).
If nv=0, odedef will never be called and the NAG defined null void function pointer, NULLFN, can be supplied in the call to d03plc.
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 ucpt[], 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: 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 Vit, 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.it, 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 point, ξi, for i=1,2,,nxi.
8: ucp[npde×nxi] const double Input
Note: the i,jth element of the matrix is stored in ucp[j-1×npde+i-1].
On entry: if nxi>0, ucp[j-1×npde+i-1] contains the value of Uix,t at the coupling point x=ξj, for i=1,2,,npde and j=1,2,,nxi.
9: ucpx[npde×nxi] const double Input
Note: the i,jth element of the matrix is stored in ucpx[j-1×npde+i-1].
On entry: if nxi>0, ucpx[j-1×npde+i-1] contains the value of Uix,t x at the coupling point x=ξj, for i=1,2,,npde and j=1,2,,nxi.
10: ucpt[npde×nxi] const double Input
Note: the i,jth element of the matrix is stored in ucpt[j-1×npde+i-1].
On entry: if nxi>0, ucpt[j-1×npde+i-1] contains the value of Ui t at the coupling point x=ξj, for i=1,2,,npde and j=1,2,,nxi.
11: r[nv] double Output
On exit: r[i-1] must contain the ith component of R, for i=1,2,,nv, 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: ires Integer * 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, d03plc returns to the calling function with the error indicator set to fail.code= NE_FAILED_DERIV.
13: comm Nag_Comm *
Pointer to structure of type Nag_Comm; the following members are relevant to odedef.
userdouble *
iuserInteger *
pPointer 
The type Pointer will be void *. Before calling d03plc you may allocate memory and initialize these pointers with various quantities for use by odedef when called from d03plc (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 d03plc. If your code inadvertently does return any NaNs or infinities, d03plc is likely to produce unexpected results.
12: nxi Integer Input
On entry: the number of ODE/PDE coupling points.
Constraints:
  • if nv=0, nxi=0;
  • if nv>0, nxi0.
13: 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: x[0]xi[0]<xi[1]<<xi[nxi-1]x[npts-1].
14: neqn Integer Input
On entry: the number of ODEs in the time direction.
Constraint: neqn=npde×npts+nv.
15: 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.
16: 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.
17: itol Integer Input
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, the error test to be satisfied is ei<1.0. itol indicates to d03plc 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: 1itol4.
18: norm Nag_NormType Input
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.
19: 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 (nv=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.
20: 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 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, there is an automatic switch to the modified Newton iteration. The default value is algopt[2]=1.0.
algopt[3]
Specifies whether or not the Petzold error test is to be employed. The Petzold error test results in extra overhead but is more suitable when algebraic equations are present, such as Pi,j=0.0, for j=1,2,,npde, for some i or when there is no V.it dependence in the coupled ODE system. If algopt[3]=1.0, the Petzold test is used. If algopt[3]=2.0, the Petzold test is not used. The default value is algopt[3]=1.0.
If algopt[0]=1.0, algopt[i-1], for i=5,6,7, are not used.
algopt[4]
Specifies the value of Theta to be used in the Theta integration method. 0.51algopt[4]0.99. The default value is algopt[4]=0.55.
algopt[5]
Specifies what method is to be used to solve the system of nonlinear equations arising on each step of the Theta method. If algopt[5]=1.0, a modified Newton iteration is used and if algopt[5]=2.0, a functional iteration method is used. The default value is algopt[5]=1.0.
algopt[6]
Specifies whether or not the integrator is allowed to switch automatically between modified Newton and functional iteration methods in order to be more efficient. If algopt[6]=1.0, switching is allowed and if algopt[6]=2.0, switching is not allowed. The default value is algopt[6]=1.0.
algopt[10]
Specifies a point in the time direction, tcrit, beyond which integration must not be attempted. The use of tcrit is described under the argument itask. If algopt[0]0.0, a value of 0.0 for algopt[10], say, should be specified even if itask subsequently specifies that tcrit will not be used.
algopt[11]
Specifies the minimum absolute step size to be allowed in the time integration. If this option is not required, algopt[11] should be set to 0.0.
algopt[12]
Specifies the maximum absolute step size to be allowed in the time integration. If this option is not required, algopt[12] should be set to 0.0.
algopt[13]
Specifies the initial step size to be attempted by the integrator. If algopt[13]=0.0, the initial step size is calculated internally.
algopt[14]
Specifies the maximum number of steps to be attempted by the integrator in any one call. If algopt[14]=0.0, no limit is imposed.
algopt[22]
Specifies what method is to be used to solve the nonlinear equations at the initial point to initialize the values of U, Ut, V and V.. If algopt[22]=1.0, a modified Newton iteration is used and if algopt[22]=2.0, functional iteration is used. The default value is algopt[22]=1.0.
algopt[28] and algopt[29] are used only for the sparse matrix algebra option, 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, 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.
21: 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.
22: 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, lrsave3mlu+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=npde2npts+6nxi+3npde+26+nxi+nv+7npts+2;
if nv>0​ and ​nxi=0,
nwkres=npde2npts+3npde+32+nv+7npts+3;
if nv=0,
nwkres=npde2npts+3npde+32+7npts+4.
Where lenode is defined by
if the BDF method is used,
lenode=6+intalgopt[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.
23: isave[lisave] Integer Communication 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.
24: 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.
25: itask Integer Input
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.
26: itrace Integer Input
On entry: the level of trace information required from d03plc 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.
27: 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.
28: 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 d03plc.
Constraint: ind=0 or 1.
On exit: ind=1.
29: comm Nag_Comm *
The NAG communication argument (see Section 3.1.1 in the Introduction to the NAG Library CL Interface).
30: 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.
31: 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, 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.
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. Incorrect specification of boundary conditions may also result in this error.
Underlying ODE solver cannot make further progress from the point ts with the supplied values of atol and rtol. ts=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, itask=value.
Constraint: itask=1, 2, 3, 4 or 5.
On entry, itol=value.
Constraint: itol=1, 2, 3 or 4.
On entry, npde=value.
Constraint: npde1.
On entry, npts=value.
Constraint: npts3.
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_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, x[i-1]=value, j=value and x[j-1]=value.
Constraint: x[0]<x[1]<<x[npts-1].
On entry, i=value, xi[i]=value and xi[i-1]=value.
Constraint: xi[i]>xi[i-1].
NE_NOT_WRITE_FILE
Cannot open file value for writing.
NE_REAL_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_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.
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

d03plc 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

d03plc is threaded by NAG for parallel execution in multithreaded implementations of the NAG Library.
d03plc 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

d03plc 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 schemes 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 and on the accuracy requested. For a given system and a fixed accuracy 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 first-order system with coupled ODEs arising from the use of the characteristic equations for the numerical boundary conditions.
The PDEs are
U1 t + U1 x +2 U2 x =0, U2 t +2 U1 x + U2 x =0,  
for x0,1 and t0.
The PDEs have an exact solution given by
U1x,t=fx-3t+gx+t, U2x,t=fx-3t-gx+t,  
where fz=expπzsin2πz, gz=exp-2πzcos2πz.
The initial conditions are given by the exact solution.
The characteristic variables are W1=U1-U2 and W2=U1+U2, corresponding to the characteristics given by dx/dt=-1 and dx/dt=3 respectively. Hence we require a physical boundary condition for W2 at the left-hand boundary and for W1 at the right-hand boundary (corresponding to the incoming characteristics), and a numerical boundary condition for W1 at the left-hand boundary and for W2 at the right-hand boundary (outgoing characteristics).
The physical boundary conditions are obtained from the exact solution, and the numerical boundary conditions are supplied in the form of the characteristic equations for the outgoing characteristics, that is
W1 t - W1 x =0  
at the left-hand boundary, and
W2 t + 3 W2 x =0  
at the right-hand boundary.
In order to specify these boundary conditions, two ODE variables V1 and V2 are introduced, defined by
V1t=W10,t=U10,t-U20,t, V2t=W21,t=U11,t+U21,t.  
The coupling points are therefore at x=0 and x=1.
The numerical boundary conditions are now
V.1- W1 x =0  
at the left-hand boundary, and
V.2+ 3 W2 x =0  
at the right-hand boundary.
The spatial derivatives are evaluated at the appropriate boundary points in bndary using one-sided differences (into the domain and therefore consistent with the characteristic directions).
The numerical flux is calculated using Roe's approximate Riemann solver (see Section 3 for details), giving
F^=12 3U1L-U1R+3U2L+U2R 3U1L+U1R+3U2L-U2R .  
Example 2 (ex2)
This example is the standard shock-tube test problem proposed by Sod (1978) for the Euler equations of gas dynamics. The problem models the flow of a gas in a long tube following the sudden breakdown of a diaphragm separating two initial gas states at different pressures and densities. There is an exact solution to this problem which is not included explicitly as the calculation is quite lengthy. The PDEs are
ρ t + m x = 0, m t + x m2ρ+γ- 1 e-m22ρ = 0, e t + x m eρ+mργ- 1 e-m22ρ = 0,  
where ρ is the density; m is the momentum, such that m=ρu, where u is the velocity; e is the specific energy; and γ is the (constant) ratio of specific heats. The pressure p is given by
p=γ-1 e-ρu22 .  
The solution domain is 0x1 for 0<t0.2, with the initial discontinuity at x=0.5, and initial conditions
ρx,0=1, mx,0=0, ex,0=2.5,   for ​x<0.5, ρx,0=0.125, mx,0=0, ex,0=0.25,   for ​x>0.5.  
The solution is uniform and constant at both boundaries for the spatial domain and time of integration stated, and hence the physical and numerical boundary conditions are indistinguishable and are both given by the initial conditions above. The evaluation of the numerical flux for the Euler equations is not trivial; the Roe algorithm given in Section 3 cannot be used directly as the Jacobian is nonlinear. However, an algorithm is available using the argument-vector method (see Roe (1981)), and this is provided in the utility function d03puc. An alternative Approxiate Riemann Solver using Osher's scheme is provided in d03pvc. Either d03puc or d03pvc can be called from numflx.

10.1 Program Text

Program Text (d03plce.c)

10.2 Program Data

None.

10.3 Program Results

Program Results (d03plce.r)
GnuplotProduced by GNUPLOT 5.0 patchlevel 0 Example Program 1 First-order System with Coupled ODEs Solution U(1,x,t) U(1,x,t) gnuplot_plot_1 gnuplot_plot_2 0 0.1 0.2 0.3 0.4 0.5 0.6 Time 0 0.2 0.4 0.6 0.8 1 x −8 −4 0
GnuplotProduced by GNUPLOT 5.0 patchlevel 0 Example Program 2 Shock Tube Test Problem of Euler Equations in Gas Dynamics DENSITY Density gnuplot_plot_1 gnuplot_plot_2 0 0.02 0.04 0.06 0.08 0.1 0.12 0.14 0.16 0.18 0.2 0.22 Time 0 0.1 0.2 0.3 0.4 0.5 0.6 0.7 0.8 0.9 1 x 0.2 0.4 0.6 0.8 1
GnuplotProduced by GNUPLOT 5.0 patchlevel 0 Shock Tube Test Problem of Euler Equations in Gas Dynamics VELOCITY Velocity gnuplot_plot_1 gnuplot_plot_2 0 0.02 0.04 0.06 0.08 0.1 0.12 0.14 0.16 0.18 0.2 0.22 Time 0 0.1 0.2 0.3 0.4 0.5 0.6 0.7 0.8 0.9 1 x 0 0.3 0.6 0.9
GnuplotProduced by GNUPLOT 5.0 patchlevel 0 Shock Tube Test Problem of Euler Equations in Gas Dynamics PRESSURE Pressure gnuplot_plot_1 gnuplot_plot_2 0 0.02 0.04 0.06 0.08 0.1 0.12 0.14 0.16 0.18 0.2 0.22 Time 0 0.1 0.2 0.3 0.4 0.5 0.6 0.7 0.8 0.9 1 x 0.1 0.2 0.3 0.4 0.5 0.6 0.7 0.8 0.9 1