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

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

NAG Toolbox: nag_pde_1d_parab_convdiff (d03pf)

Purpose

nag_pde_1d_parab_convdiff (d03pf) integrates a system of linear or nonlinear convection-diffusion equations in one space dimension, with optional source terms. 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 PDEs to a system of ordinary differential equations (ODEs), and the resulting system is solved using a backward differentiation formula (BDF) method.

Syntax

[ts, u, rsave, isave, ind, ifail] = d03pf(ts, tout, pdedef, numflx, bndary, u, x, acc, tsmax, rsave, isave, itask, itrace, ind, 'npde', npde, 'npts', npts)
[ts, u, rsave, isave, ind, ifail] = nag_pde_1d_parab_convdiff(ts, tout, pdedef, numflx, bndary, u, x, acc, tsmax, rsave, isave, itask, itrace, ind, 'npde', npde, 'npts', npts)
Note: the interface to this routine has changed since earlier releases of the toolbox:
Mark 22: lrsave, lisave have been removed from the interface
.

Description

nag_pde_1d_parab_convdiff (d03pf) integrates the system of convection-diffusion equations in conservative form:
npde
Pi,j(Uj)/(t) + (Fi)/(x) = Ci(Di)/(x) + Si,
j = 1
j=1npdePi,j Uj t + Fi x =Ci Di x +Si,
(1)
or the hyperbolic convection-only system:
(Ui)/(t) + (Fi)/(x) = 0,
Ui t + Fi x =0,
(2)
for i = 1,2,,npde,  axb,  tt0i=1,2,,npde,  axb,  tt0, where the vector UU is the set of solution values
U (x,t) = [ U1 (x,t) ,, Unpde (x,t) ]T .
U (x,t) = [ U 1 (x,t) ,, U npde (x,t) ] T .
The functions Pi,jPi,j, FiFi, CiCi and SiSi depend on xx, tt and UU; and DiDi depends on xx, tt, UU and UxUx, where UxUx is the spatial derivative of UU. Note that Pi,jPi,j, FiFi, CiCi and SiSi must not depend on any space derivatives; and none of the functions may depend on time derivatives. In terms of conservation laws, FiFi, (CiDi)/(x) CiDi x  and SiSi are the convective flux, diffusion and source terms respectively.
The integration in time is from t0t0 to touttout, over the space interval axbaxb, where a = x1a=x1 and b = xnptsb=xnpts are the leftmost and rightmost points of a user-defined mesh x1,x2,,xnptsx1,x2,,xnpts. The initial values of the functions U(x,t)U(x,t) must be given at t = t0t=t0.
The PDEs are approximated by a system of ODEs in time for the values of UiUi at mesh points using a spatial discretization method similar to the central-difference scheme used in nag_pde_1d_parab_fd (d03pc), nag_pde_1d_parab_dae_fd (d03ph) and nag_pde_1d_parab_remesh_fd (d03pp), but with the flux FiFi 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, iF^i say, must be calculated by you in terms of the left and right values of the solution vector UU (denoted by ULUL and URUR respectively), at each mid-point of the mesh xj1 / 2 = (xj1 + xj) / 2xj-1/2=(xj-1+xj)/2, for j = 2,3,,nptsj=2,3,,npts. The left and right values are calculated by nag_pde_1d_parab_convdiff (d03pf) 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 iF^i is derived from the solution of the Riemann problem given by
(Ui)/(t) + (Fi)/(y) = 0,
Ui t + Fi y =0,
(3)
where y = xxj1 / 2y=x-xj-1/2, i.e., y = 0y=0 corresponds to x = xj1 / 2x=xj-1/2, with discontinuous initial values U = ULU=UL for y < 0y<0 and U = URU=UR for y > 0y>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,jPi,j, CiCi, DiDi and SiSi. 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 = 0Ut+Fx=0 or equivalently Ut + AUx = 0Ut+AUx=0. Provided the system is linear in UU, i.e., the Jacobian matrix AA does not depend on UU, the numerical flux F^ is given by
npde
= (1/2)(FL + FR)(1/2)αk|λk|ek,
k = 1
F^=12 (FL+FR)-12k=1npdeαk|λk|ek,
(4)
where FLFL (FRFR) is the flux FF calculated at the left (right) value of UU, denoted by ULUL (URUR); the λkλk are the eigenvalues of AA; the ekek are the right eigenvectors of AA; and the αkαk are defined by
npde
URUL = αkek.
k = 1
UR-UL=k=1npdeαkek.
(5)
An example is given in Section [Example].
If the system is nonlinear, Roe's scheme requires that a linearized Jacobian is found (see Roe (1981)).
The functions Pi,jPi,j, CiCi, DiDi and SiSi (but not FiFi) must be specified in a pdedef. The numerical flux iF^i must be supplied in a separate numflx. For problems in the form (2), the actual argument nag_pde_1d_parab_convdiff_sample_pdedef (d03pfp) may be used for pdedef. nag_pde_1d_parab_convdiff_sample_pdedef (d03pfp) is included in the NAG Toolbox and sets the matrix with entries Pi,jPi,j to the identity matrix, and the functions CiCi, DiDi and SiSi to zero.
The boundary condition specification has sufficient flexibility to allow for different types of problems. For second-order problems, i.e., DiDi depending on UxUx, 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 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 only linear extrapolation is allowed in this function (for greater flexibility the function nag_pde_1d_parab_convdiff_dae (d03pl) should be used). 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. Examples can be found in Section [Example].
The boundary conditions must be specified in bndary in the form
GiL(x,t,U) = 0  at ​x = a,  i = 1,2,,npde,
GiL(x,t,U)=0  at ​x=a,  i=1,2,,npde,
(6)
at the left-hand boundary, and
GiR(x,t,U) = 0  at ​x = b,  i = 1,2,,npde,
GiR(x,t,U)=0  at ​x=b,  i=1,2,,npde,
(7)
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 UU 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 problem is subject to the following restrictions:
(i) Pi,jPi,j, FiFi, CiCi and SiSi must not depend on any space derivatives;
(ii) Pi,jPi,j, FiFi, CiCi, DiDi and SiSi must not depend on any time derivatives;
(iii) t0 < toutt0<tout, so that integration is in the forward direction;
(iv) The evaluation of the terms Pi,jPi,j, CiCi, DiDi and SiSi 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,,xnptsx1,x2,,xnpts;
(v) At least one of the functions Pi,jPi,j must be nonzero so that there is a time derivative present in the PDE problem.
In total there are npde × nptsnpde×npts ODEs in the time direction. This system is then integrated forwards in time using a BDF method.
For further details of the algorithm, see Pennington and Berzins (1994) and the references therein.

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

Parameters

Compulsory Input Parameters

