File: Example-programs-for-B_002dsplines.html

package info (click to toggle)
gsl-ref-html 1.10-1
  • links: PTS
  • area: main
  • in suites: lenny
  • size: 4,496 kB
  • ctags: 2,955
  • sloc: makefile: 33
file content (181 lines) | stat: -rw-r--r-- 6,507 bytes parent folder | download
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
<html lang="en">
<head>
<title>Example programs for B-splines - GNU Scientific Library -- Reference Manual</title>
<meta http-equiv="Content-Type" content="text/html">
<meta name="description" content="GNU Scientific Library -- Reference Manual">
<meta name="generator" content="makeinfo 4.8">
<link title="Top" rel="start" href="index.html#Top">
<link rel="up" href="Basis-Splines.html" title="Basis Splines">
<link rel="prev" href="Evaluation-of-B_002dspline-basis-functions.html" title="Evaluation of B-spline basis functions">
<link rel="next" href="References-and-Further-Reading.html" title="References and Further Reading">
<link href="http://www.gnu.org/software/texinfo/" rel="generator-home" title="Texinfo Homepage">
<!--
Copyright (C) 1996, 1997, 1998, 1999, 2000, 2001, 2002, 2003, 2004, 2005, 2006, 2007 The GSL Team.

Permission is granted to copy, distribute and/or modify this document
under the terms of the GNU Free Documentation License, Version 1.2 or
any later version published by the Free Software Foundation; with the
Invariant Sections being ``GNU General Public License'' and ``Free Software
Needs Free Documentation'', the Front-Cover text being ``A GNU Manual'',
and with the Back-Cover Text being (a) (see below).  A copy of the
license is included in the section entitled ``GNU Free Documentation
License''.

(a) The Back-Cover Text is: ``You have freedom to copy and modify this
GNU Manual, like GNU software.''-->
<meta http-equiv="Content-Style-Type" content="text/css">
<style type="text/css"><!--
  pre.display { font-family:inherit }
  pre.format  { font-family:inherit }
  pre.smalldisplay { font-family:inherit; font-size:smaller }
  pre.smallformat  { font-family:inherit; font-size:smaller }
  pre.smallexample { font-size:smaller }
  pre.smalllisp    { font-size:smaller }
  span.sc    { font-variant:small-caps }
  span.roman { font-family:serif; font-weight:normal; } 
  span.sansserif { font-family:sans-serif; font-weight:normal; } 
--></style>
</head>
<body>
<div class="node">
<p>
<a name="Example-programs-for-B-splines"></a>
<a name="Example-programs-for-B_002dsplines"></a>
Next:&nbsp;<a rel="next" accesskey="n" href="References-and-Further-Reading.html">References and Further Reading</a>,
Previous:&nbsp;<a rel="previous" accesskey="p" href="Evaluation-of-B_002dspline-basis-functions.html">Evaluation of B-spline basis functions</a>,
Up:&nbsp;<a rel="up" accesskey="u" href="Basis-Splines.html">Basis Splines</a>
<hr>
</div>

<h3 class="section">38.5 Example programs for B-splines</h3>

<p><a name="index-basis-splines_002c-examples-2434"></a>
The following program computes a linear least squares fit to data using
cubic B-spline basis functions with uniform breakpoints. The data is
generated from the curve y(x) = \cos(x) \exp(-0.1 x) on
[0, 15] with gaussian noise added.

