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

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

NAG Toolbox: nag_ode_bvp_ps_lin_solve (d02ue)

 Contents

    1  Purpose
    2  Syntax
    7  Accuracy
    9  Example

Purpose

nag_ode_bvp_ps_lin_solve (d02ue) finds the solution of a linear constant coefficient boundary value problem by using the Chebyshev integration formulation on a Chebyshev Gauss–Lobatto grid.

Syntax

[bmat, f, uc, resid, ifail] = d02ue(n, a, b, c, bmat, y, bvec, f, 'm', m)
[bmat, f, uc, resid, ifail] = nag_ode_bvp_ps_lin_solve(n, a, b, c, bmat, y, bvec, f, 'm', m)

Description

nag_ode_bvp_ps_lin_solve (d02ue) solves the constant linear coefficient ordinary differential problem
j=0 m fj+1 dju dxj = fx ,  x a,b  
subject to a set of m linear constraints at points yia,b , for i=1,2,,m:
j=0 m B i,j+1 dju dxj x=yi = βi ,  
where 1m4, B is an m×m+1 matrix of constant coefficients and βi are constants. The points yi are usually either a or b.
The function fx is supplied as an array of Chebyshev coefficients cj, j=0,1,,n for the function discretized on n+1 Chebyshev Gauss–Lobatto points (as returned by nag_ode_bvp_ps_lin_cgl_grid (d02uc)); the coefficients are normally obtained by a previous call to nag_ode_bvp_ps_lin_coeffs (d02ua). The solution and its derivatives (up to order m) are returned, in the form of their Chebyshev series representation, as arrays of Chebyshev coefficients; subsequent calls to nag_ode_bvp_ps_lin_cgl_vals (d02ub) will return the corresponding function and derivative values at the Chebyshev Gauss–Lobatto discretization points on a,b. Function and derivative values can be obtained on any uniform grid over the same range a,b by calling the interpolation function nag_ode_bvp_ps_lin_grid_vals (d02uw).

References

Clenshaw C W (1957) The numerical solution of linear differential equations in Chebyshev series Proc. Camb. Phil. Soc. 53 134–149
Coutsias E A, Hagstrom T and Torres D (1996) An efficient spectral method for ordinary differential equations with rational function coefficients Mathematics of Computation 65(214) 611–635
Greengard L (1991) Spectral integration and two-point boundary value problems SIAM J. Numer. Anal. 28(4) 1071–80
Lundbladh A, Hennigson D S and Johannson A V (1992) An efficient spectral integration method for the solution of the Navier–Stokes equations Technical report FFA–TN 1992–28 Aeronautical Research Institute of Sweden
Muite B K (2010) A numerical comparison of Chebyshev methods for solving fourth-order semilinear initial boundary value problems Journal of Computational and Applied Mathematics 234(2) 317–342

Parameters

Compulsory Input Parameters

1:     n int64int32nag_int scalar
n, where the number of grid points is n+1.
Constraint: n8 and n is even.
2:     a – double scalar
a, the lower bound of domain a,b.
Constraint: a<b.
3:     b – double scalar
b, the upper bound of domain a,b.
Constraint: b>a.
4:     cn+1 – double array
The Chebyshev coefficients cj, j=0,1,,n, for the right hand side of the boundary value problem. Usually these are obtained by a previous call of nag_ode_bvp_ps_lin_coeffs (d02ua).
5:     bmatmm+1 – double array
bmatij+1 must contain the coefficients Bi,j+1, for i=1,2,,m and j=0,1,,m, in the problem formulation of Description.
6:     ym – double array
The points, yi, for i=1,2,,m, where the boundary conditions are discretized.
7:     bvecm – double array
The values, βi, for i=1,2,,m, in the formulation of the boundary conditions given in Description.
8:     fm+1 – double array
The coefficients, fj, for j=1,2,,m+1, in the formulation of the linear boundary value problem given in Description. The highest order term, fm+1, needs to be nonzero to have a well posed problem.

Optional Input Parameters

1:     m int64int32nag_int scalar
Default: the first dimension of the array bmat and the dimension of the arrays y, bvec. (An error is raised if these dimensions are not equal.)
The order, m, of the boundary value problem to be solved.
Constraint: 1m4.

Output Parameters