1:     ts – double scalar
The initial value of the independent variable tt.
Constraint: ts < toutts<tout.
2:     tout – double scalar
The final value of tt to which the integration is to be carried out.
3:     pdedef – function handle or string containing name of m-file
pdedef must evaluate the functions Pi,jPi,j, CiCi, DiDi and SiSi which partially define the system of PDEs. Pi,jPi,j, CiCi and SiSi may depend on xx, tt and UU; DiDi may depend on xx, tt, UU and UxUx. pdedef is called approximately midway between each pair of mesh points in turn by nag_pde_1d_parab_convdiff (d03pf). The actual argument nag_pde_1d_parab_convdiff_sample_pdedef (d03pfp) may be used for pdedef for problems in the form (2). (nag_pde_1d_parab_convdiff_sample_pdedef (d03pfp) is included in the NAG Toolbox.)
[p, c, d, s, ires] = pdedef(npde, t, x, u, ux, ires)

Input Parameters

1:     npde – int64int32nag_int scalar
The number of PDEs in the system.
2:     t – double scalar
The current value of the independent variable tt.
3:     x – double scalar
The current value of the space variable xx.
4:     u(npde) – double array
u(i)ui contains the value of the component Ui(x,t)Ui(x,t), for i = 1,2,,npdei=1,2,,npde.
5:     ux(npde) – double array
ux(i)uxi contains the value of the component (Ui(x,t))/(x) Ui(x,t) x , for i = 1,2,,npdei=1,2,,npde.
6:     ires – int64int32nag_int scalar
Set to 1​ or ​1-1​ or ​1.

Output Parameters

1:     p(npde,npde) – double array
p(i,j)pij must be set to the value of Pi,j(x,t,U)Pi,j(x,t,U), for i = 1,2,,npdei=1,2,,npde and j = 1,2,,npdej=1,2,,npde.
2:     c(npde) – double array
c(i)ci must be set to the value of Ci(x,t,U)Ci(x,t,U), for i = 1,2,,npdei=1,2,,npde.
3:     d(npde) – double array
d(i)di must be set to the value of Di(x,t,U,Ux)Di(x,t,U,Ux), for i = 1,2,,npdei=1,2,,npde.
4:     s(npde) – double array
s(i)si must be set to the value of Si(x,t,U)Si(x,t,U), for i = 1,2,,npdei=1,2,,npde.
5:     ires – int64int32nag_int scalar
Should usually remain unchanged. However, you may set ires to force the integration function to take certain actions as described below:
ires = 2ires=2
Indicates to the integrator that control should be passed back immediately to the calling function with the error indicator set to ifail = 6ifail=6.
ires = 3ires=3
Indicates to the integrator that the current time step should be abandoned and a smaller time step used instead. You may wish to set ires = 3ires=3 when a physically meaningless input or output value has been generated. If you consecutively set ires = 3ires=3, then nag_pde_1d_parab_convdiff (d03pf) returns to the calling function with the error indicator set to ifail = 4ifail=4.
4:     numflx – function handle or string containing name of m-file
numflx must supply the numerical flux for each PDE given the left and right values of the solution vector uu. numflx is called approximately midway between each pair of mesh points in turn by nag_pde_1d_parab_convdiff (d03pf).
[flux, ires] = numflx(npde, t, x, uleft, uright, ires)

Input Parameters

1:     npde – int64int32nag_int scalar
The number of PDEs in the system.
2:     t – double scalar
The current value of the independent variable tt.
3:     x – double scalar
The current value of the space variable xx.
4:     uleft(npde) – double array
uleft(i)ulefti contains the left value of the component Ui(x)Ui(x), for i = 1,2,,npdei=1,2,,npde.
5:     uright(npde) – double array
uright(i)urighti contains the right value of the component Ui(x)Ui(x), for i = 1,2,,npdei=1,2,,npde.
6:     ires – int64int32nag_int scalar
Set to 1​ or ​1-1​ or ​1.

Output Parameters

1:     flux(npde) – double array
flux(i)fluxi must be set to the numerical flux iF^i, for i = 1,2,,npdei=1,2,,npde.
2:     ires – int64int32nag_int scalar
Should usually remain unchanged. However, you may set ires to force the integration function to take certain actions as described below:
ires = 2ires=2
Indicates to the integrator that control should be passed back immediately to the calling function with the error indicator set to ifail = 6ifail=6.
ires = 3ires=3
Indicates to the integrator that the current time step should be abandoned and a smaller time step used instead. You may wish to set ires = 3ires=3 when a physically meaningless input or output value has been generated. If you consecutively set ires = 3ires=3, then nag_pde_1d_parab_convdiff (d03pf) returns to the calling function with the error indicator set to ifail = 4ifail=4.
5:     bndary – function handle or string containing name of m-file
bndary must evaluate the functions GiLGiL and GiRGiR which describe the physical and numerical boundary conditions, as given by (6) and (7).
[g, ires] = bndary(npde, npts, t, x, u, ibnd, ires)

Input Parameters

1:     npde – int64int32nag_int scalar
The number of PDEs in the system.
2:     npts – int64int32nag_int scalar
The number of mesh points in the interval [a,b][a,b].
3:     t – double scalar
The current value of the independent variable tt.
4:     x(npts) – double array
The mesh points in the spatial direction. x(1)x1 corresponds to the left-hand boundary, aa, and x(npts)xnpts corresponds to the right-hand boundary, bb.
5:     u(npde,33) – double array
Contains the value of solution components in the boundary region.
If ibnd = 0ibnd=0, u(i,j)uij contains the value of the component Ui(x ,t)Ui(x ,t) at x = x(j)x=xj, for i = 1,2,,npdei=1,2,,npde and j = 1,2,3j=1,2,3.
If ibnd0ibnd0, u(i,j)uij contains the value of the component Ui(x,t)Ui(x,t) at x = x(nptsj + 1)x=xnpts-j+1, for i = 1,2,,npdei=1,2,,npde and j = 1,2,3j=1,2,3.
6:     ibnd – int64int32nag_int scalar
Specifies which boundary conditions are to be evaluated.
ibnd = 0ibnd=0
bndary must evaluate the left-hand boundary condition at x = ax=a.
ibnd0ibnd0
bndary must evaluate the right-hand boundary condition at x = bx=b.
7:     ires – int64int32nag_int scalar
Set to 1​ or ​1-1​ or ​1.

Output Parameters

