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
|
<!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, 2014, 2015, 2016 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 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 the freedom to copy and modify this
GNU Manual." -->
<!-- Created by GNU Texinfo 5.1, http://www.gnu.org/software/texinfo/ -->
<head>
<title>GNU Scientific Library – Reference Manual: Multi-parameter regression</title>
<meta name="description" content="GNU Scientific Library – Reference Manual: Multi-parameter regression">
<meta name="keywords" content="GNU Scientific Library – Reference Manual: Multi-parameter regression">
<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="Least_002dSquares-Fitting.html#Least_002dSquares-Fitting" rel="up" title="Least-Squares Fitting">
<link href="Regularized-regression.html#Regularized-regression" rel="next" title="Regularized regression">
<link href="Linear-regression-without-a-constant-term.html#Linear-regression-without-a-constant-term" rel="previous" title="Linear regression without a constant term">
<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="Multi_002dparameter-regression"></a>
<div class="header">
<p>
Next: <a href="Regularized-regression.html#Regularized-regression" accesskey="n" rel="next">Regularized regression</a>, Previous: <a href="Linear-regression.html#Linear-regression" accesskey="p" rel="previous">Linear regression</a>, Up: <a href="Least_002dSquares-Fitting.html#Least_002dSquares-Fitting" accesskey="u" rel="up">Least-Squares Fitting</a> [<a href="Function-Index.html#Function-Index" title="Index" rel="index">Index</a>]</p>
</div>
<hr>
<a name="Multi_002dparameter-regression-1"></a>
<h3 class="section">38.3 Multi-parameter regression</h3>
<a name="index-multi_002dparameter-regression"></a>
<a name="index-fits_002c-multi_002dparameter-linear"></a>
<p>This section describes routines which perform least squares fits
to a linear model by minimizing the cost function
</p><div class="example">
<pre class="example">\chi^2 = \sum_i w_i (y_i - \sum_j X_ij c_j)^2 = || y - Xc ||_W^2
</pre></div>
<p>where <em>y</em> is a vector of <em>n</em> observations, <em>X</em> is an
<em>n</em>-by-<em>p</em> matrix of predictor variables, <em>c</em>
is a vector of the <em>p</em> unknown best-fit parameters to be estimated,
and <em>||r||_W^2 = r^T W r</em>.
The matrix <em>W = </em> diag<em>(w_1,w_2,...,w_n)</em>
defines the weights or uncertainties of the observation vector.
</p>
<p>This formulation can be used for fits to any number of functions and/or
variables by preparing the <em>n</em>-by-<em>p</em> matrix <em>X</em>
appropriately. For example, to fit to a <em>p</em>-th order polynomial in
<var>x</var>, use the following matrix,
</p>
<div class="example">
<pre class="example">X_{ij} = x_i^j
</pre></div>
<p>where the index <em>i</em> runs over the observations and the index
<em>j</em> runs from 0 to <em>p-1</em>.
</p>
<p>To fit to a set of <em>p</em> sinusoidal functions with fixed frequencies
<em>\omega_1</em>, <em>\omega_2</em>, …, <em>\omega_p</em>, use,
</p>
<div class="example">
<pre class="example">X_{ij} = sin(\omega_j x_i)
</pre></div>
<p>To fit to <em>p</em> independent variables <em>x_1</em>, <em>x_2</em>, …,
<em>x_p</em>, use,
</p>
<div class="example">
<pre class="example">X_{ij} = x_j(i)
</pre></div>
<p>where <em>x_j(i)</em> is the <em>i</em>-th value of the predictor variable
<em>x_j</em>.
</p>
<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 <em>X</em>.
</p>
<p>These functions are declared in the header file <samp>gsl_multifit.h</samp>.
</p>
<dl>
<dt><a name="index-gsl_005fmultifit_005flinear_005falloc"></a>Function: <em>gsl_multifit_linear_workspace *</em> <strong>gsl_multifit_linear_alloc</strong> <em>(const size_t <var>n</var>, const size_t <var>p</var>)</em></dt>
<dd><a name="index-gsl_005fmultifit_005flinear_005fworkspace"></a>
<p>This function allocates a workspace for fitting a model to a maximum of <var>n</var>
observations using a maximum of <var>p</var> parameters. The user may later supply
a smaller least squares system if desired. The size of the workspace is
<em>O(np + p^2)</em>.
</p></dd></dl>
<dl>
<dt><a name="index-gsl_005fmultifit_005flinear_005ffree"></a>Function: <em>void</em> <strong>gsl_multifit_linear_free</strong> <em>(gsl_multifit_linear_workspace * <var>work</var>)</em></dt>
<dd><p>This function frees the memory associated with the workspace <var>w</var>.
</p></dd></dl>
<dl>
<dt><a name="index-gsl_005fmultifit_005flinear_005fsvd"></a>Function: <em>int</em> <strong>gsl_multifit_linear_svd</strong> <em>(const gsl_matrix * <var>X</var>, gsl_multifit_linear_workspace * <var>work</var>)</em></dt>
<dd><p>This function performs a singular value decomposition of the
matrix <var>X</var> and stores the SVD factors internally in <var>work</var>.
</p></dd></dl>
<dl>
<dt><a name="index-gsl_005fmultifit_005flinear_005fbsvd"></a>Function: <em>int</em> <strong>gsl_multifit_linear_bsvd</strong> <em>(const gsl_matrix * <var>X</var>, gsl_multifit_linear_workspace * <var>work</var>)</em></dt>
<dd><p>This function performs a singular value decomposition of the
matrix <var>X</var> and stores the SVD factors internally in <var>work</var>.
The matrix <var>X</var> is first balanced by applying column scaling
factors to improve the accuracy of the singular values.
</p></dd></dl>
<dl>
<dt><a name="index-gsl_005fmultifit_005flinear"></a>Function: <em>int</em> <strong>gsl_multifit_linear</strong> <em>(const gsl_matrix * <var>X</var>, const gsl_vector * <var>y</var>, gsl_vector * <var>c</var>, gsl_matrix * <var>cov</var>, double * <var>chisq</var>, gsl_multifit_linear_workspace * <var>work</var>)</em></dt>
<dd><p>This function computes the best-fit parameters <var>c</var> of the model
<em>y = X c</em> for the observations <var>y</var> and the matrix of
predictor variables <var>X</var>, using the preallocated workspace provided
in <var>work</var>. The <em>p</em>-by-<em>p</em> variance-covariance matrix of the model parameters
<var>cov</var> is set to <em>\sigma^2 (X^T X)^{-1}</em>, where <em>\sigma</em> is
the standard deviation of the fit residuals.
The sum of squares of the residuals from the best-fit,
<em>\chi^2</em>, is returned in <var>chisq</var>. If the coefficient of
determination is desired, it can be computed from the expression
<em>R^2 = 1 - \chi^2 / TSS</em>, where the total sum of squares (TSS) of
the observations <var>y</var> may be computed from <code>gsl_stats_tss</code>.
</p>
<p>The best-fit is found by singular value decomposition of the matrix
<var>X</var> using the modified Golub-Reinsch SVD algorithm, 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.
</p></dd></dl>
<dl>
<dt><a name="index-gsl_005fmultifit_005flinear_005ftsvd"></a>Function: <em>int</em> <strong>gsl_multifit_linear_tsvd</strong> <em>(const gsl_matrix * <var>X</var>, const gsl_vector * <var>y</var>, const double <var>tol</var>, gsl_vector * <var>c</var>, gsl_matrix * <var>cov</var>, double * <var>chisq</var>, size_t * <var>rank</var>, gsl_multifit_linear_workspace * <var>work</var>)</em></dt>
<dd><p>This function computes the best-fit parameters <var>c</var> of the model
<em>y = X c</em> for the observations <var>y</var> and the matrix of
predictor variables <var>X</var>, using a truncated SVD expansion.
Singular values which satisfy <em>s_i \le tol \times s_0</em>
are discarded from the fit, where <em>s_0</em> is the largest singular value.
The <em>p</em>-by-<em>p</em> variance-covariance matrix of the model parameters
<var>cov</var> is set to <em>\sigma^2 (X^T X)^{-1}</em>, where <em>\sigma</em> is
the standard deviation of the fit residuals.
The sum of squares of the residuals from the best-fit,
<em>\chi^2</em>, is returned in <var>chisq</var>. The effective rank
(number of singular values used in solution) is returned in <var>rank</var>.
If the coefficient of
determination is desired, it can be computed from the expression
<em>R^2 = 1 - \chi^2 / TSS</em>, where the total sum of squares (TSS) of
the observations <var>y</var> may be computed from <code>gsl_stats_tss</code>.
</p></dd></dl>
<dl>
<dt><a name="index-gsl_005fmultifit_005fwlinear"></a>Function: <em>int</em> <strong>gsl_multifit_wlinear</strong> <em>(const gsl_matrix * <var>X</var>, const gsl_vector * <var>w</var>, const gsl_vector * <var>y</var>, gsl_vector * <var>c</var>, gsl_matrix * <var>cov</var>, double * <var>chisq</var>, gsl_multifit_linear_workspace * <var>work</var>)</em></dt>
<dd><p>This function computes the best-fit parameters <var>c</var> of the weighted
model <em>y = X c</em> for the observations <var>y</var> with weights <var>w</var>
and the matrix of predictor variables <var>X</var>, using the preallocated
workspace provided in <var>work</var>. The <em>p</em>-by-<em>p</em> covariance matrix of the model
parameters <var>cov</var> is computed as <em>(X^T W X)^{-1}</em>. The weighted
sum of squares of the residuals from the best-fit, <em>\chi^2</em>, is
returned in <var>chisq</var>. If the coefficient of determination is
desired, it can be computed from the expression <em>R^2 = 1 - \chi^2
/ WTSS</em>, where the weighted total sum of squares (WTSS) of the
observations <var>y</var> may be computed from <code>gsl_stats_wtss</code>.
</p></dd></dl>
<dl>
<dt><a name="index-gsl_005fmultifit_005fwlinear_005ftsvd"></a>Function: <em>int</em> <strong>gsl_multifit_wlinear_tsvd</strong> <em>(const gsl_matrix * <var>X</var>, const gsl_vector * <var>w</var>, const gsl_vector * <var>y</var>, const double <var>tol</var>, gsl_vector * <var>c</var>, gsl_matrix * <var>cov</var>, double * <var>chisq</var>, size_t * <var>rank</var>, gsl_multifit_linear_workspace * <var>work</var>)</em></dt>
<dd><p>This function computes the best-fit parameters <var>c</var> of the weighted
model <em>y = X c</em> for the observations <var>y</var> with weights <var>w</var>
and the matrix of predictor variables <var>X</var>, using a truncated SVD expansion.
Singular values which satisfy <em>s_i \le tol \times s_0</em>
are discarded from the fit, where <em>s_0</em> is the largest singular value.
The <em>p</em>-by-<em>p</em> covariance matrix of the model
parameters <var>cov</var> is computed as <em>(X^T W X)^{-1}</em>. The weighted
sum of squares of the residuals from the best-fit, <em>\chi^2</em>, is
returned in <var>chisq</var>. The effective rank of the system (number of
singular values used in the solution) is returned in <var>rank</var>.
If the coefficient of determination is
desired, it can be computed from the expression <em>R^2 = 1 - \chi^2
/ WTSS</em>, where the weighted total sum of squares (WTSS) of the
observations <var>y</var> may be computed from <code>gsl_stats_wtss</code>.
</p></dd></dl>
<dl>
<dt><a name="index-gsl_005fmultifit_005flinear_005fest"></a>Function: <em>int</em> <strong>gsl_multifit_linear_est</strong> <em>(const gsl_vector * <var>x</var>, const gsl_vector * <var>c</var>, const gsl_matrix * <var>cov</var>, double * <var>y</var>, double * <var>y_err</var>)</em></dt>
<dd><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 <em>y = x.c</em>
at the point <var>x</var>.
</p></dd></dl>
<dl>
<dt><a name="index-gsl_005fmultifit_005flinear_005fresiduals"></a>Function: <em>int</em> <strong>gsl_multifit_linear_residuals</strong> <em>(const gsl_matrix * <var>X</var>, const gsl_vector * <var>y</var>, const gsl_vector * <var>c</var>, gsl_vector * <var>r</var>)</em></dt>
<dd><p>This function computes the vector of residuals <em>r = y - X c</em> for
the observations <var>y</var>, coefficients <var>c</var> and matrix of predictor
variables <var>X</var>.
</p></dd></dl>
<dl>
<dt><a name="index-gsl_005fmultifit_005flinear_005frank"></a>Function: <em>size_t</em> <strong>gsl_multifit_linear_rank</strong> <em>(const double <var>tol</var>, const gsl_multifit_linear_workspace * <var>work</var>)</em></dt>
<dd><p>This function returns the rank of the matrix <em>X</em> which must first have its
singular value decomposition computed. The rank is computed by counting the number
of singular values <em>\sigma_j</em> which satisfy <em>\sigma_j > tol \times \sigma_0</em>,
where <em>\sigma_0</em> is the largest singular value.
</p></dd></dl>
<hr>
<div class="header">
<p>
Next: <a href="Regularized-regression.html#Regularized-regression" accesskey="n" rel="next">Regularized regression</a>, Previous: <a href="Linear-regression.html#Linear-regression" accesskey="p" rel="previous">Linear regression</a>, Up: <a href="Least_002dSquares-Fitting.html#Least_002dSquares-Fitting" accesskey="u" rel="up">Least-Squares Fitting</a> [<a href="Function-Index.html#Function-Index" title="Index" rel="index">Index</a>]</p>
</div>
</body>
</html>
|