<pre class="example"><pre class="verbatim">     #include &lt;stdio.h>
     #include &lt;stdlib.h>
     #include &lt;math.h>
     #include &lt;gsl/gsl_bspline.h>
     #include &lt;gsl/gsl_multifit.h>
     #include &lt;gsl/gsl_rng.h>
     #include &lt;gsl/gsl_randist.h>
     
     /* number of data points to fit */
     #define N        200
     
     /* number of fit coefficients */
     #define NCOEFFS  8
     
     /* nbreak = ncoeffs + 2 - k = ncoeffs - 2 since k = 4 */
     #define NBREAK   (NCOEFFS - 2)
     
     int
     main (void)
     {
       const size_t n = N;
       const size_t ncoeffs = NCOEFFS;
       const size_t nbreak = NBREAK;
       size_t i, j;
       gsl_bspline_workspace *bw;
       gsl_vector *B;
       double dy;
       gsl_rng *r;
       gsl_vector *c, *w;
       gsl_vector *x, *y;
       gsl_matrix *X, *cov;
       gsl_multifit_linear_workspace *mw;
       double chisq;
     
       gsl_rng_env_setup();
       r = gsl_rng_alloc(gsl_rng_default);
     
       /* allocate a cubic bspline workspace (k = 4) */
       bw = gsl_bspline_alloc(4, nbreak);
       B = gsl_vector_alloc(ncoeffs);
     
       x = gsl_vector_alloc(n);
       y = gsl_vector_alloc(n);
       X = gsl_matrix_alloc(n, ncoeffs);
       c = gsl_vector_alloc(ncoeffs);
       w = gsl_vector_alloc(n);
       cov = gsl_matrix_alloc(ncoeffs, ncoeffs);
       mw = gsl_multifit_linear_alloc(n, ncoeffs);
     
       printf("#m=0,S=0\n");
       /* this is the data to be fitted */
       for (i = 0; i &lt; n; ++i)
         {
           double sigma;
           double xi = (15.0/(N-1)) * i;
           double yi = cos(xi) * exp(-0.1 * xi);
     
           sigma = 0.1;
           dy = gsl_ran_gaussian(r, sigma);
           yi += dy;
     
           gsl_vector_set(x, i, xi);
           gsl_vector_set(y, i, yi);
           gsl_vector_set(w, i, 1.0 / (sigma*sigma));
     
           printf("%f %f\n", xi, yi);
         }
     
       /* use uniform breakpoints on [0, 15] */
       gsl_bspline_knots_uniform(0.0, 15.0, bw);
     
       /* construct the fit matrix X */
       for (i = 0; i &lt; n; ++i)
         {
           double xi = gsl_vector_get(x, i);
     
           /* compute B_j(xi) for all j */
           gsl_bspline_eval(xi, B, bw);
     
           /* fill in row i of X */
           for (j = 0; j &lt; ncoeffs; ++j)
             {
               double Bj = gsl_vector_get(B, j);
               gsl_matrix_set(X, i, j, Bj);
             }
         }
     
       /* do the fit */
       gsl_multifit_wlinear(X, w, y, c, cov, &amp;chisq, mw);
     
       /* output the smoothed curve */
       {
         double xi, yi, yerr;
     
         printf("#m=1,S=0\n");
         for (xi = 0.0; xi &lt; 15.0; xi += 0.1)
           {
             gsl_bspline_eval(xi, B, bw);
             gsl_multifit_linear_est(B, c, cov, &amp;yi, &amp;yerr);
             printf("%f %f\n", xi, yi);
           }
       }
     
       gsl_rng_free(r);
       gsl_bspline_free(bw);
       gsl_vector_free(B);
       gsl_vector_free(x);
       gsl_vector_free(y);
       gsl_matrix_free(X);
       gsl_vector_free(c);
       gsl_vector_free(w);
       gsl_matrix_free(cov);
       gsl_multifit_linear_free(mw);
     
       return 0;
     } /* main() */
</pre></pre>
   <p>The output can be plotted with <span class="sc">gnu</span> <code>graph</code>.

<pre class="example">     $ ./a.out &gt; bspline.dat
     $ graph -T ps -X x -Y y -x 0 15 -y -1 1.3 &lt; bspline.dat &gt; bspline.ps
</pre>
<hr>The GNU Scientific Library - a free numerical library licensed under the GNU GPL<br>Back to the <a href="/software/gsl/">GNU Scientific Library Homepage</a></body></html>