1:     g(npde) – double array
g(i)gi must contain the iith component of either gLgL or gRgR in (6) and (7), depending on the value of ibnd, for i = 1,2,,npdei=1,2,,npde.
2:     ires – int64int32nag_int scalar
Should usually remain unchanged. However, you may set ires to force the integration function to take certain actions as described below:
ires = 2ires=2
Indicates to the integrator that control should be passed back immediately to the calling function with the error indicator set to ifail = 6ifail=6.
ires = 3ires=3
Indicates to the integrator that the current time step should be abandoned and a smaller time step used instead. You may wish to set ires = 3ires=3 when a physically meaningless input or output value has been generated. If you consecutively set ires = 3ires=3, then nag_pde_1d_parab_convdiff (d03pf) returns to the calling function with the error indicator set to ifail = 4ifail=4.
6:     u(npde,npts) – double array
npde, the first dimension of the array, must satisfy the constraint npde1npde1.
u(i,j)uij must contain the initial value of Ui(x,t)Ui(x,t) at x = x(j)x=xj and t = tst=ts, for i = 1,2,,npdei=1,2,,npde and j = 1,2,,nptsj=1,2,,npts.
7:     x(npts) – double array
npts, the dimension of the array, must satisfy the constraint npts3npts3.
The mesh points in the space direction. x(1)x1 must specify the left-hand boundary, aa, and x(npts)xnpts must specify the right-hand boundary, bb.
Constraint: x(1) < x(2) < < x(npts)x1<x2<<xnpts.
8:     acc(22) – double array
The components of acc contain the relative and absolute error tolerances used in the local error test in the time integration.
If E(i,j)E(i,j) is the estimated error for UiUi at the jjth mesh point, the error test is
E(i,j) = acc(1) × u(i,j) + acc(2).
E(i,j)=acc1×uij+acc2.
Constraint: acc(1)acc1 and acc(2)0.0acc20.0 (but not both zero).
9:     tsmax – double scalar
The maximum absolute step size to be allowed in the time integration. If tsmax = 0.0tsmax=0.0 then no maximum is imposed.
Constraint: tsmax0.0tsmax0.0.
10:   rsave(lrsave) – double array
lrsave, the dimension of the array, must satisfy the constraint lrsave(11 + 9 × npde) × npde × npts + (32 + 3 × npde) × npde + 7 ×  npts + 54lrsave(11+9×npde)×npde×npts+(32+3×npde)×npde+7×npts+54.
If ind = 0ind=0, rsave need not be set on entry.
If ind = 1ind=1, rsave must be unchanged from the previous call to the function because it contains required information about the iteration.
11:   isave(lisave) – int64int32nag_int array
lisave, the dimension of the array, must satisfy the constraint lisavenpde × npts + 24lisavenpde×npts+24.
If ind = 0ind=0, isave need not be set on entry.
If ind = 1ind=1, isave must be unchanged from the previous call to the function because it contains required information about the iteration. In particular:
isave(1)isave1
Contains the number of steps taken in time.
isave(2)isave2
Contains the number of residual evaluations of the resulting ODE system used. One such evaluation involves computing the PDE functions at all the mesh points, as well as one evaluation of the functions in the boundary conditions.
isave(3)isave3
Contains the number of Jacobian evaluations performed by the time integrator.
isave(4)isave4
Contains the order of the last backward differentiation formula method used.
isave(5)isave5
Contains the number of Newton iterations performed by the time integrator. Each iteration involves an ODE residual evaluation followed by a back-substitution using the LULU decomposition of the Jacobian matrix.
12:   itask – int64int32nag_int scalar
The task to be performed by the ODE integrator.
itask = 1itask=1
Normal computation of output values uu at t = toutt=tout (by overshooting and interpolating).
itask = 2itask=2
Take one step in the time direction and return.
itask = 3itask=3
Stop at first internal integration point at or beyond t = toutt=tout.
Constraint: itask = 1itask=1, 22 or 33.
13:   itrace – int64int32nag_int scalar
The level of trace information required from nag_pde_1d_parab_convdiff (d03pf) and the underlying ODE solver. itrace may take the value 1-1, 00, 11, 22 or 33.
itrace = 1itrace=-1
No output is generated.
itrace = 0itrace=0
Only warning messages from the PDE solver are printed on the current error message unit (see nag_file_set_unit_error (x04aa)).
itrace > 0itrace>0
Output from the underlying ODE solver is printed on the current advisory message unit (see nag_file_set_unit_advisory (x04ab)). This output contains details of Jacobian entries, the nonlinear iteration and the time integration during the computation of the ODE system.
If itrace < 1itrace<-1, then 1-1 is assumed and similarly if itrace > 3itrace>3, then 33 is assumed.
The advisory messages are given in greater detail as itrace increases. You are advised to set itrace = 0itrace=0, unless you are experienced with sub-chapter D02M–N.
14:   ind – int64int32nag_int scalar
Indicates whether this is a continuation call or a new integration.
ind = 0ind=0
Starts or restarts the integration in time.
ind = 1ind=1
Continues the integration after an earlier exit from the function. In this case, only the parameters tout and ifail should be reset between calls to nag_pde_1d_parab_convdiff (d03pf).
Constraint: ind = 0ind=0 or 11.

Optional Input Parameters

1:     npde – int64int32nag_int scalar
Default: The first dimension of the array u.
The number of PDEs to be solved.
Constraint: npde1npde1.
2:     npts – int64int32nag_int scalar
Default: The dimension of the array x and the second dimension of the array u. (An error is raised if these dimensions are not equal.)
The number of mesh points in the interval [a,b][a,b].
Constraint: npts3npts3.

Input Parameters Omitted from the MATLAB Interface

lrsave lisave

Output Parameters

1:     ts – double scalar
The value of tt corresponding to the solution values in u. Normally ts = toutts=tout.
2:     u(npde,npts) – double array
u(i,j)uij will contain the computed solution Ui(x,t)Ui(x,t) at x = x(j)x=xj and t = tst=ts, for i = 1,2,,npdei=1,2,,npde and j = 1,2,,nptsj=1,2,,npts.
3:     rsave(lrsave) – double array
If ind = 1ind=1, rsave must be unchanged from the previous call to the function because it contains required information about the iteration.
4:     isave(lisave) – int64int32nag_int array
If ind = 1ind=1, isave must be unchanged from the previous call to the function because it contains required information about the iteration. In particular:
isave(1)isave1
Contains the number of steps taken in time.
isave(2)isave2
Contains the number of residual evaluations of the resulting ODE system used. One such evaluation involves computing the PDE functions at all the mesh points, as well as one evaluation of the functions in the boundary conditions.
isave(3)isave3
Contains the number of Jacobian evaluations performed by the time integrator.
isave(4)isave4
Contains the order of the last backward differentiation formula method used.
isave(5)isave5
Contains the number of Newton iterations performed by the time integrator. Each iteration involves an ODE residual evaluation followed by a back-substitution using the LULU decomposition of the Jacobian matrix.
5:     ind – int64int32nag_int scalar
ind = 1ind=1.
6:     ifail – int64int32nag_int scalar
ifail = 0ifail=0 unless the function detects an error (see [Error Indicators and Warnings]).

Error Indicators and Warnings

Errors or warnings detected by the function:

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

  ifail = 1ifail=1
