File: Example-programs-for-Nonlinear-Least_002dSquares-Fitting.html

package info (click to toggle)
gsl-ref-html 1.16-1
  • links: PTS
  • area: non-free
  • in suites: jessie, jessie-kfreebsd, stretch
  • size: 5,816 kB
  • ctags: 4,130
  • sloc: makefile: 35
file content (322 lines) | stat: -rw-r--r-- 11,925 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
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
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
<!DOCTYPE html PUBLIC "-//W3C//DTD HTML 4.01 Transitional//EN" "http://www.w3.org/TR/html4/loose.dtd">
<html>
<!-- Copyright (C) 1996, 1997, 1998, 1999, 2000, 2001, 2002, 2003, 2004, 2005, 2006, 2007, 2008, 2009, 2010, 2011, 2012, 2013 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.3 or
any later version published by the Free Software Foundation; with no
Invariant Sections and no cover texts.  A copy of the license is
included in the section entitled "GNU Free Documentation License". -->
<!-- Created by GNU Texinfo 5.1, http://www.gnu.org/software/texinfo/ -->
<head>
<title>GNU Scientific Library &ndash; Reference Manual: Example programs for Nonlinear Least-Squares Fitting</title>

<meta name="description" content="GNU Scientific Library &ndash; Reference Manual: Example programs for Nonlinear Least-Squares Fitting">
<meta name="keywords" content="GNU Scientific Library &ndash; Reference Manual: Example programs for Nonlinear Least-Squares Fitting">
<meta name="resource-type" content="document">
<meta name="distribution" content="global">
<meta name="Generator" content="makeinfo">
<meta http-equiv="Content-Type" content="text/html; charset=utf-8">
<link href="index.html#Top" rel="start" title="Top">
<link href="Function-Index.html#Function-Index" rel="index" title="Function Index">
<link href="Nonlinear-Least_002dSquares-Fitting.html#Nonlinear-Least_002dSquares-Fitting" rel="up" title="Nonlinear Least-Squares Fitting">
<link href="References-and-Further-Reading-for-Nonlinear-Least_002dSquares-Fitting.html#References-and-Further-Reading-for-Nonlinear-Least_002dSquares-Fitting" rel="next" title="References and Further Reading for Nonlinear Least-Squares Fitting">
<link href="Computing-the-covariance-matrix-of-best-fit-parameters.html#Computing-the-covariance-matrix-of-best-fit-parameters" rel="previous" title="Computing the covariance matrix of best fit parameters">
<style type="text/css">
<!--
a.summary-letter {text-decoration: none}
blockquote.smallquotation {font-size: smaller}
div.display {margin-left: 3.2em}
div.example {margin-left: 3.2em}
div.indentedblock {margin-left: 3.2em}
div.lisp {margin-left: 3.2em}
div.smalldisplay {margin-left: 3.2em}
div.smallexample {margin-left: 3.2em}
div.smallindentedblock {margin-left: 3.2em; font-size: smaller}
div.smalllisp {margin-left: 3.2em}
kbd {font-style:oblique}
pre.display {font-family: inherit}
pre.format {font-family: inherit}
pre.menu-comment {font-family: serif}
pre.menu-preformatted {font-family: serif}
pre.smalldisplay {font-family: inherit; font-size: smaller}
pre.smallexample {font-size: smaller}
pre.smallformat {font-family: inherit; font-size: smaller}
pre.smalllisp {font-size: smaller}
span.nocodebreak {white-space:nowrap}
span.nolinebreak {white-space:nowrap}
span.roman {font-family:serif; font-weight:normal}
span.sansserif {font-family:sans-serif; font-weight:normal}
ul.no-bullet {list-style: none}
-->
</style>


</head>

<body lang="en" bgcolor="#FFFFFF" text="#000000" link="#0000FF" vlink="#800080" alink="#FF0000">
<a name="Example-programs-for-Nonlinear-Least_002dSquares-Fitting"></a>
<div class="header">
<p>
Next: <a href="References-and-Further-Reading-for-Nonlinear-Least_002dSquares-Fitting.html#References-and-Further-Reading-for-Nonlinear-Least_002dSquares-Fitting" accesskey="n" rel="next">References and Further Reading for Nonlinear Least-Squares Fitting</a>, Previous: <a href="Computing-the-covariance-matrix-of-best-fit-parameters.html#Computing-the-covariance-matrix-of-best-fit-parameters" accesskey="p" rel="previous">Computing the covariance matrix of best fit parameters</a>, Up: <a href="Nonlinear-Least_002dSquares-Fitting.html#Nonlinear-Least_002dSquares-Fitting" accesskey="u" rel="up">Nonlinear Least-Squares Fitting</a> &nbsp; [<a href="Function-Index.html#Function-Index" title="Index" rel="index">Index</a>]</p>
</div>
<hr>
<a name="Examples-29"></a>
<h3 class="section">38.11 Examples</h3>

