File: Multi_002dparameter-fitting.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 (149 lines) | stat: -rw-r--r-- 8,778 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
<html lang="en">
<head>
<title>Multi-parameter fitting - 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="Least_002dSquares-Fitting.html" title="Least-Squares Fitting">
<link rel="prev" href="Linear-fitting-without-a-constant-term.html" title="Linear fitting without a constant term">
<link rel="next" href="Fitting-Examples.html" title="Fitting Examples">
<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="Multi-parameter-fitting"></a>
<a name="Multi_002dparameter-fitting"></a>
Next:&nbsp;<a rel="next" accesskey="n" href="Fitting-Examples.html">Fitting Examples</a>,
Previous:&nbsp;<a rel="previous" accesskey="p" href="Linear-fitting-without-a-constant-term.html">Linear fitting without a constant term</a>,
Up:&nbsp;<a rel="up" accesskey="u" href="Least_002dSquares-Fitting.html">Least-Squares Fitting</a>
<hr>
</div>

<h3 class="section">36.4 Multi-parameter fitting</h3>

<p><a name="index-multi_002dparameter-regression-2381"></a><a name="index-fits_002c-multi_002dparameter-linear-2382"></a>The functions described in this section perform least-squares fits to a
general linear model, y = X c where y is a vector of
n observations, X is an n by p matrix of
predictor variables, and the elements of the vector c are the p unknown best-fit parameters which are to be estimated.  The chi-squared value is given by <!-- {$\chi^2 = \sum_i w_i (y_i - \sum_j X_{ij} c_j)^2$} -->
\chi^2 = \sum_i w_i (y_i - \sum_j X_{ij} c_j)^2.

   <p>This formulation can be used for fits to any number of functions and/or
variables by preparing the n-by-p matrix X
appropriately.  For example, to fit to a p-th order polynomial in
<var>x</var>, use the following matrix,

<pre class="example">     X_{ij} = x_i^j
</pre>
   <p class="noindent">where the index i runs over the observations and the index
j runs from 0 to p-1.

   <p>To fit to a set of p sinusoidal functions with fixed frequencies
\omega_1, \omega_2, <small class="dots">...</small>, \omega_p, use,

<pre class="example">     X_{ij} = sin(\omega_j x_i)
</pre>
   <p class="noindent">To fit to p independent variables x_1, x_2, <small class="dots">...</small>,
x_p, use,

<pre class="example">     X_{ij} = x_j(i)
</pre>
   <p class="noindent">where x_j(i) is the i-th value of the predictor variable
x_j.

   <p>The functions described in this section are declared in the header file
<samp><span class="file">gsl_multifit.h</span></samp>.

   <p>The solution of the general linear least-squares system requires an
additional working space for intermediate results, such as the singular
value decomposition of the matrix X.

<div class="defun">
&mdash; Function: gsl_multifit_linear_workspace * <b>gsl_multifit_linear_alloc</b> (<var>size_t n, size_t p</var>)<var><a name="index-gsl_005fmultifit_005flinear_005falloc-2383"></a></var><br>
<blockquote><p>This function allocates a workspace for fitting a model to <var>n</var>
observations using <var>p</var> parameters. 
</p></blockquote></div>

<div class="defun">
&mdash; Function: void <b>gsl_multifit_linear_free</b> (<var>gsl_multifit_linear_workspace * work</var>)<var><a name="index-gsl_005fmultifit_005flinear_005ffree-2384"></a></var><br>
<blockquote><p>This function frees the memory associated with the workspace <var>w</var>. 
</p></blockquote></div>

<div class="defun">
&mdash; Function: int <b>gsl_multifit_linear</b> (<var>const gsl_matrix * X, const gsl_vector * y, gsl_vector * c, gsl_matrix * cov, double * chisq, gsl_multifit_linear_workspace * work</var>)<var><a name="index-gsl_005fmultifit_005flinear-2385"></a></var><br>
&mdash; Function: int <b>gsl_multifit_linear_svd</b> (<var>const gsl_matrix * X, const gsl_vector * y, double tol, size_t * rank, gsl_vector * c, gsl_matrix * cov, double * chisq, gsl_multifit_linear_workspace * work</var>)<var><a name="index-gsl_005fmultifit_005flinear_005fsvd-2386"></a></var><br>
<blockquote><p>These functions compute the best-fit parameters <var>c</var> of the model
y = X c for the observations <var>y</var> and the matrix of predictor
variables <var>X</var>.  The variance-covariance matrix of the model
parameters <var>cov</var> is estimated from the scatter of the observations
about the best-fit.  The sum of squares of the residuals from the
best-fit, \chi^2, is returned in <var>chisq</var>.

        <p>The best-fit is found by singular value decomposition of the matrix
<var>X</var> using the preallocated workspace provided in <var>work</var>. The
modified Golub-Reinsch SVD algorithm is used, with column scaling to
improve the accuracy of the singular values. Any components which have
zero singular value (to machine precision) are discarded from the fit. 
In the second form of the function the components are discarded if the
ratio of singular values s_i/s_0 falls below the user-specified
tolerance <var>tol</var>, and the effective rank is returned in <var>rank</var>. 
</p></blockquote></div>

<div class="defun">
&mdash; Function: int <b>gsl_multifit_wlinear</b> (<var>const gsl_matrix * X, const gsl_vector * w, const gsl_vector * y, gsl_vector * c, gsl_matrix * cov, double * chisq, gsl_multifit_linear_workspace * work</var>)<var><a name="index-gsl_005fmultifit_005fwlinear-2387"></a></var><br>
&mdash; Function: int <b>gsl_multifit_wlinear_svd</b> (<var>const gsl_matrix * X, const gsl_vector * w, const gsl_vector * y, double tol, size_t * rank, gsl_vector * c, gsl_matrix * cov, double * chisq, gsl_multifit_linear_workspace * work</var>)<var><a name="index-gsl_005fmultifit_005fwlinear_005fsvd-2388"></a></var><br>
<blockquote>
<p>This function computes the best-fit parameters <var>c</var> of the weighted
model y = X c for the observations <var>y</var> with weights <var>w</var>
and the matrix of predictor variables <var>X</var>.  The covariance matrix of
the model parameters <var>cov</var> is computed with the given weights.  The
weighted sum of squares of the residuals from the best-fit,
\chi^2, is returned in <var>chisq</var>.

        <p>The best-fit is found by singular value decomposition of the matrix
<var>X</var> using the preallocated workspace provided in <var>work</var>. Any
components which have zero singular value (to machine precision) are
discarded from the fit.  In the second form of the function the
components are discarded if the ratio of singular values s_i/s_0
falls below the user-specified tolerance <var>tol</var>, and the effective
rank is returned in <var>rank</var>. 
</p></blockquote></div>

<div class="defun">
&mdash; Function: int <b>gsl_multifit_linear_est</b> (<var>const gsl_vector * x, const gsl_vector * c, const gsl_matrix * cov, double * y, double * y_err</var>)<var><a name="index-gsl_005fmultifit_005flinear_005fest-2389"></a></var><br>
<blockquote><p>This function uses the best-fit multilinear regression coefficients
<var>c</var> and their covariance matrix
<var>cov</var> to compute the fitted function value
<var>y</var> and its standard deviation <var>y_err</var> for the model y = x.c
at the point <var>x</var>. 
</p></blockquote></div>

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