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