<p>The following example program fits a weighted exponential model with
background to experimental data, <em>Y = A \exp(-\lambda t) + b</em>. The
first part of the program sets up the functions <code>expb_f</code> and
<code>expb_df</code> to calculate the model and its Jacobian.  The appropriate
fitting function is given by,
where we have chosen <em>t_i = i</em>.  The Jacobian matrix <em>J</em> is
the derivative of these functions with respect to the three parameters
(<em>A</em>, <em>\lambda</em>, <em>b</em>).  It is given by,
where <em>x_0 = A</em>, <em>x_1 = \lambda</em> and <em>x_2 = b</em>.
</p>
<div class="example">
<pre class="verbatim">/* expfit.c -- model functions for exponential + background */

struct data {
  size_t n;
  double * y;
  double * sigma;
};

int
expb_f (const gsl_vector * x, void *data, 
        gsl_vector * f)
{
  size_t n = ((struct data *)data)-&gt;n;
  double *y = ((struct data *)data)-&gt;y;
  double *sigma = ((struct data *) data)-&gt;sigma;

  double A = gsl_vector_get (x, 0);
  double lambda = gsl_vector_get (x, 1);
  double b = gsl_vector_get (x, 2);

  size_t i;

  for (i = 0; i &lt; n; i++)
    {
      /* Model Yi = A * exp(-lambda * i) + b */
      double t = i;
      double Yi = A * exp (-lambda * t) + b;
      gsl_vector_set (f, i, (Yi - y[i])/sigma[i]);
    }

  return GSL_SUCCESS;
}

int
expb_df (const gsl_vector * x, void *data, 
         gsl_matrix * J)
{
  size_t n = ((struct data *)data)-&gt;n;
  double *sigma = ((struct data *) data)-&gt;sigma;

  double A = gsl_vector_get (x, 0);
  double lambda = gsl_vector_get (x, 1);

  size_t i;

  for (i = 0; i &lt; n; i++)
    {
      /* Jacobian matrix J(i,j) = dfi / dxj, */
      /* where fi = (Yi - yi)/sigma[i],      */
      /*       Yi = A * exp(-lambda * i) + b  */
      /* and the xj are the parameters (A,lambda,b) */
      double t = i;
      double s = sigma[i];
      double e = exp(-lambda * t);
      gsl_matrix_set (J, i, 0, e/s); 
      gsl_matrix_set (J, i, 1, -t * A * e/s);
      gsl_matrix_set (J, i, 2, 1/s);
    }
  return GSL_SUCCESS;
}

int
expb_fdf (const gsl_vector * x, void *data,
          gsl_vector * f, gsl_matrix * J)
{
  expb_f (x, data, f);
  expb_df (x, data, J);

  return GSL_SUCCESS;
}
</pre></div>

<p>The main part of the program sets up a Levenberg-Marquardt solver and
some simulated random data. The data uses the known parameters
(1.0,5.0,0.1) combined with Gaussian noise (standard deviation = 0.1)
over a range of 40 timesteps. The initial guess for the parameters is
chosen as (0.0, 1.0, 0.0).
</p>
<div class="example">
<pre class="verbatim">#include &lt;stdlib.h&gt;
#include &lt;stdio.h&gt;
#include &lt;gsl/gsl_rng.h&gt;
#include &lt;gsl/gsl_randist.h&gt;
#include &lt;gsl/gsl_vector.h&gt;
#include &lt;gsl/gsl_blas.h&gt;
#include &lt;gsl/gsl_multifit_nlin.h&gt;

#include &quot;expfit.c&quot;

#define N 40

void print_state (size_t iter, gsl_multifit_fdfsolver * s);

