Chebyshev Approximations

This chapter describes routines for computing Chebyshev approximations to univariate functions. A Chebyshev approximation is a truncation of the series f(x) = \sum c_n T_n(x), where the Chebyshev polynomials T_n(x) = \cos(n \arccos x) provide an orthogonal basis of polynomials on the interval [-1,1] with the weight function 1 / \sqrt{1-x^2}. The first few Chebyshev polynomials are, T_0(x) = 1, T_1(x) = x, T_2(x) = 2 x^2 - 1. For further information see Abramowitz & Stegun, Chapter 22.

The functions described in this chapter are declared in the header file gsl_chebyshev.h.

Definitions

type gsl_cheb_series

A Chebyshev series is stored using the following structure:

typedef struct
{
  double * c;   /* coefficients  c[0] .. c[order] */
  int order;    /* order of expansion             */
  double a;     /* lower interval point           */
  double b;     /* upper interval point           */
  ...
} gsl_cheb_series

The approximation is made over the range [a,b] using order + 1 terms, including the coefficient c[0]. The series is computed using the following convention,

f(x) = {c_0 \over 2} + \sum_{n=1} c_n T_n(x)

which is needed when accessing the coefficients directly.

Creation and Calculation of Chebyshev Series

gsl_cheb_series *gsl_cheb_alloc(const size_t n)

This function allocates space for a Chebyshev series of order n and returns a pointer to a new gsl_cheb_series struct.

void gsl_cheb_free(gsl_cheb_series *cs)

This function frees a previously allocated Chebyshev series cs.

int gsl_cheb_init(gsl_cheb_series *cs, const gsl_function *f, const double a, const double b)

This function computes the Chebyshev approximation cs for the function f over the range (a,b) to the previously specified order. The computation of the Chebyshev approximation is an O(n^2) process, and requires n function evaluations.

Auxiliary Functions

The following functions provide information about an existing Chebyshev series.

size_t gsl_cheb_order(const gsl_cheb_series *cs)

This function returns the order of Chebyshev series cs.

size_t gsl_cheb_size(const gsl_cheb_series *cs)
double *gsl_cheb_coeffs(const gsl_cheb_series *cs)

These functions return the size of the Chebyshev coefficient array c[] and a pointer to its location in memory for the Chebyshev series cs.

Chebyshev Series Evaluation

double gsl_cheb_eval(const gsl_cheb_series *cs, double x)

This function evaluates the Chebyshev series cs at a given point x.

int gsl_cheb_eval_err(const gsl_cheb_series *cs, const double x, double *result, double *abserr)

This function computes the Chebyshev series cs at a given point x, estimating both the series result and its absolute error abserr. The error estimate is made from the first neglected term in the series.

double gsl_cheb_eval_n(const gsl_cheb_series *cs, size_t order, double x)

This function evaluates the Chebyshev series cs at a given point x, to (at most) the given order order.

int gsl_cheb_eval_n_err(const gsl_cheb_series *cs, const size_t order, const double x, double *result, double *abserr)

This function evaluates a Chebyshev series cs at a given point x, estimating both the series result and its absolute error abserr, to (at most) the given order order. The error estimate is made from the first neglected term in the series.

Derivatives and Integrals

The following functions allow a Chebyshev series to be differentiated or integrated, producing a new Chebyshev series. Note that the error estimate produced by evaluating the derivative series will be underestimated due to the contribution of higher order terms being neglected.

int gsl_cheb_calc_deriv(gsl_cheb_series *deriv, const gsl_cheb_series *cs)

This function computes the derivative of the series cs, storing the derivative coefficients in the previously allocated deriv. The two series cs and deriv must have been allocated with the same order.

int gsl_cheb_calc_integ(gsl_cheb_series *integ, const gsl_cheb_series *cs)

This function computes the integral of the series cs, storing the integral coefficients in the previously allocated integ. The two series cs and integ must have been allocated with the same order. The lower limit of the integration is taken to be the left hand end of the range a.

Examples

The following example program computes Chebyshev approximations to a step function. This is an extremely difficult approximation to make, due to the discontinuity, and was chosen as an example where approximation error is visible. For smooth functions the Chebyshev approximation converges extremely rapidly and errors would not be visible.

#include <stdio.h>
#include <gsl/gsl_math.h>
#include <gsl/gsl_chebyshev.h>

double
f (double x, void *p)
{
  (void)(p); /* avoid unused parameter warning */

  if (x < 0.5)
    return 0.25;
  else
    return 0.75;
}

int
main (void)
{
  int i, n = 10000;

  gsl_cheb_series *cs = gsl_cheb_alloc (40);

  gsl_function F;

  F.function = f;
  F.params = 0;

  gsl_cheb_init (cs, &F, 0.0, 1.0);

  for (i = 0; i < n; i++)
    {
      double x = i / (double)n;
      double r10 = gsl_cheb_eval_n (cs, 10, x);
      double r40 = gsl_cheb_eval (cs, x);
      printf ("%g %g %g %g\n",
              x, GSL_FN_EVAL (&F, x), r10, r40);
    }

  gsl_cheb_free (cs);

  return 0;
}

Fig. 26 shows output from the program with the original function, 10-th order approximation and 40-th order approximation, all sampled at intervals of 0.001 in x.

_images/cheb.png

Fig. 26 Chebyshev approximations to a step function

References and Further Reading

The following paper describes the use of Chebyshev series,

  • R. Broucke, “Ten Subroutines for the Manipulation of Chebyshev Series [C1] (Algorithm 446)”. Communications of the ACM 16(4), 254–256 (1973)