PDF version (NAG web site
, 64-bit version, 64-bit version)
NAG Toolbox: nag_pde_1d_parab_convdiff_dae (d03pl)
Purpose
nag_pde_1d_parab_convdiff_dae (d03pl) 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.
Syntax
[
ts,
u,
rsave,
isave,
ind,
ifail] = d03pl(
npde,
ts,
tout,
pdedef,
numflx,
bndary,
u,
x,
ncode,
odedef,
xi,
rtol,
atol,
itol,
norm_p,
laopt,
algopt,
rsave,
isave,
itask,
itrace,
ind, 'npts',
npts, 'nxi',
nxi, 'neqn',
neqn)
[
ts,
u,
rsave,
isave,
ind,
ifail] = nag_pde_1d_parab_convdiff_dae(
npde,
ts,
tout,
pdedef,
numflx,
bndary,
u,
x,
ncode,
odedef,
xi,
rtol,
atol,
itol,
norm_p,
laopt,
algopt,
rsave,
isave,
itask,
itrace,
ind, 'npts',
npts, 'nxi',
nxi, 'neqn',
neqn)
Note: the interface to this routine has changed since earlier releases of the toolbox:
At Mark 22: |
lrsave and lisave were removed from the interface |
Description
nag_pde_1d_parab_convdiff_dae (d03pl) integrates the system of convection-diffusion equations in conservative form:
or the hyperbolic convection-only system:
for
, where the vector
is the set of PDE solution values
The optional coupled ODEs are of the general form
where the vector
is the set of ODE solution values
denotes its derivative with respect to time, and
is the spatial derivative of
.
In
(1),
,
and
depend on
,
,
and
;
depends on
,
,
,
and
; and
depends on
,
,
,
and
linearly on
. Note that
,
,
and
must not depend on any space derivatives, and
,
,
and
must not depend on any time derivatives. In terms of conservation laws,
,
and
are the convective flux, diffusion and source terms respectively.
In
(3),
represents a vector of
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.
,
and
are the functions
,
and
evaluated at these coupling points. Each
may depend only linearly on time derivatives. Hence
(3) may be written more precisely as
where
,
is a vector of length
ncode,
is an
ncode by
ncode matrix,
is an
ncode by
matrix and the entries in
,
and
may depend on
,
,
,
and
. In practice you only need to supply a vector of information to define the ODEs and not the matrices
,
and
. (See
Arguments for the specification of
odedef.)
The integration in time is from to , over the space interval , where and are the leftmost and rightmost points of a user-defined mesh . The initial values of the functions and must be given at .
The PDEs are approximated by a system of ODEs in time for the values of
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
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,
say, must be calculated by you in terms of the
left and
right values of the solution vector
(denoted by
and
respectively), at each mid-point of the mesh
, for
. The left and right values are calculated by
nag_pde_1d_parab_convdiff_dae (d03pl) 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
is derived from the solution of the Riemann problem given by
where
, i.e.,
corresponds to
, with discontinuous initial values
for
and
for
, using an
approximate Riemann solver. This applies for either of the systems
(1) or
(2); the numerical flux is independent of the functions
,
,
and
. 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
or equivalently
. Provided the system is linear in
, i.e., the Jacobian matrix
does not depend on
, the numerical flux
is given by
where
(
) is the flux
calculated at the left (right) value of
, denoted by
(
); the
are the eigenvalues of
; the
are the right eigenvectors of
; and the
are defined by
An example is given in
Example and in the
nag_pde_1d_parab_convdiff (d03pf) documentation.
If the system is nonlinear, Roe's scheme requires that a linearized Jacobian is found (see
Roe (1981)).
The functions
,
,
and
(but
not
) must be specified in
pdedef. The numerical flux
must be supplied in a separate
numflx. For problems in the form
(2), the actual argument
nag_pde_1d_parab_convdiff_dae_sample_pdedef (d03plp) may be used for
pdedef.
nag_pde_1d_parab_convdiff_dae_sample_pdedef (d03plp) is included in the NAG Toolbox and sets the matrix with entries
to the identity matrix, and the functions
,
and
to zero.
The boundary condition specification has sufficient flexibility to allow for different types of problems. For second-order problems, i.e.,
depending on
, 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
Example and in the
nag_pde_1d_parab_convdiff (d03pf) documentation.
The boundary conditions must be specified in
bndary in the form
at the left-hand boundary, and
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
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
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:
(i) |
In (1), , for , may only appear linearly in the functions , for , with a similar restriction for and ; |
(ii) |
, , and must not depend on any space derivatives; and , , and must not depend on any time derivatives; |
(iii) |
, so that integration is in the forward direction; |
(iv) |
The evaluation of the terms , , and 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 ; |
(v) |
At least one of the functions must be nonzero so that there is a time derivative present in the PDE problem. |
In total there are
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.
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
Parameters
Compulsory Input Parameters
- 1:
– int64int32nag_int scalar
-
The number of PDEs to be solved.
Constraint:
.
- 2:
– double scalar
-
The initial value of the independent variable .
Constraint:
.
- 3:
– double scalar
-
The final value of to which the integration is to be carried out.
- 4:
– function handle or string containing name of m-file
-
pdedef must evaluate the functions
,
,
and
which partially define the system of PDEs.
and
may depend on
,
,
and
;
may depend on
,
,
,
and
; and
may depend on
,
,
,
and linearly on
.
pdedef is called approximately midway between each pair of mesh points in turn by
nag_pde_1d_parab_convdiff_dae (d03pl). The actual argument
nag_pde_1d_parab_convdiff_dae_sample_pdedef (d03plp) may be used for
pdedef for problems in the form
(2). (
nag_pde_1d_parab_convdiff_dae_sample_pdedef (d03plp) is included in the NAG Toolbox.)
[p, c, d, s, ires] = pdedef(npde, t, x, u, ux, ncode, v, vdot, ires)
Input Parameters
- 1:
– int64int32nag_int scalar
-
The number of PDEs in the system.
- 2:
– double scalar
-
The current value of the independent variable .
- 3:
– double scalar
-
The current value of the space variable .
- 4:
– double array
-
contains the value of the component , for .
- 5:
– double array
-
contains the value of the component , for .
- 6:
– int64int32nag_int scalar
-
The number of coupled ODEs in the system.
- 7:
– double array
-
If , contains the value of the component , for .
- 8:
– double array
-
If
,
contains the value of component
, for
.
Note:
, for , may only appear linearly in
, for .
- 9:
– int64int32nag_int scalar
-
Set to .
Output Parameters
- 1:
– double array
-
must be set to the value of , for and .
- 2:
– double array
-
must be set to the value of , for .
- 3:
– double array
-
must be set to the value of , for .
- 4:
– double array
-
must be set to the value of , for .
- 5:
– int64int32nag_int scalar
-
Should usually remain unchanged. However, you may set
ires to force the integration function to take certain actions as described below:
- Indicates to the integrator that control should be passed back immediately to the calling function with the error indicator set to .
- Indicates to the integrator that the current time step should be abandoned and a smaller time step used instead. You may wish to set when a physically meaningless input or output value has been generated. If you consecutively set , then nag_pde_1d_parab_convdiff_dae (d03pl) returns to the calling function with the error indicator set to .
- 5:
– 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
.
numflx is called approximately midway between each pair of mesh points in turn by
nag_pde_1d_parab_convdiff_dae (d03pl).
[flux, ires] = numflx(npde, t, x, ncode, v, uleft, uright, ires)
Input Parameters
- 1:
– int64int32nag_int scalar
-
The number of PDEs in the system.
- 2:
– double scalar
-
The current value of the independent variable .
- 3:
– double scalar
-
The current value of the space variable .
- 4:
– int64int32nag_int scalar
-
The number of coupled ODEs in the system.
- 5:
– double array
-
If , contains the value of the component , for .
- 6:
– double array
-
contains the left value of the component , for .
- 7:
– double array
-
contains the right value of the component , for .
- 8:
– int64int32nag_int scalar
-
Set to .
Output Parameters
- 1:
– double array
-
must be set to the numerical flux , for .
- 2:
– int64int32nag_int scalar
-
Should usually remain unchanged. However, you may set
ires to force the integration function to take certain actions as described below:
- Indicates to the integrator that control should be passed back immediately to the calling function with the error indicator set to .
- Indicates to the integrator that the current time step should be abandoned and a smaller time step used instead. You may wish to set when a physically meaningless input or output value has been generated. If you consecutively set , then nag_pde_1d_parab_convdiff_dae (d03pl) returns to the calling function with the error indicator set to .
- 6:
– function handle or string containing name of m-file
-
bndary must evaluate the functions
and
which describe the physical and numerical boundary conditions, as given by
(8) and
(9).
[g, ires] = bndary(npde, npts, t, x, u, ncode, v, vdot, ibnd, ires)
Input Parameters
- 1:
– int64int32nag_int scalar
-
The number of PDEs in the system.
- 2:
– int64int32nag_int scalar
-
The number of mesh points in the interval .
- 3:
– double scalar
-
The current value of the independent variable .
- 4:
– double array
-
The mesh points in the spatial direction. corresponds to the left-hand boundary, , and corresponds to the right-hand boundary, .
- 5:
– double array
-
contains the value of the component
at
, for
and
.
Note: if banded matrix algebra is to be used then the functions and may depend on the value of at the boundary point and the two adjacent points only.
- 6:
– int64int32nag_int scalar
-
The number of coupled ODEs in the system.
- 7:
– double array
-
If , contains the value of the component , for .
- 8:
– double array
-
If
,
contains the value of component
, for
.
Note:
, for , may only appear linearly in
and , for .
- 9:
– int64int32nag_int scalar
-
Specifies which boundary conditions are to be evaluated.
- bndary must evaluate the left-hand boundary condition at .
- bndary must evaluate the right-hand boundary condition at .
- 10:
– int64int32nag_int scalar
-
Set to .
Output Parameters
- 1:
– double array
-
must contain the
th component of either
or
in
(8) and
(9), depending on the value of
ibnd, for
.
- 2:
– int64int32nag_int scalar
-
Should usually remain unchanged. However, you may set
ires to force the integration function to take certain actions as described below:
- Indicates to the integrator that control should be passed back immediately to the calling function with the error indicator set to .
- Indicates to the integrator that the current time step should be abandoned and a smaller time step used instead. You may wish to set when a physically meaningless input or output value has been generated. If you consecutively set , then nag_pde_1d_parab_convdiff_dae (d03pl) returns to the calling function with the error indicator set to .
- 7:
– double array
-
The initial values of the dependent variables defined as follows:
-
contain , for and , and
-
contain , for .
- 8:
– double array
-
The mesh points in the space direction. must specify the left-hand boundary, , and must specify the right-hand boundary, .
Constraint:
.
- 9:
– int64int32nag_int scalar
-
The number of coupled ODE components.
Constraint:
.
- 10:
– function handle or string containing name of m-file
-
odedef must evaluate the functions
, which define the system of ODEs, as given in
(4).
If you wish to compute the solution of a system of PDEs only (i.e.,
),
odedef must be the string
nag_pde_1d_parab_dae_keller_remesh_fd_dummy_odedef (d03pek). (
nag_pde_1d_parab_dae_keller_remesh_fd_dummy_odedef (d03pek) is included in the NAG Toolbox.)
[r, ires] = odedef(npde, t, ncode, v, vdot, nxi, xi, ucp, ucpx, ucpt, ires)
Input Parameters
- 1:
– int64int32nag_int scalar
-
The number of PDEs in the system.
- 2:
– double scalar
-
The current value of the independent variable .
- 3:
– int64int32nag_int scalar
-
The number of coupled ODEs in the system.
- 4:
– double array
-
If , contains the value of the component , for .
- 5:
– double array
-
If , contains the value of component , for .
- 6:
– int64int32nag_int scalar
-
The number of ODE/PDE coupling points.
- 7:
– double array
-
If , contains the ODE/PDE coupling point, , for .
- 8:
– double array
-
The second dimension of the array
ucp must be at least
.
If , contains the value of at the coupling point , for and .
- 9:
– double array
-
The second dimension of the array
ucpx must be at least
.
If , contains the value of at the coupling point , for and .
- 10:
– double array
-
The second dimension of the array
ucpt must be at least
.
If , contains the value of at the coupling point , for and .
- 11:
– int64int32nag_int scalar
-
The form of
that must be returned in the array
r.
- Equation (10) must be used.
- Equation (11) must be used.
Output Parameters
- 1:
– double array
-
must contain the
th component of
, for
, where
is defined as
or
The definition of
is determined by the input value of
ires.
- 2:
– int64int32nag_int scalar
-
Should usually remain unchanged. However, you may reset
ires to force the integration function to take certain actions, as described below:
- Indicates to the integrator that control should be passed back immediately to the calling (sub)routine with the error indicator set to .
- Indicates to the integrator that the current time step should be abandoned and a smaller time step used instead. You may wish to set when a physically meaningless input or output value has been generated. If you consecutively set , then nag_pde_1d_parab_convdiff_dae (d03pl) returns to the calling function with the error indicator set to .
- 11:
– double array
-
The dimension of the array
xi
must be at least
, for , must be set to the ODE/PDE coupling points.
Constraint:
.
- 12:
– double array
-
The dimension of the array
rtol
must be at least
if
or
and at least
if
or
The relative local error tolerance.
Constraint:
for all relevant .
- 13:
– double array
-
The dimension of the array
atol
must be at least
if
or
and at least
if
or
The absolute local error tolerance.
Constraint:
for all relevant
.
Note: corresponding elements of
rtol and
atol cannot both be
.
- 14:
– int64int32nag_int scalar
-
A value to indicate the form of the local error test.
If
is the estimated local error for
, for
, and
denotes the norm, then the error test to be satisfied is
.
itol indicates to
nag_pde_1d_parab_convdiff_dae (d03pl) whether to interpret either or both of
rtol and
atol as a vector or scalar in the formation of the weights
used in the calculation of the norm (see the description of
norm_p):
itol | rtol | atol | |
1 | scalar | scalar | |
2 | scalar | vector | |
3 | vector | scalar | |
4 | vector | vector | |
Constraint:
.
- 15:
– string (length ≥ 1)
-
The type of norm to be used.
- Averaged norm.
- Averaged norm.
If
denotes the norm of the vector
u of length
neqn, then for the averaged
norm
and for the averaged
norm
See the description of
itol for the formulation of the weight vector
.
Constraint:
or .
- 16:
– string (length ≥ 1)
-
The type of matrix algebra required.
- Full matrix methods to be used.
- Banded matrix methods to be used.
- Sparse matrix methods to be used.
Constraint:
, or .
Note: you are recommended to use the banded option when no coupled ODEs are present (). 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.
- 17:
– double array
-
May be set to control various options available in the integrator. If you wish to employ all the default options, then
should be set to
. 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:
- Selects the ODE integration method to be used. If , a BDF method is used and if , a Theta method is used. The default is .
If , then
, for , are not used.
- Specifies the maximum order of the BDF integration formula to be used. may be , , , or . The default value is .
- Specifies what method is to be used to solve the system of nonlinear equations arising on each step of the BDF method. If a modified Newton iteration is used and if a functional iteration method is used. If functional iteration is selected and the integrator encounters difficulty, then there is an automatic switch to the modified Newton iteration. The default value is .
- 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
, for , for some or when there is no dependence in the coupled ODE system. If , then the Petzold test is used. If , then the Petzold test is not used. The default value is .
If , then
, for , are not used.
- Specifies the value of Theta to be used in the Theta integration method. . The default value is .
- Specifies what method is to be used to solve the system of nonlinear equations arising on each step of the Theta method. If , a modified Newton iteration is used and if , a functional iteration method is used. The default value is .
- 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 , then switching is allowed and if , then switching is not allowed. The default value is .
- Specifies a point in the time direction, , beyond which integration must not be attempted. The use of is described under the argument itask. If , a value of for , say, should be specified even if itask subsequently specifies that will not be used.
- Specifies the minimum absolute step size to be allowed in the time integration. If this option is not required, should be set to .
- Specifies the maximum absolute step size to be allowed in the time integration. If this option is not required, should be set to .
- Specifies the initial step size to be attempted by the integrator. If , then the initial step size is calculated internally.
- Specifies the maximum number of steps to be attempted by the integrator in any one call. If , then no limit is imposed.
- Specifies what method is to be used to solve the nonlinear equations at the initial point to initialize the values of , , and . If , a modified Newton iteration is used and if , functional iteration is used. The default value is .
and are used only for the sparse matrix algebra option, i.e., .
- Governs the choice of pivots during the decomposition of the first Jacobian matrix. It should lie in the range , with smaller values biasing the algorithm towards maintaining sparsity at the expense of numerical stability. If lies outside the range then the default value is used. If the functions regard the Jacobian matrix as numerically singular, then increasing towards may help, but at the cost of increased fill-in. The default value is .
- Used as the relative pivot threshold during subsequent Jacobian decompositions (see ) below which an internal error is invoked. must be greater than zero, otherwise the default value is used. If is greater than 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 ). The default value is .
- 18:
– double array
-
If
,
rsave need not be set on entry.
If
,
rsave must be unchanged from the previous call to the function because it contains required information about the iteration.
- 19:
– int64int32nag_int array
-
If
,
isave need not be set.
If
,
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:
- Contains the number of steps taken in time.
- 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.
- Contains the number of Jacobian evaluations performed by the time integrator.
- Contains the order of the BDF method last used in the time integration, if applicable. When the Theta method is used, contains no useful information.
- 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 decomposition of the Jacobian matrix.
- 20:
– int64int32nag_int scalar
-
The task to be performed by the ODE integrator.
- Normal computation of output values u at (by overshooting and interpolating).
- Take one step in the time direction and return.
- Stop at first internal integration point at or beyond .
- Normal computation of output values u at but without overshooting where is described under the argument algopt.
- Take one step in the time direction and return, without passing , where is described under the argument algopt.
Constraint:
, , , or .
- 21:
– int64int32nag_int scalar
-
The level of trace information required from
nag_pde_1d_parab_convdiff_dae (d03pl) and the underlying ODE solver.
itrace may take the value
,
,
,
or
.
- No output is generated.
- Only warning messages from the PDE solver are printed on the current error message unit (see nag_file_set_unit_error (x04aa)).
- 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 , then is assumed and similarly if , then is assumed.
The advisory messages are given in greater detail as
itrace increases. You are advised to set
, unless you are experienced with
Sub-chapter D02M–N.
- 22:
– int64int32nag_int scalar
-
Indicates whether this is a continuation call or a new integration.
- Starts or restarts the integration in time.
- Continues the integration after an earlier exit from the function. In this case, only the arguments tout and ifail should be reset between calls to nag_pde_1d_parab_convdiff_dae (d03pl).
Constraint:
or .
Optional Input Parameters
- 1:
– int64int32nag_int scalar
-
Default:
the dimension of the array
x.
The number of mesh points in the interval .
Constraint:
.
- 2:
– int64int32nag_int scalar
-
Default:
the dimension of the array
xi.
The number of ODE/PDE coupling points.
Constraints:
- if , ;
- if , .
- 3:
– int64int32nag_int scalar
-
Default:
the dimension of the array
u.
The number of ODEs in the time direction.
Constraint:
.
Output Parameters
- 1:
– double scalar
-
The value of
corresponding to the solution values in
u. Normally
.
- 2:
– double array
-
The computed solution , for and , and
, for , all evaluated at .
- 3:
– double array
-
If
,
rsave must be unchanged from the previous call to the function because it contains required information about the iteration.
- 4:
– int64int32nag_int array
-
If
,
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:
- Contains the number of steps taken in time.
- 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.
- Contains the number of Jacobian evaluations performed by the time integrator.
- Contains the order of the BDF method last used in the time integration, if applicable. When the Theta method is used, contains no useful information.
- 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 decomposition of the Jacobian matrix.
- 5:
– int64int32nag_int scalar
-
.
- 6:
– int64int32nag_int scalar
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.
-
-
On entry, | , |
or | is too small, |
or | , , , or , |
or | at least one of the coupling points defined in array xi is outside the interval [], |
or | the coupling points are not in strictly increasing order, |
or | , |
or | , |
or | , 'B' or 'S', |
or | , , or , |
or | or , |
or | mesh points are badly ordered, |
or | lrsave or lisave are too small, |
or | ncode and nxi are incorrectly defined, |
or | on initial entry to nag_pde_1d_parab_convdiff_dae (d03pl), |
or | , |
or | an element of rtol or , |
or | corresponding elements of rtol and atol are both , |
or | or . |
- W
-
The underlying ODE solver cannot make any further progress, with the values of
atol and
rtol, across the integration range from the current point
. The components of
u contain the computed values at the current point
.
- W
-
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 . The problem may have a singularity, or the error requirement may be inappropriate. Incorrect specification of boundary conditions may also result in this error.
-
-
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
in one of
pdedef,
numflx,
bndary or
odedef, when the residual in the underlying ODE solver was being evaluated. Incorrect specification of boundary conditions may also result in this error.
-
-
In solving the ODE system, a singular Jacobian has been encountered. Check the problem formulation.
- W
-
When evaluating the residual in solving the ODE system,
ires was set to
in at least one of
pdedef,
numflx,
bndary or
odedef. Integration was successful as far as
.
-
-
The values of
atol and
rtol are so small that the function is unable to start the integration in time.
-
-
In either,
pdedef,
numflx,
bndary or
odedef,
ires was set to an invalid value.
- (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 arguments and array dimensions. Setting
may provide more information. If the problem persists, contact
NAG.
- W
-
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
or
.)
-
-
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 ). If using the sparse matrix algebra option, the values of and may be inappropriate.
-
-
In solving the ODE system, the maximum number of steps specified in has been taken.
- W
-
Some error weights
became zero during the time integration (see the description of
itol). Pure relative error control (
) was requested on a variable (the
th) which has become zero. The integration was successful as far as
.
-
-
One or more of the functions , or was detected as depending on time derivatives, which is not permissible.
-
-
When using the sparse option, the value of lisave or lrsave was not sufficient (more detailed information may be directed to the current error message unit).
-
An unexpected error has been triggered by this routine. Please
contact
NAG.
-
Your licence key may have expired or may not have been installed correctly.
-
Dynamic memory allocation failed.
Accuracy
nag_pde_1d_parab_convdiff_dae (d03pl) 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.
Further Comments
nag_pde_1d_parab_convdiff_dae (d03pl) 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 . 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.
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
for
and
.
The PDEs have an exact solution given by
where
,
.
The initial conditions are given by the exact solution.
The characteristic variables are and , corresponding to the characteristics given by and respectively. Hence we require a physical boundary condition for at the left-hand boundary and for at the right-hand boundary (corresponding to the incoming characteristics), and a numerical boundary condition for at the left-hand boundary and for 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
at the left-hand boundary, and
at the right-hand boundary.
In order to specify these boundary conditions, two ODE variables
and
are introduced, defined by
The coupling points are therefore at
and
.
The numerical boundary conditions are now
at the left-hand boundary, and
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
Description for details), giving
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
where
is the density;
is the momentum, such that
, where
is the velocity;
is the specific energy; and
is the (constant) ratio of specific heats. The pressure
is given by
The solution domain is
for
, with the initial discontinuity at
, and initial conditions
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
Description 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
nag_pde_1d_parab_euler_roe (d03pu). An alternative Approxiate Riemann Solver using Osher's scheme is provided in
nag_pde_1d_parab_euler_osher (d03pv). Either
nag_pde_1d_parab_euler_roe (d03pu) or
nag_pde_1d_parab_euler_osher (d03pv) can be called from
numflx.
Open in the MATLAB editor:
d03pl_example
function d03pl_example
fprintf('d03pl example results\n\n');
fprintf('example 1\n\n');
ex1;
fprintf('\nexample 2\n\n');
ex2;
function ex1
npde = int64(2);
npts = int64(141);
ncode = int64(2);
neqn = npde*npts+ncode;
ts = 0;
tout = 0.5;
x = zeros(npts,1);
xi = [0; 1];
rtol = [0.00025];
atol = [1e-05];
itol = int64(1);
norm_p = '1';
laopt = 'S';
algopt = zeros(30,1);
algopt(1) = 1;
algopt(29) = 0.1;
algopt(30) = 1.1;
rsave = zeros(11000, 1);
isave = zeros(15700, 1, 'int64');
itask = int64(1);
itrace = int64(0);
ind = int64(0);
dx = 1/(double(npts)-1);
x = [0:dx:1];
u = exact(ts, npde, x, npts);
u = reshape(u, [2*npts,1]);
u(neqn-1) = u(1) - u(2);
u(neqn) = u(neqn-2) + u(neqn - 3);
[ts, u, rsave, isave, ind, ifail] = ...
d03pl( ...
npde, ts, tout, @pdedef, @numflx1, @bndary1, u, x, ncode, ...
@odedef, xi, rtol, atol, itol, norm_p, laopt, ...
algopt, rsave, isave, itask, itrace, ind);
xout = x(1:20:npts);
nop = length(xout);
uout(1, :) = u(1:20*npde:npts*npde);
uout(2, :) = u(2:20*npde:npts*npde);
ue = exact(tout, npde, xout, nop);
fprintf('\n%8s%12s%12s%12s%12s\n','x','U_1','u_1','U_2','u_2');
for i=1:nop
fprintf('%12.4f%12.4f%12.4f%12.4f%12.4f\n', ...
xout(i), uout(1,i), ue(1,i), uout(2,i), ue(2,i));
end
fig1 = figure;
plot(x,u(1:2:npts*npde),x,u(2:2:npts*npde));
title('Convection-diffusion solution at t=0.5');
xlabel('x');
ylabel('u_1, u_2');
legend('u_1','u_2','location','NorthWest');
function ex2
npde = int64(3);
npts = int64(141);
ncode = int64(0);
nxi = int64(0);
neqn = npde*npts+ncode;
ts = 0;
x = zeros(npts,1);
xi = [];
rtol = [0.0005];
atol = [0.005];
itol = int64(1);
norm_p = '2';
laopt = 'B';
algopt = zeros(30,1);
algopt(1) = 2;
algopt(6) = 2;
algopt(7) = 2;
algopt(13) = 0.005;
rsave = zeros(21000, 1);
isave = zeros(25700, 1, 'int64');
itask = int64(1);
itrace = int64(0);
ind = int64(0);
dx = 1/(double(npts)-1);
x = [0:dx:1];
u = uvinit(x);
for j=1:2
tout = 0.1*j;
[ts, u, rsave, isave, ind, ifail] = ...
d03pl( ...
npde, ts, tout, @pdedef, @numflx2, @bndary2, u, x, ncode, ...
@odedef, xi, rtol, atol, itol, norm_p, laopt, ...
algopt, rsave, isave, itask, itrace, ind,'nxi',nxi);
density = u(1,:);
velocity = u(2,:)./density;
pressure = 0.4*density.*(u(3,:)./density-velocity.^2/2);
if j==1
fig2 = figure;
else
fig3 = figure;
end
plot(x,density,x,velocity,x,pressure);
text = sprintf('Shock tube problem, t = %3.1f',tout);
title(text);
xlabel('x');
legend('density','velocity','pressure')
end
nsteps = 50*((isave(1)+25)/50);
nfuncs = 50*((isave(2)+25)/50);
njacs = isave(3);
niters = isave(5);
fprintf('\n Number of time steps (nearest 50) = %6d\n',nsteps);
fprintf(' Number of function evaluations (nearest 50) = %6d\n',nfuncs);
fprintf(' Number of Jacobian evaluations (nearest 1) = %6d\n',njacs);
fprintf(' Number of iterations (nearest 1) = %6d\n',niters);
function [g, ires] = bndary1(npde, npts, t, x, u, ncode, v, vdot, ibnd, ires)
g = zeros(npde, 1);
np1 = npts-1;
if (ibnd == 0)
ue = exact(t,npde,x(1),1);
g(1) = u(1,1) + u(2,1) - ue(1,1) - ue(2,1);
dudx = (u(1,2)-u(2,2)-u(1,1)+u(2,1))/(x(2)-x(1));
g(2) = vdot(1) - dudx;
else
ue = exact(t,npde,x(npts),1);
g(1) = u(1,npts) - u(2,npts) - ue(1,1) + ue(2,1);
dudx = (u(1,npts)+u(2,npts)-u(1,np1)-u(2,np1))/(x(npts)-x(np1));
g(2) = vdot(2) + 3*dudx;
end
function [u] = exact(t,npde,x,npts)
u = zeros(npde,npts);
pi2 = 2*pi;
for i = 1:double(npts)
f = exp(pi*(x(i)-3*t))*sin(pi2*(x(i)-3*t));
g = exp(-pi2*(x(i)+t))*cos(pi2*(x(i)+t));
u(1,i) = f + g;
u(2,i) = f - g;
end
function [flux, ires] = numflx1(npde, t, x, ncode, v, uleft, uright, ires)
flux = zeros(npde, 1);
flux(1) = (3*uleft(1)-uright(1)+3*uleft(2)+uright(2))/2;
flux(2) = (3*uleft(1)+uright(1)+3*uleft(2)-uright(2))/2;
function [r, ires] = odedef(npde, t, ncode, v, vdot, nxi, xi, ucp, ucpx, ...
ucpt, ires)
r = zeros(ncode, 1);
if (ires ~= -1)
r(1) = v(1) - ucp(1,1) + ucp(2,1);
r(2) = v(2) - ucp(1,2) - ucp(2,2);
end
function [p, c, d, s, ires] = pdedef(npde, t, x, u, ux, ncode, v, vdot, ires)
p = eye(npde);
c = ones(npde,1);
d = zeros(npde, 1);
s = d;
function [g, ires] = bndary2(npde, npts, t, x, u, ncode, v, vdot, ibnd, ires)
if (ibnd == 0)
g(1) = u(1,1) - 1;
g(2) = u(2,1);
g(3) = u(3,1) - 2.5;
else
g(1) = u(1,npts) - 1/8;
g(2) = u(2,npts);
g(3) = u(3,npts) - 1/4;
end
function [flux, ires] = numflx2(npde, t, x, ncode, v, uleft, uright, ires)
gamma = 1.4;
[flux, ifail] = d03pu( ...
uleft, uright, gamma);
function [u] = uvinit(x)
n = size(x,2);
u = zeros(3,n);
for i = 1:n
if x(i)<1/2
u(1,i) = 1;
u(3,i) = 2.5;
elseif x(i)== 1/2
u(1,i) = (1+1/8)/2;
u(3,i) = (2.5+1/4)/2;
else
u(1,i) = 1/8;
u(3,i) = 1/4;
end
end
d03pl example results
example 1
x U_1 u_1 U_2 u_2
0.0000 -0.0432 -0.0432 0.0432 0.0432
0.1429 -0.0221 -0.0220 0.0001 -0.0000
0.2857 -0.0199 -0.0199 -0.0231 -0.0231
0.4286 -0.0123 -0.0123 -0.0175 -0.0176
0.5714 0.0247 0.0245 0.0226 0.0224
0.7143 0.0832 0.0827 0.0830 0.0825
0.8571 0.1041 0.1036 0.1043 0.1039
1.0000 -0.0007 -0.0001 -0.0005 0.0001
example 2
Number of time steps (nearest 50) = 200
Number of function evaluations (nearest 50) = 450
Number of Jacobian evaluations (nearest 1) = 1
Number of iterations (nearest 1) = 2
PDF version (NAG web site
, 64-bit version, 64-bit version)
© The Numerical Algorithms Group Ltd, Oxford, UK. 2009–2015