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
  
     | 
    
      <html lang="en">
<head>
<title>Linear fitting without a constant term - 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-regression.html" title="Linear regression">
<link rel="next" href="Multi_002dparameter-fitting.html" title="Multi-parameter fitting">
<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="Linear-fitting-without-a-constant-term"></a>
Next: <a rel="next" accesskey="n" href="Multi_002dparameter-fitting.html">Multi-parameter fitting</a>,
Previous: <a rel="previous" accesskey="p" href="Linear-regression.html">Linear regression</a>,
Up: <a rel="up" accesskey="u" href="Least_002dSquares-Fitting.html">Least-Squares Fitting</a>
<hr>
</div>
<h3 class="section">36.3 Linear fitting without a constant term</h3>
<p>The functions described in this section can be used to perform
least-squares fits to a straight line model without a constant term,
Y = c_1 X.
<div class="defun">
— Function: int <b>gsl_fit_mul</b> (<var>const double * x, const size_t xstride, const double * y, const size_t ystride, size_t n, double * c1, double * cov11, double * sumsq</var>)<var><a name="index-gsl_005ffit_005fmul-2378"></a></var><br>
<blockquote><p>This function computes the best-fit linear regression coefficient
<var>c1</var> of the model Y = c_1 X for the datasets (<var>x</var>,
<var>y</var>), two vectors of length <var>n</var> with strides <var>xstride</var> and
<var>ystride</var>.  The errors on <var>y</var> are assumed unknown so the
variance of the parameter <var>c1</var> is estimated from
the scatter of the points around the best-fit line and returned via the
parameter <var>cov11</var>.  The sum of squares of the residuals from the
best-fit line is returned in <var>sumsq</var>. 
</p></blockquote></div>
<div class="defun">
— Function: int <b>gsl_fit_wmul</b> (<var>const double * x, const size_t xstride, const double * w, const size_t wstride, const double * y, const size_t ystride, size_t n, double * c1, double * cov11, double * sumsq</var>)<var><a name="index-gsl_005ffit_005fwmul-2379"></a></var><br>
<blockquote><p>This function computes the best-fit linear regression coefficient
<var>c1</var> of the model Y = c_1 X for the weighted datasets
(<var>x</var>, <var>y</var>), two vectors of length <var>n</var> with strides
<var>xstride</var> and <var>ystride</var>.  The vector <var>w</var>, of length <var>n</var>
and stride <var>wstride</var>, specifies the weight of each datapoint. The
weight is the reciprocal of the variance for each datapoint in <var>y</var>.
        <p>The variance of the parameter <var>c1</var> is computed using the weights
and returned via the parameter <var>cov11</var>.  The weighted sum of
squares of the residuals from the best-fit line, \chi^2, is
returned in <var>chisq</var>. 
</p></blockquote></div>
<div class="defun">
— Function: int <b>gsl_fit_mul_est</b> (<var>double x, double c1, double cov11, double * y, double * y_err</var>)<var><a name="index-gsl_005ffit_005fmul_005fest-2380"></a></var><br>
<blockquote><p>This function uses the best-fit linear regression coefficient <var>c1</var>
and its covariance <var>cov11</var> to compute the fitted function
<var>y</var> and its standard deviation <var>y_err</var> for the model Y =
c_1 X 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>
 
     |