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_coeffs (d02ua)

 Contents

    1  Purpose
    2  Syntax
    7  Accuracy
    9  Example

Purpose

nag_ode_bvp_ps_lin_coeffs (d02ua) obtains the Chebyshev coefficients of a function discretized on Chebyshev Gauss–Lobatto points. The set of discretization points on which the function is evaluated is usually obtained by a previous call to nag_ode_bvp_ps_lin_cgl_grid (d02uc).

Syntax

[c, ifail] = d02ua(n, f)
[c, ifail] = nag_ode_bvp_ps_lin_coeffs(n, f)

Description

nag_ode_bvp_ps_lin_coeffs (d02ua) computes the coefficients cj, for j=1,2,,n+1, of the interpolating Chebyshev series
12 c1 T0 x- + c2 T1 x- + c3T2 x- ++ cn+1 Tn x- ,  
which interpolates the function fx evaluated at the Chebyshev Gauss–Lobatto points
x-r = - cos r-1 π/n ,   r=1,2,,n+1 .  
Here Tjx- denotes the Chebyshev polynomial of the first kind of degree j with argument x- defined on -1,1. In terms of your original variable, x say, the input values at which the function values are to be provided are
xr = - 12 b - a cos πr-1 /n + 1 2 b + a ,   r=1,2,,n+1 , ​  
where b and a are respectively the upper and lower ends of the range of x over which the function is required.

References

Canuto C (1988) Spectral Methods in Fluid Dynamics 502 Springer
Canuto C, Hussaini M Y, Quarteroni A and Zang T A (2006) Spectral Methods: Fundamentals in Single Domains Springer
Trefethen L N (2000) Spectral Methods in MATLAB SIAM

Parameters

Compulsory Input Parameters

1:     n int64int32nag_int scalar
n, where the number of grid points is n+1. This is also the largest order of Chebyshev polynomial in the Chebyshev series to be computed.
Constraint: n>0 and n is even.
2:     fn+1 – double array
The function values fxr, for r=1,2,,n+1.

Optional Input Parameters

None.

Output Parameters

1:     cn+1 – double array
The Chebyshev coefficients, cj, for j=1,2,,n+1.
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
Constraint: n>1.
Constraint: n is even.
   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 Chebyshev coefficients computed should be accurate to within a small multiple of machine precision.

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).

Example

See Example in nag_ode_bvp_ps_lin_solve (d02ue).
function d02ua_example


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

% On [0,4], Solve u + ux = f0(x); u(0) = 0
% where f) is such that u(x) = sin(10xcos^2x).
% Set up Chebyshev grid on [a,b]
a = 0;
b = 4;
n = int64(128);
[x, ifail] = d02uc(n, a, b);

% Get Chebyshev coeficients on grid for f0(x) = 1.
z = 10*x.*cos(x).^2;
f0 = sin(z) - 10*cos(x).^2.*cos(10*x.*cos(x).^2).*(2*x.*tan(x)-1);
[f0_c, ifail] = d02ua(n, f0);

% Set up problem definition for  u + ux = f0 [(1 1).(u ux) = f0]
f = [1, 1];
% subject to u(a) = 0  [(1 0).(u ux)(a) = 0]
y = [a];
B = [1, 0];
beta = 0;

% Solve in coefficient space using f0_c for rhs.
[B, f, uc, resid, ifail] = d02ue(...
                                 n, a, b, f0_c, B, y, beta, f);

% Transform solution and derivative back to real space.
[u,  ifail] = d02ub(...
                    n, a, b, int64(0), uc(:, 1));
[ux, ifail] = d02ub(...
                    n, a, b, int64(1), uc(:, 2));

maxerr = max(abs(u - sin(z)));
fprintf('With n = %4d, maximum error in solution = %13.2e\n',n,maxerr);

% Plot solution
fig1 = figure;
plot(x,u,x,ux);
title('Solution of u + u_x = f_0 on [0,4]');
xlabel('x');
ylabel('u(x) and u_x(x)');
legend('u','u_x','Location','Northwest');

d02ua example results

With n =  128, maximum error in solution =      5.33e-09
d02ua_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