~ubuntu-branches/ubuntu/dapper/gsl-ref-html/dapper

1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
<HTML>
<HEAD>
<!-- This HTML file has been created by texi2html 1.54+ (gsl)
     from /home/bjg/gsl.redhat/doc/gsl-ref.texi on 14 September 2005 -->

<TITLE>GNU Scientific Library -- Reference Manual - Chebyshev Approximations</TITLE>
<link href="gsl-ref_29.html" rel=Next>
<link href="gsl-ref_27.html" rel=Previous>
<link href="gsl-ref_toc.html" rel=ToC>

</HEAD>
<BODY>
<p>Go to the <A HREF="gsl-ref_1.html">first</A>, <A HREF="gsl-ref_27.html">previous</A>, <A HREF="gsl-ref_29.html">next</A>, <A HREF="gsl-ref_50.html">last</A> section, <A HREF="gsl-ref_toc.html">table of contents</A>.
<P><HR><P>


<H1><A NAME="SEC411" HREF="gsl-ref_toc.html#TOC411">Chebyshev Approximations</A></H1>
<P>
<A NAME="IDX1993"></A>
<A NAME="IDX1994"></A>
<A NAME="IDX1995"></A>

</P>
<P>
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 &#38; Stegun, Chapter 22. 

</P>
<P>
The functions described in this chapter are declared in the header file
<TT>`gsl_chebyshev.h'</TT>.

</P>



<H2><A NAME="SEC412" HREF="gsl-ref_toc.html#TOC412">Definitions</A></H2>

<P>
A Chebyshev series  is stored using the following structure,

</P>

<PRE>
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_struct
</PRE>

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

</P>

<PRE>
f(x) = (c_0 / 2) + \sum_{n=1} c_n T_n(x)
</PRE>

<P>
which is needed when accessing the coefficients directly.

</P>


<H2><A NAME="SEC413" HREF="gsl-ref_toc.html#TOC413">Creation and Calculation of Chebyshev Series</A></H2>

<P>
<DL>
<DT><U>Function:</U> gsl_cheb_series * <B>gsl_cheb_alloc</B> <I>(const size_t <VAR>n</VAR>)</I>
<DD><A NAME="IDX1996"></A>
This function allocates space for a Chebyshev series of order <VAR>n</VAR>
and returns a pointer to a new <CODE>gsl_cheb_series</CODE> struct.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> void <B>gsl_cheb_free</B> <I>(gsl_cheb_series * <VAR>cs</VAR>)</I>
<DD><A NAME="IDX1997"></A>
This function frees a previously allocated Chebyshev series <VAR>cs</VAR>.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_cheb_init</B> <I>(gsl_cheb_series * <VAR>cs</VAR>, const gsl_function * <VAR>f</VAR>, const double <VAR>a</VAR>, const double <VAR>b</VAR>)</I>
<DD><A NAME="IDX1998"></A>
This function computes the Chebyshev approximation <VAR>cs</VAR> for the
function <VAR>f</VAR> 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.
</DL>

</P>


<H2><A NAME="SEC414" HREF="gsl-ref_toc.html#TOC414">Chebyshev Series Evaluation</A></H2>

<P>
<DL>
<DT><U>Function:</U> double <B>gsl_cheb_eval</B> <I>(const gsl_cheb_series * <VAR>cs</VAR>, double <VAR>x</VAR>)</I>
<DD><A NAME="IDX1999"></A>
This function evaluates the Chebyshev series <VAR>cs</VAR> at a given point <VAR>x</VAR>.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_cheb_eval_err</B> <I>(const gsl_cheb_series * <VAR>cs</VAR>, const double <VAR>x</VAR>, double * <VAR>result</VAR>, double * <VAR>abserr</VAR>)</I>
<DD><A NAME="IDX2000"></A>
This function computes the Chebyshev series <VAR>cs</VAR> at a given point
<VAR>x</VAR>, estimating both the series <VAR>result</VAR> and its absolute error
<VAR>abserr</VAR>.  The error estimate is made from the first neglected term
in the series.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> double <B>gsl_cheb_eval_n</B> <I>(const gsl_cheb_series * <VAR>cs</VAR>, size_t <VAR>order</VAR>, double <VAR>x</VAR>)</I>
<DD><A NAME="IDX2001"></A>
This function evaluates the Chebyshev series <VAR>cs</VAR> at a given point
<VAR>n</VAR>, to (at most) the given order <VAR>order</VAR>.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_cheb_eval_n_err</B> <I>(const gsl_cheb_series * <VAR>cs</VAR>, const size_t <VAR>order</VAR>, const double <VAR>x</VAR>, double * <VAR>result</VAR>, double * <VAR>abserr</VAR>)</I>
<DD><A NAME="IDX2002"></A>
This function evaluates a Chebyshev series <VAR>cs</VAR> at a given point
<VAR>x</VAR>, estimating both the series <VAR>result</VAR> and its absolute error
<VAR>abserr</VAR>, to (at most) the given order <VAR>order</VAR>.  The error
estimate is made from the first neglected term in the series.
</DL>

</P>



<H2><A NAME="SEC415" HREF="gsl-ref_toc.html#TOC415">Derivatives and Integrals</A></H2>

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

</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_cheb_calc_deriv</B> <I>(gsl_cheb_series * <VAR>deriv</VAR>, const gsl_cheb_series * <VAR>cs</VAR>)</I>
<DD><A NAME="IDX2003"></A>
This function computes the derivative of the series <VAR>cs</VAR>, storing
the derivative coefficients in the previously allocated <VAR>deriv</VAR>.
The two series <VAR>cs</VAR> and <VAR>deriv</VAR> must have been allocated with
the same order.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_cheb_calc_integ</B> <I>(gsl_cheb_series * <VAR>integ</VAR>, const gsl_cheb_series * <VAR>cs</VAR>)</I>
<DD><A NAME="IDX2004"></A>
This function computes the integral of the series <VAR>cs</VAR>, storing the
integral coefficients in the previously allocated <VAR>integ</VAR>.  The two
series <VAR>cs</VAR> and <VAR>integ</VAR> 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 <VAR>a</VAR>.
</DL>

</P>


<H2><A NAME="SEC416" HREF="gsl-ref_toc.html#TOC416">Examples</A></H2>

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

</P>

<PRE>
#include &#60;stdio.h&#62;
#include &#60;gsl/gsl_math.h&#62;
#include &#60;gsl/gsl_chebyshev.h&#62;

double
f (double x, void *p)
{
  if (x &#60; 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, &#38;F, 0.0, 1.0);

  for (i = 0; i &#60; 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 (&#38;F, x), r10, r40);
    }

  gsl_cheb_free (cs);

  return 0;
}
</PRE>

<P>
The output from the program gives the original function, 10-th order
approximation and 40-th order approximation, all sampled at intervals of
0.001 in x.

</P>



<H2><A NAME="SEC417" HREF="gsl-ref_toc.html#TOC417">References and Further Reading</A></H2>

<P>
The following paper describes the use of Chebyshev series,

</P>

<UL>
<LI>

R. Broucke, "Ten Subroutines for the Manipulation of Chebyshev Series
[C1] (Algorithm 446)". <CITE>Communications of the ACM</CITE> 16(4), 254--256
(1973)
</UL>

<P><HR><P>
<p>Go to the <A HREF="gsl-ref_1.html">first</A>, <A HREF="gsl-ref_27.html">previous</A>, <A HREF="gsl-ref_29.html">next</A>, <A HREF="gsl-ref_50.html">last</A> section, <A HREF="gsl-ref_toc.html">table of contents</A>.
</BODY>
</HTML>