int
main (void)
{
  const gsl_multifit_fdfsolver_type *T;
  gsl_multifit_fdfsolver *s;
  int status;
  unsigned int i, iter = 0;
  const size_t n = N;
  const size_t p = 3;

  gsl_matrix *covar = gsl_matrix_alloc (p, p);
  double y[N], sigma[N];
  struct data d = { n, y, sigma};
  gsl_multifit_function_fdf f;
  double x_init[3] = { 1.0, 0.0, 0.0 };
  gsl_vector_view x = gsl_vector_view_array (x_init, p);
  const gsl_rng_type * type;
  gsl_rng * r;

  gsl_rng_env_setup();

  type = gsl_rng_default;
  r = gsl_rng_alloc (type);

  f.f = &amp;expb_f;
  f.df = &amp;expb_df;
  f.fdf = &amp;expb_fdf;
  f.n = n;
  f.p = p;
  f.params = &amp;d;

  /* This is the data to be fitted */

  for (i = 0; i &lt; n; i++)
    {
      double t = i;
      y[i] = 1.0 + 5 * exp (-0.1 * t) 
                 + gsl_ran_gaussian (r, 0.1);
      sigma[i] = 0.1;
      printf (&quot;data: %u %g %g\n&quot;, i, y[i], sigma[i]);
    };

  T = gsl_multifit_fdfsolver_lmsder;
  s = gsl_multifit_fdfsolver_alloc (T, n, p);
  gsl_multifit_fdfsolver_set (s, &amp;f, &amp;x.vector);

  print_state (iter, s);

  do
    {
      iter++;
      status = gsl_multifit_fdfsolver_iterate (s);

      printf (&quot;status = %s\n&quot;, gsl_strerror (status));

      print_state (iter, s);

      if (status)
        break;

      status = gsl_multifit_test_delta (s-&gt;dx, s-&gt;x,
                                        1e-4, 1e-4);
    }
  while (status == GSL_CONTINUE &amp;&amp; iter &lt; 500);

  gsl_multifit_covar (s-&gt;J, 0.0, covar);

#define FIT(i) gsl_vector_get(s-&gt;x, i)
#define ERR(i) sqrt(gsl_matrix_get(covar,i,i))

  { 
    double chi = gsl_blas_dnrm2(s-&gt;f);
    double dof = n - p;
    double c = GSL_MAX_DBL(1, chi / sqrt(dof)); 

    printf(&quot;chisq/dof = %g\n&quot;,  pow(chi, 2.0) / dof);

    printf (&quot;A      = %.5f +/- %.5f\n&quot;, FIT(0), c*ERR(0));
    printf (&quot;lambda = %.5f +/- %.5f\n&quot;, FIT(1), c*ERR(1));
    printf (&quot;b      = %.5f +/- %.5f\n&quot;, FIT(2), c*ERR(2));
  }

  printf (&quot;status = %s\n&quot;, gsl_strerror (status));

  gsl_multifit_fdfsolver_free (s);
  gsl_matrix_free (covar);
  gsl_rng_free (r);
  return 0;
}

void
print_state (size_t iter, gsl_multifit_fdfsolver * s)
{
  printf (&quot;iter: %3u x = % 15.8f % 15.8f % 15.8f &quot;
          &quot;|f(x)| = %g\n&quot;,
          iter,
          gsl_vector_get (s-&gt;x, 0), 
          gsl_vector_get (s-&gt;x, 1),
          gsl_vector_get (s-&gt;x, 2), 
          gsl_blas_dnrm2 (s-&gt;f));
}
</pre></div>

<p>The iteration terminates when the change in x is smaller than 0.0001, as
both an absolute and relative change.  Here are the results of running
the program:
</p>
<div class="smallexample">
<pre class="smallexample">iter: 0 x=1.00000000 0.00000000 0.00000000 |f(x)|=117.349
status=success
iter: 1 x=1.64659312 0.01814772 0.64659312 |f(x)|=76.4578
status=success
iter: 2 x=2.85876037 0.08092095 1.44796363 |f(x)|=37.6838
status=success
iter: 3 x=4.94899512 0.11942928 1.09457665 |f(x)|=9.58079
status=success
iter: 4 x=5.02175572 0.10287787 1.03388354 |f(x)|=5.63049
status=success
iter: 5 x=5.04520433 0.10405523 1.01941607 |f(x)|=5.44398
status=success
iter: 6 x=5.04535782 0.10404906 1.01924871 |f(x)|=5.44397
chisq/dof = 0.800996
A      = 5.04536 +/- 0.06028
lambda = 0.10405 +/- 0.00316
b      = 1.01925 +/- 0.03782
status = success
</pre></div>

<p>The approximate values of the parameters are found correctly, and the
chi-squared value indicates a good fit (the chi-squared per degree of
freedom is approximately 1).  In this case the errors on the parameters
can be estimated from the square roots of the diagonal elements of the
covariance matrix.  
</p>
<p>If the chi-squared value shows a poor fit (i.e. <em>chi^2/dof &gt;&gt; 1</em>) then the error estimates obtained from the
covariance matrix will be too small.  In the example program the error estimates
are multiplied by <em>\sqrt{\chi^2/dof}</em> in this case, a common way of increasing the
errors for a poor fit.  Note that a poor fit will result from the use
an inappropriate model, and the scaled error estimates may then
be outside the range of validity for Gaussian errors.
</p>

<hr>
<div class="header">
<p>
Next: <a href="References-and-Further-Reading-for-Nonlinear-Least_002dSquares-Fitting.html#References-and-Further-Reading-for-Nonlinear-Least_002dSquares-Fitting" accesskey="n" rel="next">References and Further Reading for Nonlinear Least-Squares Fitting</a>, Previous: <a href="Computing-the-covariance-matrix-of-best-fit-parameters.html#Computing-the-covariance-matrix-of-best-fit-parameters" accesskey="p" rel="previous">Computing the covariance matrix of best fit parameters</a>, Up: <a href="Nonlinear-Least_002dSquares-Fitting.html#Nonlinear-Least_002dSquares-Fitting" accesskey="u" rel="up">Nonlinear Least-Squares Fitting</a> &nbsp; [<a href="Function-Index.html#Function-Index" title="Index" rel="index">Index</a>]</p>
</div>



</body>
</html>