PDF version (NAG web site
, 64-bit version, 64-bit version)
NAG Toolbox: nag_pde_1d_parab_dae_coll (d03pj)
Purpose
nag_pde_1d_parab_dae_coll (d03pj) integrates a system of linear or nonlinear parabolic partial differential equations (PDEs), in one space variable with scope for coupled ordinary differential equations (ODEs). The spatial discretization is performed using a Chebyshev collocation method, and the method of lines is employed to reduce the PDEs to a system of ODEs. The resulting system is solved using a backward differentiation formula (BDF) method or a Theta method (switching between Newton's method and functional iteration).
Syntax
[
ts,
u,
x,
rsave,
isave,
ind,
user,
cwsav,
lwsav,
iwsav,
rwsav,
ifail] = d03pj(
npde,
m,
ts,
tout,
pdedef,
bndary,
u,
xbkpts,
npoly,
npts,
ncode,
odedef,
xi,
uvinit,
rtol,
atol,
itol,
norm_p,
laopt,
algopt,
rsave,
isave,
itask,
itrace,
ind,
cwsav,
lwsav,
iwsav,
rwsav, 'nbkpts',
nbkpts, 'nxi',
nxi, 'neqn',
neqn, 'user',
user)
[
ts,
u,
x,
rsave,
isave,
ind,
user,
cwsav,
lwsav,
iwsav,
rwsav,
ifail] = nag_pde_1d_parab_dae_coll(
npde,
m,
ts,
tout,
pdedef,
bndary,
u,
xbkpts,
npoly,
npts,
ncode,
odedef,
xi,
uvinit,
rtol,
atol,
itol,
norm_p,
laopt,
algopt,
rsave,
isave,
itask,
itrace,
ind,
cwsav,
lwsav,
iwsav,
rwsav, 'nbkpts',
nbkpts, 'nxi',
nxi, 'neqn',
neqn, 'user',
user)
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_dae_coll (d03pj) integrates the system of parabolic-elliptic equations and coupled ODEs
where
(1) defines the PDE part and
(2) generalizes the coupled ODE part of the problem.
In
(1),
and
depend on
,
,
,
, and
;
depends on
,
,
,
,
and
linearly on
. The vector
is the set of PDE solution values
and the vector
is the partial derivative with respect to
. Note that
,
and
must not depend on
. The vector
is the set of ODE solution values
and
denotes its derivative with respect to time.
In
(2),
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 some of the PDE spatial mesh points.
,
,
,
and
are the functions
,
,
,
and
evaluated at these coupling points. Each
may only depend linearly on time derivatives. Hence the equation
(2) 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 need only 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 of a user-defined set of break-points . The coordinate system in space is defined by the value of ; for Cartesian coordinates, for cylindrical polar coordinates and for spherical polar coordinates.
The PDE system which is defined by the functions
,
and
must be specified in
pdedef.
The initial values of the functions
and
must be given at
. These values are calculated in
uvinit.
The functions
which may be thought of as fluxes, are also used in the definition of the boundary conditions. The boundary conditions must have the form
where
or
. The functions
may only depend
linearly on
.
The boundary conditions must be specified in
bndary.
The algebraic-differential equation system which is defined by the functions
must be specified in
odedef. You must also specify the coupling points
in the array
xi. Thus, the problem is subject to the following restrictions:
(i) |
in (1), , for , may only appear linearly in the functions
, for , with a similar restriction for ; |
(ii) |
and the flux must not depend on any time derivatives; |
(iii) |
, so that integration is in the forward direction; |
(iv) |
the evaluation of the functions , and is done at both the break-points and internally selected points for each element in turn, that is , and are evaluated twice at each break-point. 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; |
(vi) |
if and , which is the left boundary point, then it must be ensured that the PDE solution is bounded at this point. This can be done either by specifying the solution at or by specifying a zero flux there, that is and . |
The parabolic equations are approximated by a system of ODEs in time for the values of
at the mesh points. This ODE system is obtained by approximating the PDE solution between each pair of break-points by a Chebyshev polynomial of degree
npoly. The interval between each pair of break-points is treated by
nag_pde_1d_parab_dae_coll (d03pj) as an element, and on this element, a polynomial and its space and time derivatives are made to satisfy the system of PDEs at
spatial points, which are chosen internally by the code and the break-points. The user-defined break-points and the internally selected points together define the mesh. The smallest value that
npoly can take is one, in which case, the solution is approximated by piecewise linear polynomials between consecutive break-points and the method is similar to an ordinary finite element method.
In total there are
mesh points in the spatial direction, and
ODEs in the time direction; one ODE at each break-point for each PDE component,
ODEs for each PDE component between each pair of break-points, and
ncode coupled ODEs. The system is then integrated forwards in time using a Backward Differentiation Formula (BDF) method or a Theta method.
References
Berzins M (1990) Developments in the NAG Library software for parabolic equations Scientific Software Systems (eds J C Mason and M G Cox) 59–72 Chapman and Hall
Berzins M and Dew P M (1991) Algorithm 690: Chebyshev polynomial software for elliptic-parabolic systems of PDEs ACM Trans. Math. Software 17 178–206
Berzins M, Dew P M and Furzeland R M (1988) Software tools for time-dependent equations in simulation and optimization of large systems Proc. IMA Conf. Simulation and Optimization (ed A J Osiadcz) 35–50 Clarendon Press, Oxford
Berzins M and Furzeland R M (1992) An adaptive theta method for the solution of stiff and nonstiff differential equations Appl. Numer. Math. 9 1–19
Zaturska N B, Drazin P G and Banks W H H (1988) On the flow of a viscous fluid driven along a channel by a suction at porous walls Fluid Dynamics Research 4
Parameters
Compulsory Input Parameters
- 1:
– int64int32nag_int scalar
-
The number of PDEs to be solved.
Constraint:
.
- 2:
– int64int32nag_int scalar
-
The coordinate system used:
- Indicates Cartesian coordinates.
- Indicates cylindrical polar coordinates.
- Indicates spherical polar coordinates.
Constraint:
, or .
- 3:
– double scalar
-
The initial value of the independent variable .
Constraint:
.
- 4:
– double scalar
-
The final value of to which the integration is to be carried out.
- 5:
– function handle or string containing name of m-file
-
pdedef must compute the functions
,
and
which define the system of PDEs. The functions may depend on
,
,
,
and
;
may depend linearly on
. The functions must be evaluated at a set of points.
[p, q, r, ires, user] = pdedef(npde, t, x, nptl, u, ux, ncode, v, vdot, ires, user)
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 array
-
Contains a set of mesh points at which , and are to be evaluated. and contain successive user-supplied break-points and the elements of the array will satisfy .
- 4:
– int64int32nag_int scalar
-
The number of points at which evaluations are required (the value of ).
- 5:
– double array
-
contains the value of the component where , for and .
- 6:
– double array
-
contains the value of the component where , for and .
- 7:
– int64int32nag_int scalar
-
The number of coupled ODEs in the system.
- 8:
– double array
-
If , contains the value of the component , for .
- 9:
– double array
-
If
,
contains the value of component
, for
.
Note:
, for , may only appear linearly in
, for .
- 10:
– int64int32nag_int scalar
-
Set to .
- 11:
– Any MATLAB object
pdedef is called from
nag_pde_1d_parab_dae_coll (d03pj) with the object supplied to
nag_pde_1d_parab_dae_coll (d03pj).
Output Parameters
- 1:
– double array
-
must be set to the value of where , for , and .
- 2:
– double array
-
must be set to the value of where , for and .
- 3:
– double array
-
must be set to the value of where , for and .
- 4:
– 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 (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_dae_coll (d03pj) returns to the calling function with the error indicator set to .
- 5:
– Any MATLAB object
- 6:
– function handle or string containing name of m-file
-
bndary must compute the functions
and
which define the boundary conditions as in equation
(4).
[beta, gamma, ires, user] = bndary(npde, t, u, ux, ncode, v, vdot, ibnd, ires, user)
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 array
-
contains the value of the component
at the boundary specified by
ibnd, for
.
- 4:
– double array
-
contains the value of the component
at the boundary specified by
ibnd, for
.
- 5:
– int64int32nag_int scalar
-
The number of coupled ODEs in the system.
- 6:
– double array
-
If , contains the value of the component , for .
- 7:
– double array
-
If
,
contains the value of component
, for
.
Note:
, for , may only appear linearly in
, for .
- 8:
– int64int32nag_int scalar
-
Specifies which boundary conditions are to be evaluated.
- bndary must set up the coefficients of the left-hand boundary, .
- bndary must set up the coefficients of the right-hand boundary, .
- 9:
– int64int32nag_int scalar
-
Set to .
- 10:
– Any MATLAB object
bndary is called from
nag_pde_1d_parab_dae_coll (d03pj) with the object supplied to
nag_pde_1d_parab_dae_coll (d03pj).
Output Parameters
- 1:
– double array
-
must be set to the value of
at the boundary specified by
ibnd, for
.
- 2:
– double array
-
must be set to the value of
at the boundary specified by
ibnd, for
.
- 3:
– 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 (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_dae_coll (d03pj) returns to the calling function with the error indicator set to .
- 4:
– Any MATLAB object
- 7:
– double array
-
If
the value of
u must be unchanged from the previous call.
- 8:
– double array
-
The values of the break-points in the space direction. must specify the left-hand boundary, , and must specify the right-hand boundary, .
Constraint:
.
- 9:
– int64int32nag_int scalar
-
The degree of the Chebyshev polynomial to be used in approximating the PDE solution between each pair of break-points.
Constraint:
.
- 10:
– int64int32nag_int scalar
-
The number of mesh points in the interval .
Constraint:
.
- 11:
– int64int32nag_int scalar
-
The number of coupled ODE components.
Constraint:
.
- 12:
– function handle or string containing name of m-file
-
odedef must evaluate the functions
, which define the system of ODEs, as given in
(3).
If you wish to compute the solution of a system of PDEs only (
),
odedef must be the string
nag_pde_1d_parab_remesh_fd_dummy_odedef (d53pck).
[f, ires, user] = odedef(npde, t, ncode, v, vdot, nxi, xi, ucp, ucpx, rcp, ucpt, ucptx, ires, user)
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 points, , 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
rcp must be at least
.
contains the value of the flux at the coupling point , for and .
- 11:
– double array
-
The second dimension of the array
ucpt must be at least
.
If , contains the value of at the coupling point , for and .
- 12:
– double array
-
The second dimension of the array
ucptx must be at least
.
contains the value of at the coupling point , for and .
- 13:
– int64int32nag_int scalar
-
The form of
that must be returned in the array
f.
- Equation (5) must be used.
- Equation (6) must be used.
- 14:
– Any MATLAB object
odedef is called from
nag_pde_1d_parab_dae_coll (d03pj) with the object supplied to
nag_pde_1d_parab_dae_coll (d03pj).
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_dae_coll (d03pj) returns to the calling function with the error indicator set to .
- 3:
– Any MATLAB object
- 13:
– double array
-
The dimension of the array
xi
must be at least
, for , must be set to the ODE/PDE coupling points.
Constraint:
.
- 14:
– function handle or string containing name of m-file
-
uvinit must compute the initial values of the PDE and the ODE components
, for
and
, and
, for
.
[u, v, user] = uvinit(npde, npts, x, ncode, user)
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 array
-
, for , contains the current values of the space variable .
- 4:
– int64int32nag_int scalar
-
The number of coupled ODEs in the system.
- 5:
– Any MATLAB object
uvinit is called from
nag_pde_1d_parab_dae_coll (d03pj) with the object supplied to
nag_pde_1d_parab_dae_coll (d03pj).
Output Parameters
- 1:
– double array
-
If , contains the value of the component , for and .
- 2:
– double array
-
contains the value of component , for .
- 3:
– Any MATLAB object
- 15:
– 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 .
- 16:
– 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
.
- 17:
– int64int32nag_int scalar
-
A value to indicate the form of the local error test.
itol indicates to
nag_pde_1d_parab_dae_coll (d03pj) whether to interpret either or both of
rtol or
atol as a vector or scalar. The error test to be satisfied is
, where
is defined as follows:
itol | rtol | atol | |
1 | scalar | scalar | |
2 | scalar | vector | |
3 | vector | scalar | |
4 | vector | vector | |
In the above, denotes the estimated local error for the th component of the coupled PDE/ODE system in time, , for .
The choice of norm used is defined by the argument
norm_p.
Constraint:
.
- 18:
– string (length ≥ 1)
-
The type of norm to be used.
- Maximum norm.
- Averaged norm.
If
denotes the norm of the vector
u of length
neqn, then for the averaged
norm
while for the maximum norm
See the description of
itol for the formulation of the weight vector
.
Constraint:
or .
- 19:
– 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 (i.e., ).
- 20:
– 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 value 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, .
- 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 this 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 .
- Is used as a relative pivot threshold during subsequent Jacobian decompositions (see ) below which an internal error is invoked. If is greater than no check is made on the pivot size, and this may be a necessary option if the Jacobian is found to be numerically singular (see ). The default value is .
- 21:
– 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.
- 22:
– int64int32nag_int array
-
If
,
isave need not be set on entry.
If
,
isave must be unchanged from the previous call to the function because it contains required information about the iteration required for subsequent calls. In particular:
- Contains the number of steps taken in time.
- 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.
- Contains the number of Jacobian evaluations performed by the time integrator.
- Contains the order of the ODE method last used in the time integration.
- 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.
- 23:
– int64int32nag_int scalar
-
Specifies the task to be performed by the ODE integrator.
- Normal computation of output values u at .
- One step 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 .
- 24:
– int64int32nag_int scalar
-
The level of trace information required from
nag_pde_1d_parab_dae_coll (d03pj) 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.
- 25:
– 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_dae_coll (d03pj).
Constraint:
or .
- 26:
– cell array of strings
-
- 27:
– logical array
-
- 28:
– int64int32nag_int array
-
- 29:
– double array
-
Optional Input Parameters
- 1:
– int64int32nag_int scalar
-
Default:
the dimension of the array
xbkpts.
The number of break-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:
.
- 4:
– Any MATLAB object
user is not used by
nag_pde_1d_parab_dae_coll (d03pj), but is passed to
pdedef,
bndary,
odedef and
uvinit. Note that for large objects it may be more efficient to use a global variable which is accessible from the m-files than to use
user.
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
, evaluated at
, as follows:
-
contain , for and , and
-
contain , for .
- 3:
– double array
-
The mesh points chosen by
nag_pde_1d_parab_dae_coll (d03pj) in the spatial direction. The values of
x will satisfy
.
- 4:
– double array
-
If
,
rsave must be unchanged from the previous call to the function because it contains required information about the iteration.
- 5:
– int64int32nag_int array
-
If
,
isave must be unchanged from the previous call to the function because it contains required information about the iteration required for subsequent calls. In particular:
- Contains the number of steps taken in time.
- 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.
- Contains the number of Jacobian evaluations performed by the time integrator.
- Contains the order of the ODE method last used in the time integration.
- 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.
- 6:
– int64int32nag_int scalar
-
.
- 7:
– Any MATLAB object
- 8:
– cell array of strings
-
- 9:
– logical array
-
- 10:
– int64int32nag_int array
-
- 11:
– double array
-
- 12:
– 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, | is too small, |
or | , , , or , |
or | , or , |
or | at least one of the coupling point in array xi is outside the interval [], |
or | , |
or | , |
or | , |
or | or , |
or | , , or , |
or | or , |
or | ncode and nxi are incorrectly defined, |
or | , |
or | , or , |
or | or , |
or | break-points are badly ordered, |
or | lrsave is too small, |
or | lisave is too small, |
or | the ODE integrator has not been correctly defined; check algopt argument, |
or | either an element of rtol or , |
or | all the elements of rtol and atol are zero. |
- 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.
-
-
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 at least
pdedef,
bndary or
odedef, when the residual in the underlying ODE solver was being evaluated.
-
-
In solving the ODE system, a singular Jacobian has been encountered. You should check your problem formulation.
- W
-
When evaluating the residual in solving the ODE system,
ires was set to
in at least
pdedef,
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 one of
pdedef,
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 error message unit).
-
-
In solving the ODE system, the maximum number of steps specified in have 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
.
-
-
The flux function 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_dae_coll (d03pj) controls the accuracy of the integration in the time direction but not the accuracy of the approximation in space. The spatial accuracy depends on both the number of mesh points and on their distribution in space. In the time integration only the local error over a single step is controlled and so the accuracy over a number of steps cannot be guaranteed. You should therefore test the effect of varying the accuracy argument
atol and
rtol.
Further Comments
The argument specification allows you to include equations with only first-order derivatives in the space direction but there is no guarantee that the method of integration will be satisfactory for such systems. The position and nature of the boundary conditions in particular are critical in defining a stable problem.
The time taken depends on the complexity of the parabolic system and on the accuracy requested.
Example
This example provides a simple coupled system of one PDE and one ODE.
for
.
The left boundary condition at
is
The right boundary condition at
is
The initial conditions at
are defined by the exact solution:
and the coupling point is at
.
Open in the MATLAB editor:
d03pj_example
function d03pj_example
fprintf('d03pj example results\n\n');
npde = int64(1);
m = int64(0);
ts = 0.0001;
tout = 0.2;
u = zeros(22,1);
xbkpts = [0:0.1:1];
npoly = int64(2);
npts = int64(21);
ncode = int64(1);
xi = [1];
rtol = [0.0001];
atol = rtol;
itol = int64(1);
normt = 'A';
laopt = 'F';
algopt = zeros(30,1);
rsave = zeros(900, 1);
isave = zeros(24, 1, 'int64');
itask = int64(1);
itrace = int64(0);
ind = int64(0);
cwsav = {''; ''; ''; ''; ''; ''; ''; ''; ''; ''};
lwsav = false(100, 1);
iwsav = zeros(505, 1, 'int64');
rwsav = zeros(1100, 1);
t = [0.1:3.1/19:3.2];
for i_t = 1:numel(t)
tout = t(i_t);
[ts, u, x, rsave, isave, ind, user, cwsav, lwsav, iwsav, rwsav, ifail] = ...
d03pj( ...
npde, m, ts, tout, @pdedef, @bndary, u, xbkpts, ...
npoly, npts, ncode, @odedef, xi, @uvinit, rtol, atol, itol, ...
normt, laopt, algopt, rsave, isave, itask, itrace, ind, ...
cwsav, lwsav, iwsav, rwsav);
if (i_t==1)
xs = reshape(x,[7,3]);
end
if mod(i_t+4,6)==0
fprintf('\nThe solution at t = %7.4f is:\n',ts);
for j = 0:2
fprintf('%10s%12s','x','u(x,t)');
end
fprintf('\n');
us = reshape(u(1:21),[7,3]);
for i = 1:7
for j = 1:3
fprintf('%12.2f%10.4f',xs(i,j),us(i,j));
end
fprintf('\n');
end
end
v(:,i_t) = u(1:(end-1));
end
fig1 = figure;
x = reshape(xs,[21,1]);
mesh(t,x,v);
xlabel('t');
ylabel('x');
zlabel('u(x,t)');
title('Coupled Parabolic PDE/ODE using Collocation and BDF');
function [p,q,r,ires,user] = pdedef(npde,t,x,nptl,u,ux,ncode,v,vdot,ires,user)
p = zeros(npde,npde,nptl);
q = zeros(npde,nptl);
r = zeros(npde,nptl);
p(1,1,:) = v(1)*v(1);
r(1,:) = ux(1,1:nptl);
q(1,1:nptl) = -x(1:nptl)'.*ux(1,1:nptl)*(v(1)*vdot(1));
function [beta,gamma,ires,user] = bndary(npde,t,u,ux,ncode,v,vdot,ibnd,...
ires,user)
beta = zeros(npde,1);
gamma = zeros(npde,1);
beta(1) = 1;
if (ibnd == 0)
gamma(1) = -v(1)*exp(t);
else
gamma(1) = -v(1)*vdot(1);
end
function [f,ires,user] = odedef(npde,t,ncode,v,vdot,nxi,xi,ucp, ...
ucpx,rcp,ucpt,ucptx,ires,user)
f = zeros(ncode,1);
if (ires == 1)
f(1) = vdot(1) - v(1)*ucp(1,1) - ucpx(1,1) - 1 - t;
elseif (ires == -1)
f(1) = vdot(1);
end
function [u,v,user] = uvinit(npde,npts,x,ncode,user)
u = zeros(npde,npts);
ts = 1e-4;
v(1) = ts;
u = exp(ts*(1-x)) - 1;
d03pj example results
The solution at t = 0.2632 is:
x u(x,t) x u(x,t) x u(x,t)
0.00 0.3014 0.35 0.1869 0.70 0.0825
0.05 0.2844 0.40 0.1714 0.75 0.0683
0.10 0.2676 0.45 0.1561 0.80 0.0543
0.15 0.2510 0.50 0.1409 0.85 0.0406
0.20 0.2347 0.55 0.1260 0.90 0.0270
0.25 0.2185 0.60 0.1113 0.95 0.0135
0.30 0.2026 0.65 0.0968 1.00 0.0003
The solution at t = 1.2421 is:
x u(x,t) x u(x,t) x u(x,t)
0.00 2.4634 0.35 1.2421 0.70 0.4514
0.05 2.2548 0.40 1.1071 0.75 0.3639
0.10 2.0588 0.45 0.9801 0.80 0.2818
0.15 1.8746 0.50 0.8609 0.85 0.2045
0.20 1.7014 0.55 0.7488 0.90 0.1319
0.25 1.5387 0.60 0.6434 0.95 0.0637
0.30 1.3858 0.65 0.5444 1.00 -0.0004
The solution at t = 2.2211 is:
x u(x,t) x u(x,t) x u(x,t)
0.00 8.2169 0.35 3.2347 0.70 0.9449
0.05 7.2478 0.40 2.7893 0.75 0.7401
0.10 6.3805 0.45 2.3908 0.80 0.5569
0.15 5.6045 0.50 2.0341 0.85 0.3929
0.20 4.9099 0.55 1.7149 0.90 0.2461
0.25 4.2885 0.60 1.4293 0.95 0.1147
0.30 3.7323 0.65 1.1736 1.00 -0.0029
The solution at t = 3.2000 is:
x u(x,t) x u(x,t) x u(x,t)
0.00 23.5336 0.35 6.9983 0.70 1.6052
0.05 19.9047 0.40 5.8146 0.75 1.2189
0.10 16.8115 0.45 4.8060 0.80 0.8898
0.15 14.1765 0.50 3.9466 0.85 0.6092
0.20 11.9308 0.55 3.2142 0.90 0.3701
0.25 10.0177 0.60 2.5902 0.95 0.1662
0.30 8.3873 0.65 2.0584 1.00 -0.0076
PDF version (NAG web site
, 64-bit version, 64-bit version)
© The Numerical Algorithms Group Ltd, Oxford, UK. 2009–2015