On entry,tstouttstout,
ortouttstout-ts is too small,
oritask = 1itask=1, 22 or 33,
ornpts < 3npts<3,
ornpde < 1npde<1,
orind0ind0 or 11,
orincorrect user-defined mesh, i.e., x(i)x(i + 1)xixi+1 for some i = 1,2,,npts1i=1,2,,npts-1,
orlrsave or lisave are too small,
orind = 1ind=1 on initial entry to nag_pde_1d_parab_convdiff (d03pf),
oracc(1)acc1 or acc(2) < 0.0acc2<0.0,
oracc(1)acc1 or acc(2)acc2 are both zero,
ortsmax < 0.0tsmax<0.0.
W ifail = 2ifail=2
The underlying ODE solver cannot make any further progress, with the values of acc, across the integration range from the current point t = tst=ts. The components of u contain the computed values at the current point t = tst=ts.
W ifail = 3ifail=3
In the underlying ODE solver, there were repeated error test failures on an attempted step, before completing the requested task, but the integration was successful as far as t = tst=ts. The problem may have a singularity, or the error requirement may be inappropriate. Incorrect specification of boundary conditions may also result in this error.
  ifail = 4ifail=4
In setting up the ODE system, the internal initialization function was unable to initialize the derivative of the ODE system. This could be due to the fact that ires was repeatedly set to 33 in one of pdedef, numflx or bndary when the residual in the underlying ODE solver was being evaluated. Incorrect specification of boundary conditions may also result in this error.
  ifail = 5ifail=5
In solving the ODE system, a singular Jacobian has been encountered. Check the problem formulation.
W ifail = 6ifail=6
When evaluating the residual in solving the ODE system, ires was set to 22 in at least one of pdedef, numflx or bndary. Integration was successful as far as t = tst=ts.
  ifail = 7ifail=7
The values of acc(1)acc1 and acc(2)acc2 are so small that the function is unable to start the integration in time.
  ifail = 8ifail=8
In either, pdedef, numflx or bndary, ires was set to an invalid value.
  ifail = 9ifail=9 (nag_ode_ivp_stiff_imp_revcom (d02nn))
A serious error has occurred in an internal call to the specified function. Check the problem specification and all parameters and array dimensions. Setting itrace = 1itrace=1 may provide more information. If the problem persists, contact NAG.
W ifail = 10ifail=10
The required task has been completed, but it is estimated that a small change in the values of acc is unlikely to produce any change in the computed solution. (Only applies when you are not operating in one step mode, that is when itask2itask2.)
  ifail = 11ifail=11
An error occurred during Jacobian formulation of the ODE system (a more detailed error description may be directed to the current advisory message unit when itrace1itrace1).
  ifail = 12ifail=12
Not applicable.
  ifail = 13ifail=13
Not applicable.
  ifail = 14ifail=14
One or more of the functions Pi,jPi,j, DiDi or CiCi was detected as depending on time derivatives, which is not permissible.

Accuracy

nag_pde_1d_parab_convdiff (d03pf) 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 components of the accuracy parameter, acc.

Further Comments

nag_pde_1d_parab_convdiff (d03pf) 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 tsmax. It is worth experimenting with this parameter, 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.

Example

For this function two examples are presented. There is a single example program for nag_pde_1d_parab_convdiff (d03pf), with a main program and the code to solve the two example problems given in Example 1 (EX1) and Example 2 (EX2).
Example 1 (EX1)
This example is a simple first-order system which illustrates the calculation of the numerical flux using Roe's approximate Riemann solver, and the specification of numerical boundary conditions using extrapolated characteristic variables. The PDEs are
(U1)/(t) + (U1)/(x) + (U2)/(x) = 0,
(U2)/(t) + 4 (U1)/(x) + (U2)/(x) = 0,
U1 t + U1 x + U2 x = 0, U2 t +4 U1 x + U2 x = 0,
for x [0,1] x [0,1]  and t 0 t 0 . The PDEs have an exact solution given by
U1 (x,t) = (1/2) {exp(x + t) + exp(x3t)} + (1/4) {sin(2π(x3t)2)sin(2π(x + t)2)} + 2 t2 2 x t ,
U2 (x,t) = exp(x3t) exp(x + t) + (1/2) {sin(2π(x3t)2) + sin(2π(x3t)2)} + x2 + 5 t2 2 x t .
U1 (x,t) = 12 { exp( x + t ) + exp(x-3t) } + 14 { sin( 2 π (x-3t) 2 ) - sin( 2 π (x+t) 2 ) } + 2 t2 - 2 x t , U2 (x,t) = exp(x-3t) - exp(x+t) + 12 { sin( 2 π (x-3t) 2 ) + sin( 2 π ( x - 3 t ) 2 ) } + x2 + 5 t2 - 2 x t .
The initial conditions are given by the exact solution. The characteristic variables are 2U1 + U22U1+U2 and 2U1U22U1-U2 corresponding to the characteristics given by dx / dt = 3dx/dt=3 and dx / dt = 1dx/dt=-1 respectively. Hence a physical boundary condition is required for 2U1 + U22U1+U2 at the left-hand boundary, and for 2U1U22U1-U2 at the right-hand boundary (corresponding to the incoming characteristics); and a numerical boundary condition is required for 2U1U22U1-U2 at the left-hand boundary, and for 2U1 + U22U1+U2 at the right-hand boundary (outgoing characteristics). The physical boundary conditions are obtained from the exact solution, and the numerical boundary conditions are calculated by linear extrapolation of the appropriate characteristic variable. The numerical flux is calculated using Roe's approximate Riemann solver: Using the notation in Section [Description], the flux vector FF and the Jacobian matrix AA are
F =
[ U1 + U2 4U1 + U2 ]
  and   A =
