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_coll_nth_comp (d02tg)

 Contents

    1  Purpose
    2  Syntax
    7  Accuracy
    9  Example

Purpose

nag_ode_bvp_coll_nth_comp (d02tg) solves a system of linear ordinary differential equations by least squares fitting of a series of Chebyshev polynomials using collocation.

Syntax

[c, ifail] = d02tg(m, l, x0, x1, k1, kp, coeff, bdyc, 'n', n)
[c, ifail] = nag_ode_bvp_coll_nth_comp(m, l, x0, x1, k1, kp, coeff, bdyc, 'n', n)

Description

nag_ode_bvp_coll_nth_comp (d02tg) calculates an approximate solution of a linear or linearized system of ordinary differential equations as a Chebyshev series. Suppose there are n differential equations for n variables y1,y2,,yn, over the range x0,x1. Let the ith equation be
j=1 mi+1 k=1 n fkjix y kj-1 x= rix  
where yk j x=djykx dxj . coeff evaluates the coefficients fkjix and the right-hand side rix for each i, 1in, at any point x. The boundary conditions may be applied either at the end points or at intermediate points; they are written in the same form as the differential equations, and specified by bdyc. For example the jth boundary condition out of those associated with the ith differential equation takes the form
j=1 li+1 k=1 n f kj ij xij ykj-1 xij=rij xij ,  
where xij lies between x0 and x1. It is assumed in this function that certain of the boundary conditions are associated with each differential equation. This is for your convenience; the grouping does not affect the results.
The degree of the polynomial solution must be the same for all variables. You specify the degree required, k1-1, and the number of collocation points, kp, in the range. The function sets up a system of linear equations for the Chebyshev coefficients, with n equations for each collocation point and one for each boundary condition. The collocation points are chosen at the extrema of a shifted Chebyshev polynomial of degree kp-1. The boundary conditions are satisfied exactly, and the remaining equations are solved by a least squares method. The result produced is a set of Chebyshev coefficients for the n functions y1,y2,,yn, with the range normalized to -1,1.
nag_fit_1dcheb_eval2 (e02ak) can be used to evaluate the components of the solution at any point on the range x0,x1 (see Example for an example). nag_fit_1dcheb_deriv (e02ah) and nag_fit_1dcheb_integ (e02aj) may be used to obtain Chebyshev series representations of derivatives and integrals (respectively) of the components of the solution.

References

Picken S M (1970) Algorithms for the solution of differential equations in Chebyshev-series by the selected points method Report Math. 94 National Physical Laboratory

Parameters

Compulsory Input Parameters

1:     mn int64int32nag_int array
mi must be set to the highest order derivative occurring in the ith equation, for i=1,2,,n.
Constraint: mi1, for i=1,2,,n.
2:     ln int64int32nag_int array
li must be set to the number of boundary conditions associated with the ith equation, for i=1,2,,n.
Constraint: li0, for i=1,2,,n.
3:     x0 – double scalar
The left-hand boundary, x0.
4:     x1 – double scalar
The right-hand boundary, x1.
Constraint: x1>x0.
5:     k1 int64int32nag_int scalar
The number of coefficients, k1, to be returned in the Chebyshev series representation of the solution (hence, the degree of the polynomial approximation is k1-1).
Constraint: k11+ max 1in mi.
6:     kp int64int32nag_int scalar
The number of collocation points to be used, kp.
Constraint: n×kpn×k1+ i=1 n li .
7:     coeff – function handle or string containing name of m-file
coeff defines the system of differential equations (see Description). It must evaluate the coefficient functions fkjix and the right-hand side function rix of the ith equation at a given point. Only nonzero entries of the array a and rhs need be specifically assigned, since all elements are set to zero by nag_ode_bvp_coll_nth_comp (d02tg) before calling coeff.
[a, rhs] = coeff(x, ii, a, ia, ia1, rhs)

Input Parameters

1:     x – double scalar
x, the point at which the functions must be evaluated.
2:     ii int64int32nag_int scalar
The equation for which the coefficients and right-hand side are to be evaluated.
3:     aiaia1 – double array
All elements of a are set to zero.
4:     ia int64int32nag_int scalar
5:     ia1 int64int32nag_int scalar
The first dimension of the array a and the second dimension of the array a.
6:     rhs – double scalar
Is set to zero.

Output Parameters

1:     aiaia1 – double array
akj must contain the value fkjix, for 1kn, 1jmi+1.
2:     rhs – double scalar
It must contain the value rix.
8:     bdyc – function handle or string containing name of m-file
bdyc defines the boundary conditions (see Description). It must evaluate the coefficient functions fkj ij and right-hand side function rij in the jth boundary condition associated with the ith equation, at the point xij at which the boundary condition is applied. Only nonzero entries of the array a and rhs need be specifically assigned, since all elements are set to zero by nag_ode_bvp_coll_nth_comp (d02tg) before calling bdyc.
[x, a, rhs] = bdyc(ii, j, a, ia, ia1, rhs)

Input Parameters

