NAG Library Routine Document
d03pdf
(dim1_parab_coll_old)
d03pda (dim1_parab_coll)
1
Purpose
d03pdf/d03pda integrates a system of linear or nonlinear parabolic partial differential equations (PDEs) in one space variable. 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 ordinary differential equations (ODEs). The resulting system is solved using a backward differentiation formula method.
d03pda is a version of
d03pdf that has additional arguments in order to make it safe for use in multithreaded applications (see
Section 5).
2
Specification
2.1
Specification for d03pdf
Fortran Interface
Subroutine d03pdf ( |
npde, m, ts, tout, pdedef, bndary, u, nbkpts, xbkpts, npoly, npts, x, uinit, acc, rsave, lrsave, isave, lisave, itask, itrace, ind, ifail) |
Integer, Intent (In) | :: | npde, m, nbkpts, npoly, npts, lrsave, lisave, itask, itrace | Integer, Intent (Inout) | :: | isave(lisave), ind, ifail | Real (Kind=nag_wp), Intent (In) | :: | tout, xbkpts(nbkpts), acc | Real (Kind=nag_wp), Intent (Inout) | :: | ts, u(npde,npts), rsave(lrsave) | Real (Kind=nag_wp), Intent (Out) | :: | x(npts) | External | :: | pdedef, bndary, uinit |
|
C Header Interface
#include <nagmk26.h>
void |
d03pdf_ (const Integer *npde, const Integer *m, double *ts, const double *tout, void (NAG_CALL *pdedef)(const Integer *npde, const double *t, const double x[], const Integer *nptl, const double u[], const double ux[], double p[], double q[], double r[], Integer *ires), void (NAG_CALL *bndary)(const Integer *npde, const double *t, const double u[], const double ux[], const Integer *ibnd, double beta[], double gamma[], Integer *ires), double u[], const Integer *nbkpts, const double xbkpts[], const Integer *npoly, const Integer *npts, double x[], void (NAG_CALL *uinit)(const Integer *npde, const Integer *npts, const double x[], double u[]), const double *acc, double rsave[], const Integer *lrsave, Integer isave[], const Integer *lisave, const Integer *itask, const Integer *itrace, Integer *ind, Integer *ifail) |
|
2.2
Specification for d03pda
Fortran Interface
Subroutine d03pda ( |
npde, m, ts, tout, pdedef, bndary, u, nbkpts, xbkpts, npoly, npts, x, uinit, acc, rsave, lrsave, isave, lisave, itask, itrace, ind, iuser, ruser, cwsav, lwsav, iwsav, rwsav, ifail) |
Integer, Intent (In) | :: | npde, m, nbkpts, npoly, npts, lrsave, lisave, itask, itrace | Integer, Intent (Inout) | :: | isave(lisave), ind, iuser(*), iwsav(505), ifail | Real (Kind=nag_wp), Intent (In) | :: | tout, xbkpts(nbkpts), acc | Real (Kind=nag_wp), Intent (Inout) | :: | ts, u(npde,npts), rsave(lrsave), ruser(*), rwsav(1100) | Real (Kind=nag_wp), Intent (Out) | :: | x(npts) | Logical, Intent (Inout) | :: | lwsav(100) | Character (80), Intent (Inout) | :: | cwsav(10) | External | :: | pdedef, bndary, uinit |
|
C Header Interface
#include <nagmk26.h>
void |
d03pda_ (const Integer *npde, const Integer *m, double *ts, const double *tout, void (NAG_CALL *pdedef)(const Integer *npde, const double *t, const double x[], const Integer *nptl, const double u[], const double ux[], double p[], double q[], double r[], Integer *ires, Integer iuser[], double ruser[]), void (NAG_CALL *bndary)(const Integer *npde, const double *t, const double u[], const double ux[], const Integer *ibnd, double beta[], double gamma[], Integer *ires, Integer iuser[], double ruser[]), double u[], const Integer *nbkpts, const double xbkpts[], const Integer *npoly, const Integer *npts, double x[], void (NAG_CALL *uinit)(const Integer *npde, const Integer *npts, const double x[], double u[], Integer iuser[], double ruser[]), const double *acc, double rsave[], const Integer *lrsave, Integer isave[], const Integer *lisave, const Integer *itask, const Integer *itrace, Integer *ind, Integer iuser[], double ruser[], char cwsav[], logical lwsav[], Integer iwsav[], double rwsav[], Integer *ifail, const Charlen length_cwsav) |
|
3
Description
d03pdf/d03pda integrates the system of parabolic equations:
where
,
and
depend on
,
,
,
and the vector
is the set of solution values
and the vector
is its partial derivative with respect to
. Note that
,
and
must not depend on
.
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 system is defined by the functions
,
and
which must be specified in
pdedef.
The initial values of the functions
must be given at
, and must be specified in
uinit.
The functions
, for
, which may be thought of as fluxes, are also used in the definition of the boundary conditions for each equation. The boundary conditions must have the form
where
or
.
The boundary conditions must be specified in
bndary. Thus, the problem is subject to the following restrictions:
(i) |
, so that integration is in the forward direction; |
(ii) |
, and the flux must not depend on any time derivatives; |
(iii) |
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 break-points ; |
(iv) |
at least one of the functions must be nonzero so that there is a time derivative present in the problem; |
(v) |
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 by either specifying the solution at or by specifying a zero flux there, that is and . See also Section 9. |
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
d03pdf/d03pda 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. In the case of just one element, the break-points are the boundaries. 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 and () ODEs for each PDE component between each pair of break-points. The system is then integrated forwards in time using a backward differentiation formula method.
4
References
Berzins M (1990) Developments in the NAG Library software for parabolic equations Scientific Software Systems (eds J C Mason and M G Cox) 59–72 Chapman and Hall
Berzins M and Dew P M (1991) Algorithm 690: Chebyshev polynomial software for elliptic-parabolic systems of PDEs ACM Trans. Math. Software 17 178–206
Zaturska N B, Drazin P G and Banks W H H (1988) On the flow of a viscous fluid driven along a channel by a suction at porous walls Fluid Dynamics Research 4
5
Arguments
- 1: – IntegerInput
-
On entry: the number of PDEs in the system to be solved.
Constraint:
.
- 2: – IntegerInput
-
On entry: the coordinate system used:
- Indicates Cartesian coordinates.
- Indicates cylindrical polar coordinates.
- Indicates spherical polar coordinates.
Constraint:
, or .
- 3: – Real (Kind=nag_wp)Input/Output
-
On entry: the initial value of the independent variable .
On exit: the value of
corresponding to the solution values in
u. Normally
.
Constraint:
.
- 4: – Real (Kind=nag_wp)Input
-
On entry: the final value of to which the integration is to be carried out.
- 5: – Subroutine, supplied by the user.External Procedure
-
pdedef must compute the values of the functions
,
and
which define the system of PDEs. The functions may depend on
,
,
and
and must be evaluated at a set of points.
The specification of
pdedef
for
d03pdf is:
Fortran Interface
Integer, Intent (In) | :: | npde, nptl | Integer, Intent (Inout) | :: | ires | Real (Kind=nag_wp), Intent (In) | :: | t, x(nptl), u(npde,nptl), ux(npde,nptl) | Real (Kind=nag_wp), Intent (Out) | :: | p(npde,npde,nptl), q(npde,nptl), r(npde,nptl) |
|
C Header Interface
#include <nagmk26.h>
void |
pdedef (const Integer *npde, const double *t, const double x[], const Integer *nptl, const double u[], const double ux[], double p[], double q[], double r[], Integer *ires) |
|
The specification of
pdedef
for
d03pda is:
Fortran Interface
Subroutine pdedef ( |
npde, t, x, nptl, u, ux, p, q, r, ires, iuser, ruser) |
Integer, Intent (In) | :: | npde, nptl | Integer, Intent (Inout) | :: | ires, iuser(*) | Real (Kind=nag_wp), Intent (In) | :: | t, x(nptl), u(npde,nptl), ux(npde,nptl) | Real (Kind=nag_wp), Intent (Inout) | :: | ruser(*) | Real (Kind=nag_wp), Intent (Out) | :: | p(npde,npde,nptl), q(npde,nptl), r(npde,nptl) |
|
C Header Interface
#include <nagmk26.h>
void |
pdedef (const Integer *npde, const double *t, const double x[], const Integer *nptl, const double u[], const double ux[], double p[], double q[], double r[], Integer *ires, Integer iuser[], double ruser[]) |
|
- 1: – IntegerInput
-
On entry: the number of PDEs in the system.
- 2: – Real (Kind=nag_wp)Input
-
On entry: the current value of the independent variable .
- 3: – Real (Kind=nag_wp) arrayInput
-
On entry: 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: – IntegerInput
-
On entry: the number of points at which evaluations are required (the value of ).
- 5: – Real (Kind=nag_wp) arrayInput
-
On entry: contains the value of the component where , for and .
- 6: – Real (Kind=nag_wp) arrayInput
-
On entry: contains the value of the component where , for and .
- 7: – Real (Kind=nag_wp) arrayOutput
-
On exit: must be set to the value of where , for , and .
- 8: – Real (Kind=nag_wp) arrayOutput
-
On exit: must be set to the value of where , for and .
- 9: – Real (Kind=nag_wp) arrayOutput
-
On exit: must be set to the value of where , for and .
- 10: – IntegerInput/Output
-
On entry: set to .
On exit: should usually remain unchanged. However, you may set
ires to force the integration routine 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 , d03pdf/d03pda returns to the calling subroutine with the error indicator set to .
- Note: the following are additional arguments for specific use with d03pda. Users of d03pdf therefore need not read the remainder of this description.
- 11: – Integer arrayUser Workspace
- 12: – Real (Kind=nag_wp) arrayUser Workspace
-
pdedef is called with the arguments
iuser and
ruser as supplied to
d03pdf/d03pda. You should use the arrays
iuser and
ruser to supply information to
pdedef.
pdedef must either be a module subprogram USEd by, or declared as EXTERNAL in, the (sub)program from which
d03pdf/d03pda is called. Arguments denoted as
Input must
not be changed by this procedure.
Note: pdedef should not return floating-point NaN (Not a Number) or infinity values, since these are not handled by
d03pdf/d03pda. If your code inadvertently
does return any NaNs or infinities,
d03pdf/d03pda is likely to produce unexpected results.
- 6: – Subroutine, supplied by the user.External Procedure
-
bndary must compute the functions
and
which define the boundary conditions as in equation
(3).
The specification of
bndary
for
d03pdf is:
Fortran Interface
Integer, Intent (In) | :: | npde, ibnd | Integer, Intent (Inout) | :: | ires | Real (Kind=nag_wp), Intent (In) | :: | t, u(npde), ux(npde) | Real (Kind=nag_wp), Intent (Out) | :: | beta(npde), gamma(npde) |
|
C Header Interface
#include <nagmk26.h>
void |
bndary (const Integer *npde, const double *t, const double u[], const double ux[], const Integer *ibnd, double beta[], double gamma[], Integer *ires) |
|
The specification of
bndary
for
d03pda is:
Fortran Interface
Integer, Intent (In) | :: | npde, ibnd | Integer, Intent (Inout) | :: | ires, iuser(*) | Real (Kind=nag_wp), Intent (In) | :: | t, u(npde), ux(npde) | Real (Kind=nag_wp), Intent (Inout) | :: | ruser(*) | Real (Kind=nag_wp), Intent (Out) | :: | beta(npde), gamma(npde) |
|
C Header Interface
#include <nagmk26.h>
void |
bndary (const Integer *npde, const double *t, const double u[], const double ux[], const Integer *ibnd, double beta[], double gamma[], Integer *ires, Integer iuser[], double ruser[]) |
|
- 1: – IntegerInput
-
On entry: the number of PDEs in the system.
- 2: – Real (Kind=nag_wp)Input
-
On entry: the current value of the independent variable .
- 3: – Real (Kind=nag_wp) arrayInput
-
On entry:
contains the value of the component
at the boundary specified by
ibnd, for
.
- 4: – Real (Kind=nag_wp) arrayInput
-
On entry:
contains the value of the component
at the boundary specified by
ibnd, for
.
- 5: – IntegerInput
-
On entry: 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, .
- 6: – Real (Kind=nag_wp) arrayOutput
-
On exit:
must be set to the value of
at the boundary specified by
ibnd, for
.
- 7: – Real (Kind=nag_wp) arrayOutput
-
On exit:
must be set to the value of
at the boundary specified by
ibnd, for
.
- 8: – IntegerInput/Output
-
On entry: set to .
On exit: should usually remain unchanged. However, you may set
ires to force the integration routine 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 , d03pdf/d03pda returns to the calling subroutine with the error indicator set to .
- Note: the following are additional arguments for specific use with d03pda. Users of d03pdf therefore need not read the remainder of this description.
- 9: – Integer arrayUser Workspace
- 10: – Real (Kind=nag_wp) arrayUser Workspace
-
bndary is called with the arguments
iuser and
ruser as supplied to
d03pdf/d03pda. You should use the arrays
iuser and
ruser to supply information to
bndary.
bndary must either be a module subprogram USEd by, or declared as EXTERNAL in, the (sub)program from which
d03pdf/d03pda is called. Arguments denoted as
Input must
not be changed by this procedure.
Note: bndary should not return floating-point NaN (Not a Number) or infinity values, since these are not handled by
d03pdf/d03pda. If your code inadvertently
does return any NaNs or infinities,
d03pdf/d03pda is likely to produce unexpected results.
- 7: – Real (Kind=nag_wp) arrayInput/Output
-
On entry: if
the value of
u must be unchanged from the previous call.
On exit: will contain the computed solution at .
- 8: – IntegerInput
-
On entry: the number of break-points in the interval .
Constraint:
.
- 9: – Real (Kind=nag_wp) arrayInput
-
On entry: the values of the break-points in the space direction. must specify the left-hand boundary, , and must specify the right-hand boundary, .
Constraint:
.
- 10: – IntegerInput
-
On entry: the degree of the Chebyshev polynomial to be used in approximating the PDE solution between each pair of break-points.
Constraint:
.
- 11: – IntegerInput
-
On entry: the number of mesh points in the interval .
Constraint:
.
- 12: – Real (Kind=nag_wp) arrayOutput
-
On exit: the mesh points chosen by
d03pdf/d03pda in the spatial direction. The values of
x will satisfy
.
- 13: – Subroutine, supplied by the user.External Procedure
-
uinit must compute the initial values of the PDE components
, for
and
.
The specification of
uinit
for
d03pdf is:
Fortran Interface
Integer, Intent (In) | :: | npde, npts | Real (Kind=nag_wp), Intent (In) | :: | x(npts) | Real (Kind=nag_wp), Intent (Out) | :: | u(npde,npts) |
|
C Header Interface
#include <nagmk26.h>
void |
uinit (const Integer *npde, const Integer *npts, const double x[], double u[]) |
|
The specification of
uinit
for
d03pda is:
Fortran Interface
Integer, Intent (In) | :: | npde, npts | Integer, Intent (Inout) | :: | iuser(*) | Real (Kind=nag_wp), Intent (In) | :: | x(npts) | Real (Kind=nag_wp), Intent (Inout) | :: | ruser(*) | Real (Kind=nag_wp), Intent (Out) | :: | u(npde,npts) |
|
C Header Interface
#include <nagmk26.h>
void |
uinit (const Integer *npde, const Integer *npts, const double x[], double u[], Integer iuser[], double ruser[]) |
|
- 1: – IntegerInput
-
On entry: the number of PDEs in the system.
- 2: – IntegerInput
-
On entry: the number of mesh points in the interval .
- 3: – Real (Kind=nag_wp) arrayInput
-
On entry: , contains the values of the th mesh point, for .
- 4: – Real (Kind=nag_wp) arrayOutput
-
On exit: must be set to the initial value , for and .
- Note: the following are additional arguments for specific use with d03pda. Users of d03pdf therefore need not read the remainder of this description.
- 5: – Integer arrayUser Workspace
- 6: – Real (Kind=nag_wp) arrayUser Workspace
-
uinit is called with the arguments
iuser and
ruser as supplied to
d03pdf/d03pda. You should use the arrays
iuser and
ruser to supply information to
uinit.
uinit must either be a module subprogram USEd by, or declared as EXTERNAL in, the (sub)program from which
d03pdf/d03pda is called. Arguments denoted as
Input must
not be changed by this procedure.
Note: uinit should not return floating-point NaN (Not a Number) or infinity values, since these are not handled by
d03pdf/d03pda. If your code inadvertently
does return any NaNs or infinities,
d03pdf/d03pda is likely to produce unexpected results.
- 14: – Real (Kind=nag_wp)Input
-
On entry: a positive quantity for controlling the local error estimate in the time integration. If
is the estimated error for
at the
th mesh point, the error test is:
Constraint:
.
- 15: – Real (Kind=nag_wp) arrayCommunication Array
-
If
,
rsave need not be set on entry.
If
,
rsave must be unchanged from the previous call to the routine because it contains required information about the iteration.
- 16: – IntegerInput
-
On entry: the dimension of the array
rsave as declared in the (sub)program from which
d03pdf/d03pda is called.
Constraint:
.
- 17: – Integer arrayCommunication Array
-
If
,
isave need not be set on entry.
If
,
isave must be unchanged from the previous call to the routine because it contains required information about the iteration. 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 last backward differentiation formula method used.
- 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 decomposition of the Jacobian matrix.
- 18: – IntegerInput
-
On entry: the dimension of the array
isave as declared in the (sub)program from which
d03pdf/d03pda is called.
Constraint:
.
- 19: – IntegerInput
-
On entry: 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 .
Constraint:
, or .
- 20: – IntegerInput
-
On entry: the level of trace information required from
d03pdf/d03pda 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 x04aaf).
- Output from the underlying ODE solver is printed on the current advisory message unit (see x04abf). This output contains details of Jacobian entries, the nonlinear iteration and the time integration during the computation of the ODE system.
If , is assumed and similarly if , 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.
- 21: – IntegerInput/Output
-
On entry: 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 routine. In this case, only the arguments tout and ifail should be reset between calls to d03pdf/d03pda.
Constraint:
or .
On exit: .
- 22: – IntegerInput/Output
-
Note: for d03pda, ifail does not occur in this position in the argument list. See the additional arguments described below.
On entry:
ifail must be set to
,
. If you are unfamiliar with this argument you should refer to
Section 3.4 in How to Use the NAG Library and its Documentation for details.
For environments where it might be inappropriate to halt program execution when an error is detected, the value
is recommended. If the output of error messages is undesirable, then the value
is recommended. Otherwise, if you are not familiar with this argument, the recommended value is
.
When the value is used it is essential to test the value of ifail on exit.
On exit:
unless the routine detects an error or a warning has been flagged (see
Section 6).
- Note: the following are additional arguments for specific use with d03pda. Users of d03pdf therefore need not read the remainder of this description.
- 22: – Integer arrayUser Workspace
-
iuser is not used by
d03pdf/d03pda, but is passed directly to
pdedef,
bndary and
uinit and may be used to pass information to these routines.
- 23: – Real (Kind=nag_wp) arrayUser Workspace
- 24: – Character(80) arrayCommunication Array
-
- 25: – Logical arrayCommunication Array
-
- 26: – Integer arrayCommunication Array
-
- 27: – Real (Kind=nag_wp) arrayCommunication Array
-
- 28: – IntegerInput/Output
-
Note: see the argument description for
ifail above.
6
Error Indicators and Warnings
If on entry
or
, explanatory error messages are output on the current error message unit (as defined by
x04aaf).
Errors or warnings detected by the routine:
-
On entry, .
Constraint: .
On entry, , , and .
Constraint: .
On entry, .
Constraint: or .
On entry, .
Constraint: , or .
On entry, .
Constraint: .
On entry, .
Constraint: .
On entry, .
Constraint: , or .
On entry, and .
Constraint: or
On entry, .
Constraint: .
On entry, .
Constraint: .
On entry, .
Constraint: .
On entry, .
Constraint: .
On entry, , and .
Constraint: .
On entry, on initial entry .
Constraint: on initial entry .
On entry, and .
Constraint: .
On entry, is too small:
and .
-
Underlying ODE solver cannot make further progress from the point
ts with the supplied value of
acc.
,
.
-
Repeated errors in an attempted step of underlying ODE solver. Integration was successful as far as
ts:
.
-
In setting up the ODE system an internal auxiliary was unable to initialize the derivative. This could be due to your setting
in
pdedef or
bndary.
-
Singular Jacobian of ODE system. Check problem formulation.
-
In evaluating residual of ODE system,
has been set in
pdedef or
bndary. Integration is successful as far as
ts:
.
-
acc was too small to start integration:
.
-
ires set to an invalid value in call to
pdedef or
bndary.
-
Serious error in internal call to an auxiliary. Increase
itrace for further details.
-
Integration completed, but a small change in
acc is unlikely to result in a changed solution.
.
-
Error during Jacobian formulation for ODE system. Increase
itrace for further details.
-
Flux function appears to depend on time derivatives.
An unexpected error has been triggered by this routine. Please
contact
NAG.
See
Section 3.9 in How to Use the NAG Library and its Documentation for further information.
Your licence key may have expired or may not have been installed correctly.
See
Section 3.8 in How to Use the NAG Library and its Documentation for further information.
Dynamic memory allocation failed.
See
Section 3.7 in How to Use the NAG Library and its Documentation for further information.
7
Accuracy
d03pdf/d03pda controls the accuracy of the integration in the time direction but not the accuracy of the approximation in space. The spatial accuracy depends on the degree of the polynomial approximation
npoly, and on both the number of break-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,
acc.
8
Parallelism and Performance
d03pdf/d03pda is threaded by NAG for parallel execution in multithreaded implementations of the NAG Library.
d03pdf/d03pda makes calls to BLAS and/or LAPACK routines, which may be threaded within the vendor library used by this implementation. Consult the documentation for the vendor library for further information.
Please consult the
X06 Chapter Introduction for information on how to control and interrogate the OpenMP environment used within this routine. Please also consult the
Users' Note for your implementation for any additional implementation-specific information.
d03pdf/d03pda is designed to solve parabolic systems (possibly including elliptic equations) with second-order derivatives in space. The argument specification allows you to include equations with only first-order derivatives in the space direction but there is no guarantee that the method of integration will be satisfactory for such systems. The position and nature of the boundary conditions in particular are critical in defining a stable problem.
The time taken depends on the complexity of the parabolic system and on the accuracy requested.
10
Example
The problem consists of a fourth-order PDE which can be written as a pair of second-order elliptic-parabolic PDEs for
and
,
where
and
. The boundary conditions are given by
The initial conditions at
are given by
The absence of boundary conditions for
does not pose any difficulties provided that the derivative flux boundary conditions are assigned to the first PDE
(4) which has the correct flux,
. The conditions on
at the boundaries are assigned to the second PDE by setting
in equation
(3) and placing the Dirichlet boundary conditions on
in the function
.
10.1
Program Text
Note: the following programs illustrate the use of d03pdf and d03pda.
Program Text (d03pdfe.f90)
Program Text (d03pdae.f90)
10.2
Program Data
Program Data (d03pdfe.d)
Program Data (d03pdae.d)
10.3
Program Results
Program Results (d03pdfe.r)
Program Results (d03pdae.r)