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