1:     ii int64int32nag_int scalar
The differential equation with which the condition is associated.
2:     j int64int32nag_int scalar
The boundary condition for which the coefficients and right-hand side are to be evaluated.
3:     aiaia1 – double array
All elements of a are set to zero.
4:     ia int64int32nag_int scalar
5:     ia1 int64int32nag_int scalar
The first dimension of the array a and the second dimension of the array a.
6:     rhs – double scalar
Is set to zero.

Output Parameters

1:     x – double scalar
xij, the value at which the boundary condition is applied.
2:     aiaia1 – double array
The value f kj ij xij , for 1kn, 1jmi+1.
3:     rhs – double scalar
The value rij xij .

Optional Input Parameters

1:     n int64int32nag_int scalar
Default: the dimension of the arrays m, l. (An error is raised if these dimensions are not equal.)
n, the number of differential equations in the system.
Constraint: n1.

Output Parameters

1:     cldcn – double array
The kth column of c contains the computed Chebyshev coefficients of the kth component of the solution, yk; that is, the computed solution is:
yk=i=1k1cikTi-1x,  1kn,  
where Tix is the Chebyshev polynomial of the first kind and denotes that the first coefficient, c1k, is halved.
2:     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:
   ifail=1
On entry,n<1,
ormi<1 for some i,
orli<0 for some i,
orx0x1,
ork1<1+mi for some i,
or n×kp<n×k1+ i=1 n li ,
orldc<k1.
   ifail=2
On entry,lw is too small (see Arguments),
orliw<n×k1.
   ifail=3
Either the boundary conditions are not linearly independent, or the rank of the matrix of equations for the coefficients is less than the number of unknowns. Increasing kp may overcome this latter problem.
   ifail=4
The least squares function nag_linsys_real_gen_lsqsol (f04am) has failed to correct the first approximate solution (see nag_linsys_real_gen_lsqsol (f04am)). Increasing kp may remove this difficulty.
   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

Estimates of the accuracy of the solution may be obtained by using the checks described in Further Comments. The Chebyshev coefficients are calculated by a stable numerical method.

Further Comments

The time taken by nag_ode_bvp_coll_nth_comp (d02tg) depends on the complexity of the system of differential equations, the degree of the polynomial solution and the number of matching points.
If the number of matching points kp is equal to the number of coefficients k1 minus the average number of boundary conditions 1n i=1 n li , then the least squares solution reduces to simple solution of linear equations and true collocation results. The accuracy of the solution may be checked by repeating the calculation with different values of k1. If the Chebyshev coefficients decrease rapidly, the size of the last two or three gives an indication of the error. If they do not decrease rapidly, it may be desirable to use a different method. Note that the Chebyshev coefficients are calculated for the range normalized to -1,1.
Generally the number of boundary conditions required is equal to the sum of the orders of the n differential equations. However, in some cases fewer boundary conditions are needed, because the assumption of a polynomial solution is equivalent to one or more boundary conditions (since it excludes singular solutions).
A system of nonlinear differential equations must be linearized before using the function. The calculation is repeated iteratively. On each iteration the linearized equation is used. In the example in Example, the y variables are to be determined at the current iteration whilst the z variables correspond to the solution determined at the previous iteration, (or the initial approximation on the first iteration). For a starting approximation, we may take, say, a linear function, and set up the appropriate Chebyshev coefficients before starting the iteration. For example, if y1=ax+b in the range x0,x1, we set B, the array of coefficients,
In some cases a better initial approximation may be needed and can be obtained by using nag_fit_1dcheb_arb (e02ad) or nag_fit_1dcheb_glp (e02af) to obtain a Chebyshev series for an approximate solution. The coefficients of the current iterate must be communicated to coeff and bdyc, e.g., via a global variable. (See Example.) The convergence of the (Newton) iteration cannot be guaranteed in general, though it is usually satisfactory from a good starting approximation.

Example

This example solves the nonlinear system
2y1+y22-1 y1+y2=0, 2y2-y1=0,  
in the range -1,1, with y1=0, y2=3, y2=0 at x=-1.
Suppose an approximate solution is z1, z2 such that y1z1, y2z2: then the first equation gives, on linearizing,
2y1+z22-1 y1+2z1z2+1 y2=2z1z22.  
The starting approximation is taken to be z1=0, z2=3. In the program below, the array B is used to hold the coefficients of the previous iterate (or of the starting approximation). We iterate until the Chebyshev coefficients converge to five figures. nag_fit_1dcheb_eval2 (e02ak) is used to calculate the solution from its Chebyshev coefficients.
function d02tg_example


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

global b x0 x1 ndegree

% Initialize variables and arrays.
n       = 2;
nderiv  = int64([1 2]);
nbc     = nderiv;
x0      = -1;
x1      = 1;
ncoeff  = int64(9);
ncolloc = int64(15);
ndegree = ncoeff-1;

% Iterate to covergence within tolerance of 1.0e0-5.
iter = 0;
emax = 1.0;
b = zeros(ncoeff, n);
b(1,2) = 6.0;
while emax >= 1.0e-05
  iter = iter + 1;

  [c, ifail] = d02tg(...
                     nderiv, nbc, x0, x1, ncoeff, ncolloc, ...
                     @coeff, @bdyc);
  % Output coefficients.
  fprintf('\n Iteration %2d. Chebyshev coefficients are:\n',iter);
  for j = 1:n
    fprintf('%8.4f',c(1:ncoeff,j));
    fprintf('\n');
  end
  
  emax = max(abs(c - b));
  b = c;