[ 1 1 4 1 ]
,
F=[ U1+U2 4U1+U2 ]  and   A=[ 1 1 4 1 ],
and the eigenvalues of AA are 33 and 1-1 with right eigenvectors 12
[]
T
[12]T and − 12
[]
T
[-12]T respectively. Using equation (4) the αkαk are given by
[ U1R − U1L U2R − U2L ]
= α1
[ 1 2 ]
+ α2
[ − 1 2 ]
,
[ U1R-U1L U2R-U2L ]=α1 [ 1 2 ]+α2 [ -1 2 ],
that is
α1 = (1/4) (2U1R2U1L + U2RU2L)   and   α2 = (1/4) (2U1R + 2U1L + U2RU2L) .
α1 = 14 ( 2 U1R - 2 U1L + U2R - U2L )   and   α2 = 14 ( -2 U1R + 2 U1L + U2R - U2L ) .
FLFL is given by
FL =
[ U1L + U2L 4U1L + U2L ]
,
FL = [ U1L+U2L 4U1L+U2L ] ,
and similarly for FRFR. From equation (4), the numerical flux vector is
F̂ = (1/2)
[ U1L + U2L + 0U1R + U2R 4U1L + U2L + 4U1R + U2R ]
− (1/2) α1 |3|
[ 1 2 ]
− (1/2) α2 | − 1|
[ − 1 2 ]
,
F^ = 12 [ U1L+U2L+0U1R+U2R 4U1L+U2L+4U1R+U2R ] - 12 α1 |3| [ 1 2 ] - 12 α2 |-1| [ -1 2 ] ,
that is
F̂ = (1/2)
[ 3U1L − 0U1R + (3/2)U2L + (1/2) U2R 6U1L + 2U1R + 3U2L − 0U2R ]
.
F^ = 12 [ 3U1L-0U1R+32U2L+12 U2R 6U1L+ 2U1R+ 3U2L-0U2R ] .
Example 2 (EX2)
This example is an advection-diffusion equation in which the flux term depends explicitly on xx:
(U)/(t) + x (U)/(x) = ε (2U)/(x2),
U t +x U x =ε 2U x2 ,
for x[1,1]x[-1,1] and 0t100t10. The parameter εε is taken to be 0.010.01. The two physical boundary conditions are U(1,t) = 3.0U(-1,t)=3.0 and U(1,t) = 5.0U(1,t)=5.0 and the initial condition is U(x,0) = x + 4U(x,0)=x+4. The integration is run to steady state at which the solution is known to be U = 4U=4 across the domain with a narrow boundary layer at both boundaries. In order to write the PDE in conservative form, a source term must be introduced, i.e.,
(U)/(t) + ((xU))/(x) = ε (2U)/(x2) + U.
U t + (xU) x =ε 2U x2 +U.
As in Example 1, the numerical flux is calculated using the Roe approximate Riemann solver. The Riemann problem to solve locally is
(U)/(t) + ((xU))/(x) = 0.
U t + (x U) x =0.
The xx in the flux term is assumed to be constant at a local level, and so using the notation in Section [Description], F = xUF=xU and A = xA=x. The eigenvalue is xx and the eigenvector (a scalar in this case) is 11. The numerical flux is therefore
F̂ =
{ xUL if ​x ≥ 0, xUR if ​x < 0.
F^={ xUL if ​x0, xUR if ​x<0.
function nag_pde_1d_parab_convdiff_example
ts = 0;
tout = 0.1;
u = [1, 1.010050167084168, 1.020201340026756, 1.030454533953517, ...
    1.040810774192388, 1.051271096376024, 1.06183654654536, ...
    1.072508181254217, 1.083287067674959, 1.09417428370521, ...
    1.105170918075648, 1.116278070458871, 1.127496851579376, ...
    1.138828383324622, 1.150273798857227, 1.161834242728283, ...
    1.17351087099181, 1.185304851320365, 1.19721736312181, ...
    1.209249597657251, 1.22140275816017, 1.233678059956743, ...
    1.246076730587381, 1.258600009929478, 1.271249150321405, ...
    1.284025416687741, 1.296930086665772, 1.309964450733247, ...
    1.323129812337437, 1.336427488025472, 1.349858807576003, ...
    1.363425114132178, 1.377127764335957, 1.39096812846378, ...
    1.404947590563594, 1.419067548593257, 1.43332941456034, ...
    1.447734614663325, 1.462284589434224, 1.476980793882643, ...
    1.49182469764127, 1.506817785112853, 1.521961555618634, ...
    1.537257523548281, 1.552707218511336, 1.568312185490169, ...
    1.584073984994482, 1.59999419321736, 1.616074402192893, ...
    1.632316219955379, 1.648721270700128, 1.665291194945886, ...
    1.682027649698886, 1.698932308618551, 1.716006862184859, ...
    1.733253017867395, 1.750672500296101, 1.768267051433735, ...
    1.786038430750073, 1.803988415397857, 1.822118800390509, ...
    1.840431398781637, 1.858928041846342, 1.877610579264343, ...
    1.896480879304951, 1.915540829013896, 1.934792334402032, ...
    1.95423732063594, 1.973877732230448, 1.993715533243082, ...
    2.013752707470477, 2.033991258646751, 2.054433210643888, ...
    2.075080607674122, 2.095935514494364, 2.117000016612675, ...
    2.138276220496818, 2.159766253784915, 2.181472265498201, ...
    2.203396426255937, 2.225540928492468, 2.247907986676472, ...
    2.270499837532406, 2.293318740264183, 2.316366976781091, ...
    2.339646851925991, 2.363160693705795, 2.386910853524276, ...
    2.41089970641721, 2.435129651289874, 2.45960311115695, ...
    2.484322533384816, 2.509290389936298, 2.534509177617855, ...
    2.559981418329271, 2.585709659315846, 2.611696473423118, ...
    2.637944459354153, 2.664456241929417, 2.691234472349262, ...
    2.718281828459045;
     0, 0.0007283184893761486, 0.002913271477003692, 0.006554836638408964, ...
    0.01165292715673064, 0.01820731731182075, 0.02621753831672977, ...
    0.03568274442576258, 0.04660154936295328, 0.05897183315779799, ...
    0.07279051952931347, 0.08805332403289762, 0.1047544732799356, ...
    0.1228863956604543, 0.1424393841471642, 0.1634012319375826, ...
    0.1857568419021818, 0.2094878110529674, 0.2345719915307735, ...
    0.2609830299326747, 0.2886898871648548, 0.3176563414121236, ...
    0.3478404772637143, 0.3791941645260716, 0.411662530786469, ...
    0.4451834323650897, 0.4796869289054211, 0.515094767500187, ...
    0.5513198829282125, 0.588265921280956, 0.6258267949789966, ...
    0.6638862779100518, 0.702317650151088, 0.7409834024559822, ...
    0.7797350113834701, 0.8184127965923141, 0.8568458724243957, ...
    0.8948522064124793, 0.9322387977664519, 0.9688019891867994, ...
    1.004327925502015, 1.038593172601059, 1.071365509904915, ...
    1.102404909163783, 1.131464711648689, 1.15829301479833, ...
    1.182634278054065, 1.204231155939365, 1.222826564387686, ...
    1.238165983870076, 1.25, 1.258087078981321, ...
    1.26219657150561, 1.26211193449322, 1.257634155413497, ...
    1.248585358827545, 1.234812569297339, 1.216191598935723, ...
    1.192631021684547, 1.164076189969761, 1.130513242775789, ...
    1.09197304750735, 1.04853501138347, 1.000330691672591, ...
    0.9475471279847715, 0.8904298142605735, 0.8292852232270232, ...
    0.764482792135553, 0.6964562757751382, 0.6257043712966427, ...
    0.5527905195293141, 0.4783417894574812, 0.403046756591692, ...
    0.3276522923413058, 0.2529591903857135, 0.1798165676349105, ...
    0.1091149918194066, 0.0417783051635255, -0.02124586596926006, ...
    -0.07899690212315057, -0.1305132427757891, -0.1748456304888814, ...
    -0.2110714466863253, -0.2383100047810946, -0.2557386309271434, ...
    -0.2626093261547741, -0.2582657667668944, -0.2421603634586167, ...
    -0.2138710646187889, -0.1731175567333843, -0.1197764858882511, ...
    -0.05389530029404987, 0.02429570418093795, 0.1143735639632994, ...
    0.2157154834078929, 0.327494747956721, 0.4486801170717872, ...
    0.5780390672923718, 0.7141452106973597, 0.855390155185901, ...
    1];
x = [0;
     0.01;
     0.02;
     0.03;
     0.04;
     0.05;
     0.06;
     0.07;
     0.08;
     0.09;
     0.1;
     0.11;
     0.12;
     0.13;
     0.14;
     0.15;
     0.16;
     0.17;
     0.18;
     0.19;
     0.2;
     0.21;
     0.22;
     0.23;
     0.24;
     0.25;
     0.26;
     0.27;
     0.28;
     0.29;
     0.3;
     0.31;
     0.32;
     0.33;
     0.34;
     0.35;
     0.36;
     0.37;
     0.38;
     0.39;
     0.4;
     0.41;
     0.42;
     0.43;
     0.44;
     0.45;
     0.46;
     0.47;
     0.48;
     0.49;
     0.5;
     0.51;
     0.52;
     0.53;
     0.54;
     0.55;
     0.56;
     0.57;
     0.58;
     0.59;
     0.6;
     0.61;
     0.62;
     0.63;
     0.64;
     0.65;
     0.66;
     0.67;
     0.68;
     0.69;
     0.7;
     0.71;
     0.72;
     0.73;
     0.74;
     0.75;
     0.76;
     0.77;
     0.78;
     0.79;
     0.8;
     0.81;
     0.82;
     0.83;
     0.84;
     0.85;
     0.86;
     0.87;
     0.88;
     0.89;
     0.9;
     0.91;
     0.92;
     0.93;
     0.94;
     0.95;
     0.96;
     0.97;
     0.98;
     0.99;
     1];
acc = [0.0001;
     1e-05];
tsmax = 0;
rsave = zeros(6695, 1);
isave = zeros(226, 1, 'int64');
itask = int64(1);
itrace = int64(0);
ind = int64(0);
[tsOut, uOut, rsaveOut, isaveOut, indOut, ifail] = ...
   nag_pde_1d_parab_convdiff(ts, tout, @pdedef, @numflx, @bndary, u, x, acc, ...
    tsmax, rsave, isave, itask, itrace, ind);
 tsOut, uOut, indOut, ifail

function [p, c, d, s, ires] = pdedef(npde, t, x, u, ux, ires)
    %  PDEDEF routine for simple hyperbolic system, equivalent to
    %  routine nag_pde_1d_parab_convdiff.
    p = zeros(npde, npde);
    c = zeros(npde, 1);
    d = zeros(npde, 1);
    s = zeros(npde, 1);
    for i = 1:double(npde)
      p(i,i) = 1;
    end
function [flux,ires] = numflx(npde, t, x, uleft, uright, ires)
      flux = zeros(npde,1);
      flux(1) = 0.5*(-uright(1)+3*uleft(1)+0.5*uright(2)+1.5*uleft(2));
      flux(2) = 0.5*(2*uright(1)+6*uleft(1)-uright(2)+3*uleft(2));

function [g, ires] = bndary(npde, npts, t, x, u, ibnd, ires)
      np1 = double(npts)-1; % npts is int64
      np2 = double(npts)-2; % npts is int64
      g = zeros(npde,1);
      if (ibnd == 0)
         ue = exact(t,npde,x(1),1);
         c = (x(2)-x(1))/(x(3)-x(2));
         exu1 = (1+c)*u(1,2) - c*u(1,3);
         exu2 = (1+c)*u(2,2) - c*u(2,3);
         g(1) = 2*u(1,1) + u(2,1) - 2*ue(1,1) - ue(2,1);
         g(2) = 2*u(1,1) - u(2,1) - 2*exu1 + exu2;
      else
         ue = exact(t,npde,x(npts),1);
         c = (x(npts)-x(np1))/(x(np1)-x(np2));
         exu1 = (1+c)*u(1,2) - c*u(1,3);
         exu2 = (1+c)*u(2,2) - c*u(2,3);
         g(1) = 2*u(1,1) - u(2,1) - 2*ue(1,1) + ue(2,1);
         g(2) = 2*u(1,1) + u(2,1) - 2*exu1 - exu2;
      end

% exact solution (for comparison and b.c. purposes)function u = exact(t,npde,x,npts)
  u = zeros(npde,npts);
  for i = 1:double(npts)
    x1 = x(i) + t;
    x2 = x(i) - 3*t;
    u(1,i)=0.5*(exp(x1)+exp(x2))+0.25*(sin(2*pi*x2^2)-sin(2*pi*x1^2))+2*t^2-2*x(i)*t;
    u(2,i)=exp(x2)-exp(x1)+0.5*(sin(2*pi*x2^2)+sin(2*pi*x1^2))+x(i)^2+5*t^2-2*x(i)*t;
  end
 

tsOut =

    0.1000


uOut =

  Columns 1 through 9

    1.0615    1.0573    1.0534    1.0492    1.0450    1.0408    1.0367    1.0325    1.0283
   -0.0155   -0.0296   -0.0435   -0.0562   -0.0679   -0.0785   -0.0879   -0.0961   -0.1032

  Columns 10 through 18

    1.0242    1.0202    1.0163    1.0123    1.0089    1.0055    1.0021    0.9990    0.9962
   -0.1091   -0.1137   -0.1171   -0.1195   -0.1202   -0.1199   -0.1186   -0.1161   -0.1125

  Columns 19 through 27

    0.9936    0.9913    0.9892    0.9876    0.9863    0.9855    0.9852    0.9853    0.9860
   -0.1077   -0.1020   -0.0953   -0.0877   -0.0793   -0.0698   -0.0597   -0.0489   -0.0377

  Columns 28 through 36

    0.9872    0.9891    0.9916    0.9949    0.9990    1.0039    1.0098    1.0167    1.0247
   -0.0260   -0.0140   -0.0018    0.0106    0.0231    0.0355    0.0478    0.0598    0.0713

  Columns 37 through 45

    1.0338    1.0441    1.0556    1.0684    1.0826    1.0982    1.1154    1.1341    1.1543
    0.0823    0.0926    0.1021    0.1106    0.1180    0.1242    0.1291    0.1325    0.1344

  Columns 46 through 54

    1.1763    1.1999    1.2253    1.2524    1.2809    1.3112    1.3433    1.3773    1.4129
    0.1346    0.1329    0.1295    0.1244    0.1175    0.1089    0.0982    0.0854    0.0706

  Columns 55 through 63

    1.4499    1.4885    1.5285    1.5698    1.6123    1.6558    1.7001    1.7452    1.7907
    0.0542    0.0361    0.0163   -0.0048   -0.0272   -0.0507   -0.0751   -0.1000   -0.1251

  Columns 64 through 72

    1.8365    1.8822    1.9278    1.9728    2.0171    2.0603    2.1023    2.1425    2.1809
   -0.1503   -0.1752   -0.1995   -0.2227   -0.2445   -0.2644   -0.2822   -0.2976   -0.3103

  Columns 73 through 81

    2.2171    2.2505    2.2813    2.3091    2.3335    2.3546    2.3723    2.3860    2.3959
   -0.3197   -0.3253   -0.3264   -0.3233   -0.3165   -0.3056   -0.2900   -0.2699   -0.2453

  Columns 82 through 90

    2.4022    2.4048    2.4035    2.3983    2.3895    2.3774    2.3623    2.3445    2.3243
   -0.2159   -0.1831   -0.1483   -0.1081   -0.0646   -0.0188    0.0281    0.0757    0.1233

  Columns 91 through 99

    2.3024    2.2793    2.2555    2.2313    2.2074    2.1844    2.1632    2.1443    2.1279
    0.1698    0.2143    0.2558    0.2934    0.3262    0.3534    0.3737    0.3862    0.3908

  Columns 100 through 101

    2.1135    2.1029
    0.3871    0.3760


indOut =

                    1


ifail =

                    0


function d03pf_example
ts = 0;
tout = 0.1;
u = [1, 1.010050167084168, 1.020201340026756, 1.030454533953517, ...
    1.040810774192388, 1.051271096376024, 1.06183654654536, ...
    1.072508181254217, 1.083287067674959, 1.09417428370521, ...
    1.105170918075648, 1.116278070458871, 1.127496851579376, ...
    1.138828383324622, 1.150273798857227, 1.161834242728283, ...
    1.17351087099181, 1.185304851320365, 1.19721736312181, ...
    1.209249597657251, 1.22140275816017, 1.233678059956743, ...
    1.246076730587381, 1.258600009929478, 1.271249150321405, ...
    1.284025416687741, 1.296930086665772, 1.309964450733247, ...
    1.323129812337437, 1.336427488025472, 1.349858807576003, ...
    1.363425114132178, 1.377127764335957, 1.39096812846378, ...
    1.404947590563594, 1.419067548593257, 1.43332941456034, ...
    1.447734614663325, 1.462284589434224, 1.476980793882643, ...
    1.49182469764127, 1.506817785112853, 1.521961555618634, ...
    1.537257523548281, 1.552707218511336, 1.568312185490169, ...
    1.584073984994482, 1.59999419321736, 1.616074402192893, ...
    1.632316219955379, 1.648721270700128, 1.665291194945886, ...
    1.682027649698886, 1.698932308618551, 1.716006862184859, ...
    1.733253017867395, 1.750672500296101, 1.768267051433735, ...
    1.786038430750073, 1.803988415397857, 1.822118800390509, ...
    1.840431398781637, 1.858928041846342, 1.877610579264343, ...
    1.896480879304951, 1.915540829013896, 1.934792334402032, ...
    1.95423732063594, 1.973877732230448, 1.993715533243082, ...
    2.013752707470477, 2.033991258646751, 2.054433210643888, ...
    2.075080607674122, 2.095935514494364, 2.117000016612675, ...
    2.138276220496818, 2.159766253784915, 2.181472265498201, ...
    2.203396426255937, 2.225540928492468, 2.247907986676472, ...
    2.270499837532406, 2.293318740264183, 2.316366976781091, ...
    2.339646851925991, 2.363160693705795, 2.386910853524276, ...
    2.41089970641721, 2.435129651289874, 2.45960311115695, ...
    2.484322533384816, 2.509290389936298, 2.534509177617855, ...
    2.559981418329271, 2.585709659315846, 2.611696473423118, ...
    2.637944459354153, 2.664456241929417, 2.691234472349262, ...
    2.718281828459045;
     0, 0.0007283184893761486, 0.002913271477003692, 0.006554836638408964, ...
    0.01165292715673064, 0.01820731731182075, 0.02621753831672977, ...
    0.03568274442576258, 0.04660154936295328, 0.05897183315779799, ...
    0.07279051952931347, 0.08805332403289762, 0.1047544732799356, ...
    0.1228863956604543, 0.1424393841471642, 0.1634012319375826, ...
    0.1857568419021818, 0.2094878110529674, 0.2345719915307735, ...
    0.2609830299326747, 0.2886898871648548, 0.3176563414121236, ...
    0.3478404772637143, 0.3791941645260716, 0.411662530786469, ...
    0.4451834323650897, 0.4796869289054211, 0.515094767500187, ...
    0.5513198829282125, 0.588265921280956, 0.6258267949789966, ...
    0.6638862779100518, 0.702317650151088, 0.7409834024559822, ...
    0.7797350113834701, 0.8184127965923141, 0.8568458724243957, ...
    0.8948522064124793, 0.9322387977664519, 0.9688019891867994, ...
    1.004327925502015, 1.038593172601059, 1.071365509904915, ...
    1.102404909163783, 1.131464711648689, 1.15829301479833, ...
    1.182634278054065, 1.204231155939365, 1.222826564387686, ...
    1.238165983870076, 1.25, 1.258087078981321, ...
    1.26219657150561, 1.26211193449322, 1.257634155413497, ...
    1.248585358827545, 1.234812569297339, 1.216191598935723, ...
    1.192631021684547, 1.164076189969761, 1.130513242775789, ...
    1.09197304750735, 1.04853501138347, 1.000330691672591, ...
    0.9475471279847715, 0.8904298142605735, 0.8292852232270232, ...
    0.764482792135553, 0.6964562757751382, 0.6257043712966427, ...
    0.5527905195293141, 0.4783417894574812, 0.403046756591692, ...
    0.3276522923413058, 0.2529591903857135, 0.1798165676349105, ...
    0.1091149918194066, 0.0417783051635255, -0.02124586596926006, ...
    -0.07899690212315057, -0.1305132427757891, -0.1748456304888814, ...
    -0.2110714466863253, -0.2383100047810946, -0.2557386309271434, ...
    -0.2626093261547741, -0.2582657667668944, -0.2421603634586167, ...
    -0.2138710646187889, -0.1731175567333843, -0.1197764858882511, ...
    -0.05389530029404987, 0.02429570418093795, 0.1143735639632994, ...
    0.2157154834078929, 0.327494747956721, 0.4486801170717872, ...
    0.5780390672923718, 0.7141452106973597, 0.855390155185901, ...
    1];
x = [0;
     0.01;
     0.02;
     0.03;
     0.04;
     0.05;
     0.06;
     0.07;
     0.08;
     0.09;
     0.1;
     0.11;
     0.12;
     0.13;
     0.14;
     0.15;
     0.16;
     0.17;
     0.18;
     0.19;
     0.2;
     0.21;
     0.22;
     0.23;
     0.24;
     0.25;
     0.26;
     0.27;
     0.28;
     0.29;
     0.3;
     0.31;
     0.32;
     0.33;
     0.34;
     0.35;
     0.36;
     0.37;
     0.38;
     0.39;
     0.4;
     0.41;
     0.42;
     0.43;
     0.44;
     0.45;
     0.46;
     0.47;
     0.48;
     0.49;
     0.5;
     0.51;
     0.52;
     0.53;
     0.54;
     0.55;
     0.56;
     0.57;
     0.58;
     0.59;
     0.6;
     0.61;
     0.62;
     0.63;
     0.64;
     0.65;
     0.66;
     0.67;
     0.68;
     0.69;
     0.7;
     0.71;
     0.72;
     0.73;
     0.74;
     0.75;
     0.76;
     0.77;
     0.78;
     0.79;
     0.8;
     0.81;
     0.82;
     0.83;
     0.84;
     0.85;
     0.86;
     0.87;
     0.88;
     0.89;
     0.9;
     0.91;
     0.92;
     0.93;
     0.94;
     0.95;
     0.96;
     0.97;
     0.98;
     0.99;
     1];
acc = [0.0001;
     1e-05];
tsmax = 0;
rsave = zeros(6695, 1);
isave = zeros(226, 1, 'int64');
itask = int64(1);
itrace = int64(0);
ind = int64(0);
[tsOut, uOut, rsaveOut, isaveOut, indOut, ifail] = ...
   d03pf(ts, tout, @pdedef, @numflx, @bndary, u, x, acc, ...
    tsmax, rsave, isave, itask, itrace, ind);
 tsOut, uOut, indOut, ifail

function [p, c, d, s, ires] = pdedef(npde, t, x, u, ux, ires)
    %  PDEDEF routine for simple hyperbolic system, equivalent to
    %  routine D03PFP
    p = zeros(npde, npde);
    c = zeros(npde, 1);
    d = zeros(npde, 1);
    s = zeros(npde, 1);
    for i = 1:double(npde)
      p(i,i) = 1;
    end
function [flux,ires] = numflx(npde, t, x, uleft, uright, ires)
      flux = zeros(npde,1);
      flux(1) = 0.5*(-uright(1)+3*uleft(1)+0.5*uright(2)+1.5*uleft(2));
      flux(2) = 0.5*(2*uright(1)+6*uleft(1)-uright(2)+3*uleft(2));

function [g, ires] = bndary(npde, npts, t, x, u, ibnd, ires)
      np1 = double(npts)-1; % npts is int64
      np2 = double(npts)-2; % npts is int64
      g = zeros(npde,1);
      if (ibnd == 0)
         ue = exact(t,npde,x(1),1);
         c = (x(2)-x(1))/(x(3)-x(2));
         exu1 = (1+c)*u(1,2) - c*u(1,3);
         exu2 = (1+c)*u(2,2) - c*u(2,3);
         g(1) = 2*u(1,1) + u(2,1) - 2*ue(1,1) - ue(2,1);
         g(2) = 2*u(1,1) - u(2,1) - 2*exu1 + exu2;
      else
         ue = exact(t,npde,x(npts),1);
         c = (x(npts)-x(np1))/(x(np1)-x(np2));
         exu1 = (1+c)*u(1,2) - c*u(1,3);
         exu2 = (1+c)*u(2,2) - c*u(2,3);
         g(1) = 2*u(1,1) - u(2,1) - 2*ue(1,1) + ue(2,1);
         g(2) = 2*u(1,1) + u(2,1) - 2*exu1 - exu2;
      end

% exact solution (for comparison and b.c. purposes)function u = exact(t,npde,x,npts)
  u = zeros(npde,npts);
  for i = 1:double(npts)
    x1 = x(i) + t;
    x2 = x(i) - 3*t;
    u(1,i)=0.5*(exp(x1)+exp(x2))+0.25*(sin(2*pi*x2^2)-sin(2*pi*x1^2))+2*t^2-2*x(i)*t;
    u(2,i)=exp(x2)-exp(x1)+0.5*(sin(2*pi*x2^2)+sin(2*pi*x1^2))+x(i)^2+5*t^2-2*x(i)*t;
  end
 

tsOut =

    0.1000


uOut =

  Columns 1 through 9

    1.0615    1.0573    1.0534    1.0492    1.0450    1.0408    1.0367    1.0325    1.0283
   -0.0155   -0.0296   -0.0435   -0.0562   -0.0679   -0.0785   -0.0879   -0.0961   -0.1032

  Columns 10 through 18

    1.0242    1.0202    1.0163    1.0123    1.0089    1.0055    1.0021    0.9990    0.9962
   -0.1091   -0.1137   -0.1171   -0.1195   -0.1202   -0.1199   -0.1186   -0.1161   -0.1125

  Columns 19 through 27

    0.9936    0.9913    0.9892    0.9876    0.9863    0.9855    0.9852    0.9853    0.9860
   -0.1077   -0.1020   -0.0953   -0.0877   -0.0793   -0.0698   -0.0597   -0.0489   -0.0377

  Columns 28 through 36

    0.9872    0.9891    0.9916    0.9949    0.9990    1.0039    1.0098    1.0167    1.0247
   -0.0260   -0.0140   -0.0018    0.0106    0.0231    0.0355    0.0478    0.0598    0.0713

  Columns 37 through 45

    1.0338    1.0441    1.0556    1.0684    1.0826    1.0982    1.1154    1.1341    1.1543
    0.0823    0.0926    0.1021    0.1106    0.1180    0.1242    0.1291    0.1325    0.1344

  Columns 46 through 54

    1.1763    1.1999    1.2253    1.2524    1.2809    1.3112    1.3433    1.3773    1.4129
    0.1346    0.1329    0.1295    0.1244    0.1175    0.1089    0.0982    0.0854    0.0706

  Columns 55 through 63

    1.4499    1.4885    1.5285    1.5698    1.6123    1.6558    1.7001    1.7452    1.7907
    0.0542    0.0361    0.0163   -0.0048   -0.0272   -0.0507   -0.0751   -0.1000   -0.1251

  Columns 64 through 72

    1.8365    1.8822    1.9278    1.9728    2.0171    2.0603    2.1023    2.1425    2.1809
   -0.1503   -0.1752   -0.1995   -0.2227   -0.2445   -0.2644   -0.2822   -0.2976   -0.3103

  Columns 73 through 81

    2.2171    2.2505    2.2813    2.3091    2.3335    2.3546    2.3723    2.3860    2.3959
   -0.3197   -0.3253   -0.3264   -0.3233   -0.3165   -0.3056   -0.2900   -0.2699   -0.2453

  Columns 82 through 90

    2.4022    2.4048    2.4035    2.3983    2.3895    2.3774    2.3623    2.3445    2.3243
   -0.2159   -0.1831   -0.1483   -0.1081   -0.0646   -0.0188    0.0281    0.0757    0.1233

  Columns 91 through 99

    2.3024    2.2793    2.2555    2.2313    2.2074    2.1844    2.1632    2.1443    2.1279
    0.1698    0.2143    0.2558    0.2934    0.3262    0.3534    0.3737    0.3862    0.3908

  Columns 100 through 101

    2.1135    2.1029
    0.3871    0.3760


indOut =

                    1


ifail =

                    0



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

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