1fe4a79e5784aff1868c09e6c3d1aa56463508df
[pspp-builds.git] / lib / linreg / linreg.c
1 /* lib/linreg/linreg.c
2
3  Copyright (C) 2005 Free Software Foundation, Inc.
4  Written by Jason H. Stover.
5
6  This program is free software; you can redistribute it and/or modify
7  it under the terms of the GNU General Public License as published by
8  the Free Software Foundation; either version 2 of the License, or (at
9  your option) any later version.
10
11  This program is distributed in the hope that it will be useful, but
12  WITHOUT ANY WARRANTY; without even the implied warranty of
13  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
14  General Public License for more details.
15
16  You should have received a copy of the GNU General Public License
17  along with this program; if not, write to the Free Software
18  Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
19  02111-1307, USA.
20 */
21
22 /*
23   Find the least-squares estimate of b for the linear model:
24
25   Y = Xb + Z
26
27   where Y is an n-by-1 column vector, X is an n-by-p matrix of 
28   independent variables, b is a p-by-1 vector of regression coefficients,
29   and Z is an n-by-1 normally-distributed random vector with independent
30   identically distributed components with mean 0.
31
32   This estimate is found via the sweep operator or singular-value
33   decomposition with gsl.
34
35
36   References:
37
38   1. Matrix Computations, third edition. GH Golub and CF Van Loan.
39   The Johns Hopkins University Press. 1996. ISBN 0-8018-5414-8.
40
41   2. Numerical Analysis for Statisticians. K Lange. Springer. 1999.
42   ISBN 0-387-94979-8.
43
44   3. Numerical Linear Algebra for Applications in Statistics. JE Gentle.
45   Springer. 1998. ISBN 0-387-98542-5.
46 */
47
48 #include "pspp_linreg.h"
49 #include <gsl/gsl_errno.h>
50 /*
51   Get the mean and standard deviation of a vector
52   of doubles via a form of the Kalman filter as
53   described on page 32 of [3].
54  */
55 static int
56 linreg_mean_std (gsl_vector_const_view v, double *mp, double *sp, double *ssp)
57 {
58   size_t i;
59   double j = 0.0;
60   double d;
61   double tmp;
62   double mean;
63   double variance;
64
65   mean = gsl_vector_get (&v.vector, 0);
66   variance = 0;
67   for (i = 1; i < v.vector.size; i++)
68     {
69       j = (double) i + 1.0;
70       tmp = gsl_vector_get (&v.vector, i);
71       d = (tmp - mean) / j;
72       mean += d;
73       variance += j * (j - 1.0) * d * d;
74     }
75   *mp = mean;
76   *sp = sqrt (variance / (j - 1.0));
77   *ssp = variance;
78
79   return GSL_SUCCESS;
80 }
81
82 /*
83   Allocate a pspp_linreg_cache and return a pointer
84   to it. n is the number of cases, p is the number of 
85   independent variables.
86  */
87 pspp_linreg_cache *
88 pspp_linreg_cache_alloc (size_t n, size_t p)
89 {
90   pspp_linreg_cache *c;
91
92   c = (pspp_linreg_cache *) malloc (sizeof (pspp_linreg_cache));
93   c->indep_means = gsl_vector_alloc (p);
94   c->indep_std = gsl_vector_alloc (p);
95   c->ssx = gsl_vector_alloc (p);        /* Sums of squares for the independent
96                                            variables.
97                                          */
98   c->ss_indeps = gsl_vector_alloc (p);  /* Sums of squares for the model 
99                                            parameters. 
100                                          */
101   c->cov = gsl_matrix_alloc (p + 1, p + 1);     /* Covariance matrix. */
102   c->n_obs = n;
103   c->n_indeps = p;
104   /*
105      Default settings.
106    */
107   c->method = PSPP_LINREG_SWEEP;
108
109   return c;
110 }
111
112 void
113 pspp_linreg_cache_free (pspp_linreg_cache * c)
114 {
115   int i;
116
117   gsl_vector_free (c->indep_means);
118   gsl_vector_free (c->indep_std);
119   gsl_vector_free (c->ss_indeps);
120   gsl_matrix_free (c->cov);
121 #if 0  
122   for (i = 1; i < c->n_coeffs; i++)
123     {
124       pspp_linreg_coeff_free (c->coeff + i);
125     }
126 #endif
127   free (c);
128 }
129
130 /*
131   Fit the linear model via least squares. All pointers passed to pspp_linreg
132   are assumed to be allocated to the correct size and initialized to the
133   values as indicated by opts. 
134  */
135 int
136 pspp_linreg (const gsl_vector * Y, const gsl_matrix * X,
137              const pspp_linreg_opts * opts, pspp_linreg_cache * cache)
138 {
139   int rc;
140   gsl_matrix *design;
141   gsl_matrix_view xtx;
142   gsl_matrix_view xm;
143   gsl_matrix_view xmxtx;
144   gsl_vector_view xty;
145   gsl_vector_view xi;
146   gsl_vector_view xj;
147   gsl_vector *param_estimates;
148
149   size_t i;
150   size_t j;
151   double tmp;
152   double m;
153   double s;
154   double ss;
155
156   if (cache == NULL)
157     {
158       return GSL_EFAULT;
159     }
160   if (opts->get_depvar_mean_std)
161     {
162       linreg_mean_std (gsl_vector_const_subvector (Y, 0, Y->size),
163                        &m, &s, &ss);
164       cache->depvar_mean = m;
165       cache->depvar_std = s;
166       cache->sst = ss;
167     }
168   for (i = 0; i < cache->n_indeps; i++)
169     {
170       if (opts->get_indep_mean_std[i])
171         {
172           linreg_mean_std (gsl_matrix_const_column (X, i), &m, &s, &ss);
173           gsl_vector_set (cache->indep_means, i, m);
174           gsl_vector_set (cache->indep_std, i, s);
175           gsl_vector_set (cache->ssx, i, ss);
176         }
177     }
178   cache->dft = cache->n_obs - 1;
179   cache->dfm = cache->n_indeps;
180   cache->dfe = cache->dft - cache->dfm;
181   cache->n_coeffs = X->size2 + 1; /* Adjust this later to allow for regression
182                                      through the origin.
183                                   */
184   if (cache->method == PSPP_LINREG_SWEEP)
185     {
186       gsl_matrix *sw;
187       /*
188          Subtract the means to improve the condition of the design
189          matrix. This requires copying X and Y. We do not divide by the
190          standard deviations of the independent variables here since doing
191          so would cause a miscalculation of the residual sums of
192          squares. Dividing by the standard deviation is done GSL's linear
193          regression functions, so if the design matrix has a poor
194          condition, use QR decomposition.
195
196          The design matrix here does not include a column for the intercept
197          (i.e., a column of 1's). If using PSPP_LINREG_QR, we need that column,
198          so design is allocated here when sweeping, or below if using QR.
199        */
200       design = gsl_matrix_alloc (X->size1, X->size2);
201       for (i = 0; i < X->size2; i++)
202         {
203           m = gsl_vector_get (cache->indep_means, i);
204           for (j = 0; j < X->size1; j++)
205             {
206               tmp = (gsl_matrix_get (X, j, i) - m);
207               gsl_matrix_set (design, j, i, tmp);
208             }
209         }
210       sw = gsl_matrix_calloc (cache->n_indeps + 1, cache->n_indeps + 1);
211       xtx = gsl_matrix_submatrix (sw, 0, 0, cache->n_indeps, cache->n_indeps);
212
213       for (i = 0; i < xtx.matrix.size1; i++)
214         {
215           tmp = gsl_vector_get (cache->ssx, i);
216           gsl_matrix_set (&(xtx.matrix), i, i, tmp);
217           xi = gsl_matrix_column (design, i);
218           for (j = (i + 1); j < xtx.matrix.size2; j++)
219             {
220               xj = gsl_matrix_column (design, j);
221               gsl_blas_ddot (&(xi.vector), &(xj.vector), &tmp);
222               gsl_matrix_set (&(xtx.matrix), i, j, tmp);
223             }
224         }
225
226       gsl_matrix_set (sw, cache->n_indeps, cache->n_indeps, cache->sst);
227       xty = gsl_matrix_column (sw, cache->n_indeps);
228       /*
229          This loop starts at 1, with i=0 outside the loop, so we can get
230          the model sum of squares due to the first independent variable.
231        */
232       xi = gsl_matrix_column (design, 0);
233       gsl_blas_ddot (&(xi.vector), Y, &tmp);
234       gsl_vector_set (&(xty.vector), 0, tmp);
235       tmp *= tmp / gsl_vector_get (cache->ssx, 0);
236       gsl_vector_set (cache->ss_indeps, 0, tmp);
237       for (i = 1; i < cache->n_indeps; i++)
238         {
239           xi = gsl_matrix_column (design, i);
240           gsl_blas_ddot (&(xi.vector), Y, &tmp);
241           gsl_vector_set (&(xty.vector), i, tmp);
242         }
243
244       /*
245          Sweep on the matrix sw, which contains XtX, XtY and YtY.
246        */
247       pspp_reg_sweep (sw);
248       cache->sse = gsl_matrix_get (sw, cache->n_indeps, cache->n_indeps);
249       cache->mse = cache->sse / cache->dfe;
250       /*
251          Get the intercept.
252        */
253       m = cache->depvar_mean;
254       for (i = 0; i < cache->n_indeps; i++)
255         {
256           tmp = gsl_matrix_get (sw, i, cache->n_indeps);
257           cache->coeff[i + 1].estimate = tmp;
258           m -= tmp * gsl_vector_get (cache->indep_means, i);
259         }
260       /*
261          Get the covariance matrix of the parameter estimates.
262          Only the upper triangle is necessary. 
263        */
264
265       /*
266          The loops below do not compute the entries related
267          to the estimated intercept.
268        */
269       for (i = 0; i < cache->n_indeps; i++)
270         for (j = i; j < cache->n_indeps; j++)
271           {
272             tmp = -1.0 * cache->mse * gsl_matrix_get (sw, i, j);
273             gsl_matrix_set (cache->cov, i + 1, j + 1, tmp);
274           }
275       /*
276          Get the covariances related to the intercept.
277        */
278       xtx = gsl_matrix_submatrix (sw, 0, 0, cache->n_indeps, cache->n_indeps);
279       xmxtx = gsl_matrix_submatrix (cache->cov, 0, 1, 1, cache->n_indeps);
280       xm = gsl_matrix_view_vector (cache->indep_means, 1, cache->n_indeps);
281       rc = gsl_blas_dsymm (CblasRight, CblasUpper, cache->mse,
282                            &xtx.matrix, &xm.matrix, 0.0, &xmxtx.matrix);
283       if (rc == GSL_SUCCESS)
284         {
285           tmp = cache->mse / cache->n_obs;
286           for (i = 1; i < 1 + cache->n_indeps; i++)
287             {
288               tmp -= gsl_matrix_get (cache->cov, 0, i)
289                 * gsl_vector_get (cache->indep_means, i - 1);
290             }
291           gsl_matrix_set (cache->cov, 0, 0, tmp);
292
293           cache->coeff[0].estimate = m;
294         }
295       else
296         {
297           fprintf (stderr, "%s:%d:gsl_blas_dsymm: %s\n",
298                    __FILE__, __LINE__, gsl_strerror (rc));
299           exit (rc);
300         }
301       gsl_matrix_free (sw);
302     }
303   else
304     {
305       /*
306          Use QR decomposition via GSL.
307        */
308
309       param_estimates = gsl_vector_alloc (1 + X->size2);
310       design = gsl_matrix_alloc (X->size1, 1 + X->size2);
311
312       for (j = 0; j < X->size1; j++)
313         {
314           gsl_matrix_set (design, j, 0, 1.0);
315           for (i = 0; i < X->size2; i++)
316             {
317               tmp = gsl_matrix_get (X, j, i);
318               gsl_matrix_set (design, j, i + 1, tmp);
319             }
320         }
321       gsl_multifit_linear_workspace *wk =
322         gsl_multifit_linear_alloc (design->size1, design->size2);
323       rc = gsl_multifit_linear (design, Y, param_estimates,
324                                 cache->cov, &(cache->sse), wk);
325       for (i = 0; i < cache->n_coeffs; i++)
326         {
327           cache->coeff[i].estimate = gsl_vector_get (param_estimates, i);
328         }
329       if (rc == GSL_SUCCESS)
330         {
331           gsl_multifit_linear_free (wk);
332           gsl_vector_free (param_estimates);
333         }
334       else
335         {
336           fprintf (stderr, "%s:%d: gsl_multifit_linear returned %d\n",
337                    __FILE__, __LINE__, rc);
338         }
339     }
340
341
342   cache->ssm = cache->sst - cache->sse;
343   /*
344      Get the remaining sums of squares for the independent
345      variables.
346    */
347   m = 0;
348   for (i = 1; i < cache->n_indeps; i++)
349     {
350       j = i - 1;
351       m += gsl_vector_get (cache->ss_indeps, j);
352       tmp = cache->ssm - m;
353       gsl_vector_set (cache->ss_indeps, i, tmp);
354     }
355
356   gsl_matrix_free (design);
357   return GSL_SUCCESS;
358 }