1:     bmatmm+1 – double array
The coefficients have been scaled to form an equivalent problem defined on the domain -1,1.
2:     fm+1 – double array
The coefficients have been scaled to form an equivalent problem defined on the domain -1,1.
3:     ucn+1m+1 – double array
The Chebyshev coefficients in the Chebyshev series representations of the solution and derivatives of the solution to the boundary value problem. The n+1 elements uc1:n+11 contain the coefficients representing the solution Uxi, for i=0,1,,n. uc1:n+1j+1 contains the coefficients representing the jth derivative of U, for j=1,2,,m.
4:     resid – double scalar
The maximum residual resulting from substituting the solution vectors returned in uc into both linear equations of Description representing the linear boundary value problem and associated boundary conditions. That is
max max i=1,m j=0 m B i,j+1 dju dxj x=yi - βi , max i=1, n+1 j=0 m f j+1 dju dxj x=xi - fx .  
5:     ifail int64int32nag_int scalar
ifail=0 unless the function detects an error (see Error Indicators and Warnings).

Error Indicators and Warnings

Errors or warnings detected by the function:

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

   ifail=1
Constraint: n is even.
Constraint: n8.
   ifail=2
Constraint: a<b.
   ifail=3
On entry, fm+1=0.0.
   ifail=6
Constraint: 1m4.
   ifail=7
Internal error while unpacking matrix during iterative refinement.
Please contact NAG.
   ifail=8
Singular matrix encountered during iterative refinement.
Please check that your system is well posed.
W  ifail=9
During iterative refinement, the maximum number of iterations was reached.
W  ifail=10
During iterative refinement, convergence was achieved, but the residual is more than 100×machine precision.
   ifail=-99
An unexpected error has been triggered by this routine. Please contact NAG.
   ifail=-399
Your licence key may have expired or may not have been installed correctly.
   ifail=-999
Dynamic memory allocation failed.

Accuracy

The accuracy should be close to machine precision for well conditioned boundary value problems.

Further Comments

The number of operations is of the order n logn  and the memory requirements are On; thus the computation remains efficient and practical for very fine discretizations (very large values of n). Collocation methods will be faster for small problems, but the method of nag_ode_bvp_ps_lin_solve (d02ue) should be faster for larger discretizations.

Example

This example solves the third-order problem 4Uxxx + 3Uxx + 2Ux + U = 2sinx - 2cosx  on -π/2,π/2  subject to the boundary conditions U -π/2 = 0 , 3Uxx -π/2 + 2Ux -π/2 + U -π/2 =2, and 3Uxx π/2 + 2Ux π/2 + U π/2 =-2  using the Chebyshev integration formulation on a Chebyshev Gauss–Lobatto grid of order 16.
function d02ue_example


fprintf('d02ue example results\n\n');

n = int64(16);
a = -pi/2;
b =  pi/2;


% Set up boundary condition on left side of domain
y = [a, b];
% Set up Dirichlet condition using exact solution at x=a.
bmat = zeros(2, 3);
bmat(1, 1:2) = [1, 1];
bmat(2, 1:2) = [1, 1];
bvec = [cos(a) - sin(a), cos(b) - sin(b)];

% Set up problem definition
f = [1, 2, 3];

% Set up solution grid
[x, ifail] = d02uc(n, a, b);

% Set up problem right hand sides for grid and transform
f0 = -2*sin(x) - 2*cos(x);
[c, ifail] = d02ua(n, f0);

% Solve in coefficient space
[bmat, f, uc, resid, ifail] = d02ue(n, a, b, c, bmat, y, bvec, f);

% Evaluate solution and derivatives on Chebyshev grid
u = zeros(n+1, 3);
for q=0:2
  [u(:, q+1),  ifail] = d02ub(n, a, b, int64(q), uc(:, q+1));
end

% Print Solution
fprintf('\nNumerical solution U and its first two derivatives\n');
fprintf('      x          U          Ux         Uxx\n');
fprintf('%10.4f %10.4f %10.4f %10.4f\n', [x u]');


d02ue example results


Numerical solution U and its first two derivatives
      x          U          Ux         Uxx
   -1.5708    -0.0000     1.0000     0.0000
   -1.5406     0.0302     0.9995    -0.0302
   -1.4512     0.1193     0.9929    -0.1193
   -1.3061     0.2616     0.9652    -0.2616
   -1.1107     0.4440     0.8960    -0.4440
   -0.8727     0.6428     0.7661    -0.6428
   -0.6011     0.8247     0.5656    -0.8247
   -0.3064     0.9534     0.3017    -0.9534
   -0.0000     1.0000     0.0000    -1.0000
    0.3064     0.9534    -0.3017    -0.9534
    0.6011     0.8247    -0.5656    -0.8247
    0.8727     0.6428    -0.7661    -0.6428
    1.1107     0.4440    -0.8960    -0.4440
    1.3061     0.2616    -0.9652    -0.2616
    1.4512     0.1193    -0.9929    -0.1193
    1.5406     0.0302    -0.9995    -0.0302
    1.5708    -0.0000    -1.0000    -0.0000

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

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