end

% Use these coefficients to evaluate the solution.
npt = 9;
fprintf('\n Solution evaluated at %d equally spaced points.\n\n', npt);
fprintf('%7s%16s\n','x','y');

% Prepare to store results for plotting.
dx = (x1-x0)/(double(npt-1));
xarray = [x0:dx:x1];
yarray = zeros(npt, n);

for ipt = 1:npt
  xx = xarray(ipt);
  fprintf('%10.4f', xx);
  for jeqn = 1:n
    [yarray(ipt, jeqn), ifail] = ...
    e02ak(...
           ndegree, x0, x1, c(:,jeqn), int64(1), xx);
    fprintf('%10.4f', yarray(ipt, jeqn));
  end
  fprintf('\n');
end

% Plot results.
fig1 = figure;
display_plot(xarray, yarray);


function [x, a, rhs] = bdyc(ii, jj, a, ia, ia1, rhs)
  % Evaluate coefficient functions and rhs of boundary conditions.

  x = -1;
  a(ii,jj) = 1;
  if (ii == 2 && jj == 1)
    rhs = 3;
  end

function [a, rhs] = coeff(x, ii, a, ia, ia1, rhs)
  % Evaluate coefficient functions and rhs of differential equations.

  global b x0 x1 ndegree

  if (ii <= 1)
    [z1, ifail] = e02ak(ndegree, x0, x1, b(:,1), int64(1), x);
    [z2, ifail] = e02ak(ndegree, x0, x1, b(:,2), int64(1), x);
    a(1,1) = z2*z2 - 1;
    a(1,2) = 2;
    a(2,1) = 2*z1*z2 + 1;
    rhs = 2*z1*z2*z2;
  else
    a(1,2) = -1;
    a(2,3) =  2;
  end

function display_plot(x, y)
  % Plot both curves.
  [haxes, hline1, hline2] = plotyy(x, y(:,1), x, y(:,2));
  % Set the axis limits and the tick specifications to beautify the plot.
  set(haxes(1), 'YLim', [-0.5 0.1]);
  set(haxes(1), 'YMinorTick', 'on');
  set(haxes(1), 'YTick', [-0.5:0.1:0]);
  set(haxes(2), 'YLim', [2.6 3.2]);
  set(haxes(2), 'YMinorTick', 'on');
  set(haxes(2), 'YTick', [2.6:0.1:3.2]);
  for iaxis = 1:2
    % These properties must be the same for both sets of axes.
    set(haxes(iaxis), 'XLim', [-1 1]);
  end
  % Add title.
  title('Solution by Chebyshev Collocation');
  % Label the axes.
  xlabel('x');
  ylabel(haxes(1), 'y_1');
  ylabel(haxes(2), 'y_2');
  % Add a legend.
  legend('y_1','y_2','Location','Best')
  % Set some features of the three lines.
  set(hline1, 'Linewidth', 0.25, 'Marker', '+', 'LineStyle', '-');
  set(hline2, 'Linewidth', 0.25, 'Marker', 'x', 'LineStyle', '--');
d02tg example results


 Iteration  1. Chebyshev coefficients are:
 -0.5659 -0.1162  0.0906 -0.0468  0.0196 -0.0069  0.0021 -0.0006  0.0001
  5.7083 -0.1642 -0.0087  0.0059 -0.0025  0.0009 -0.0003  0.0001 -0.0000

 Iteration  2. Chebyshev coefficients are:
 -0.6338 -0.1599  0.0831 -0.0445  0.0193 -0.0071  0.0023 -0.0006  0.0001
  5.6881 -0.1792 -0.0144  0.0053 -0.0023  0.0008 -0.0003  0.0001 -0.0000

 Iteration  3. Chebyshev coefficients are:
 -0.6344 -0.1604  0.0828 -0.0446  0.0193 -0.0071  0.0023 -0.0006  0.0001
  5.6880 -0.1793 -0.0145  0.0053 -0.0023  0.0008 -0.0003  0.0001 -0.0000

 Iteration  4. Chebyshev coefficients are:
 -0.6344 -0.1604  0.0828 -0.0446  0.0193 -0.0071  0.0023 -0.0006  0.0001
  5.6880 -0.1793 -0.0145  0.0053 -0.0023  0.0008 -0.0003  0.0001 -0.0000

 Solution evaluated at 9 equally spaced points.

      x               y
   -1.0000    0.0000    3.0000
   -0.7500   -0.2372    2.9827
   -0.5000   -0.3266    2.9466
   -0.2500   -0.3640    2.9032
    0.0000   -0.3828    2.8564
    0.2500   -0.3951    2.8077
    0.5000   -0.4055    2.7577
    0.7500   -0.4154    2.7064
    1.0000   -0.4255    2.6538
d02tg_fig1.png

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