8028697c6b7f3410b12927a0fe7ae1d82c5a6b4b
[pspp-builds.git] / src / math / linreg.c
1 /* PSPP - a program for statistical analysis.
2    Copyright (C) 2005 Free Software Foundation, Inc. 
3
4    This program is free software: you can redistribute it and/or modify
5    it under the terms of the GNU General Public License as published by
6    the Free Software Foundation, either version 3 of the License, or
7    (at your option) any later version.
8
9    This program is distributed in the hope that it will be useful,
10    but WITHOUT ANY WARRANTY; without even the implied warranty of
11    MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
12    GNU General Public License for more details.
13
14    You should have received a copy of the GNU General Public License
15    along with this program.  If not, see <http://www.gnu.org/licenses/>. */
16
17 #include <config.h>
18 #include <gsl/gsl_blas.h>
19 #include <gsl/gsl_cblas.h>
20 #include <gsl/gsl_errno.h>
21 #include <gsl/gsl_fit.h>
22 #include <gsl/gsl_multifit.h>
23 #include <linreg/sweep.h>
24 #include <math/coefficient.h>
25 #include <math/linreg.h>
26 #include <math/coefficient.h>
27 #include <math/design-matrix.h>
28 #include <src/data/category.h>
29 #include <src/data/variable.h>
30 #include <src/data/value.h>
31 #include <gl/xalloc.h>
32
33 /*
34   Find the least-squares estimate of b for the linear model:
35
36   Y = Xb + Z
37
38   where Y is an n-by-1 column vector, X is an n-by-p matrix of
39   independent variables, b is a p-by-1 vector of regression coefficients,
40   and Z is an n-by-1 normally-distributed random vector with independent
41   identically distributed components with mean 0.
42
43   This estimate is found via the sweep operator or singular-value
44   decomposition with gsl.
45
46
47   References:
48
49   1. Matrix Computations, third edition. GH Golub and CF Van Loan.
50   The Johns Hopkins University Press. 1996. ISBN 0-8018-5414-8.
51
52   2. Numerical Analysis for Statisticians. K Lange. Springer. 1999.
53   ISBN 0-387-94979-8.
54
55   3. Numerical Linear Algebra for Applications in Statistics. JE Gentle.
56   Springer. 1998. ISBN 0-387-98542-5.
57 */
58
59
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_, const struct variable **v)
99 {
100   const pspp_linreg_cache *c = c_;
101   const struct variable *tmp;
102   int i;
103   int j;
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   for (j = 0; j < c->n_coeffs; j++)
115     {
116       tmp = pspp_coeff_get_var (c->coeff[j], 0);
117       assert (tmp != NULL);
118       /* Repeated variables are likely to bunch together, at the end
119          of the array. */
120       i = result - 1;
121       while (i >= 0 && v[i] != tmp)
122         {
123           i--;
124         }
125       if (i < 0 && result < c->n_coeffs)
126         {
127           v[result] = tmp;
128           result++;
129         }
130     }
131   return result;
132 }
133
134 /*
135   Allocate a pspp_linreg_cache and return a pointer
136   to it. n is the number of cases, p is the number of
137   independent variables.
138  */
139 pspp_linreg_cache *
140 pspp_linreg_cache_alloc (const struct variable *depvar, const struct variable **indep_vars,
141                          size_t n, size_t p)
142 {
143   size_t i;
144   pspp_linreg_cache *c;
145
146   c = (pspp_linreg_cache *) malloc (sizeof (pspp_linreg_cache));
147   c->depvar = depvar;
148   c->indep_vars = indep_vars;
149   c->indep_means = gsl_vector_alloc (p);
150   c->indep_std = gsl_vector_alloc (p);
151   c->ssx = gsl_vector_alloc (p);        /* Sums of squares for the
152                                            independent variables.
153                                          */
154   c->ss_indeps = gsl_vector_alloc (p);  /* Sums of squares for the
155                                            model parameters.
156                                          */
157   c->n_obs = n;
158   c->n_indeps = p;
159   c->n_coeffs = 0;
160   for (i = 0; i < p; i++)
161     {
162       if (var_is_numeric (indep_vars[i]))
163         {
164           c->n_coeffs++;
165         }
166       else
167         {
168           c->n_coeffs += cat_get_n_categories (indep_vars[i]) - 1;
169         }
170     }
171
172   c->cov = gsl_matrix_alloc (c->n_coeffs + 1, c->n_coeffs + 1);
173   /*
174      Default settings.
175    */
176   c->method = PSPP_LINREG_SWEEP;
177   c->predict = pspp_linreg_predict;
178   c->residual = pspp_linreg_residual;   /* The procedure to compute my
179                                            residuals. */
180   c->get_vars = pspp_linreg_get_vars;   /* The procedure that returns
181                                            pointers to model
182                                            variables. */
183   c->resid = NULL;              /* The variable storing my residuals. */
184   c->pred = NULL;               /* The variable storing my predicted values. */
185
186   return c;
187 }
188
189 bool
190 pspp_linreg_cache_free (void *m)
191 {
192   int i;
193
194   pspp_linreg_cache *c = m;
195   if (c != NULL)
196     {
197       gsl_vector_free (c->indep_means);
198       gsl_vector_free (c->indep_std);
199       gsl_vector_free (c->ss_indeps);
200       gsl_matrix_free (c->cov);
201       gsl_vector_free (c->ssx);
202       for (i = 0; i < c->n_coeffs; i++)
203         {
204           pspp_coeff_free (c->coeff[i]);
205         }
206       free (c->coeff);
207       free (c);
208     }
209   return true;
210 }
211 static void
212 cache_init (pspp_linreg_cache *cache)
213 {
214   assert (cache != NULL);
215   cache->dft = cache->n_obs - 1;
216   cache->dfm = cache->n_indeps;
217   cache->dfe = cache->dft - cache->dfm;
218   cache->intercept = 0.0;
219 }
220
221 static void
222 post_sweep_computations (pspp_linreg_cache *cache, const struct design_matrix *dm,
223                          gsl_matrix *sw)
224 {
225   gsl_matrix *xm;
226   gsl_matrix_view xtx;
227   gsl_matrix_view xmxtx;
228   double m;
229   double tmp;
230   size_t i;
231   size_t j;
232   int rc;
233   
234   assert (sw != NULL);
235   assert (cache != NULL);
236
237   cache->sse = gsl_matrix_get (sw, cache->n_indeps, cache->n_indeps);
238   cache->mse = cache->sse / cache->dfe;
239   /*
240     Get the intercept.
241   */
242   m = cache->depvar_mean;
243   for (i = 0; i < cache->n_indeps; i++)
244     {
245       tmp = gsl_matrix_get (sw, i, cache->n_indeps);
246       cache->coeff[i]->estimate = tmp;
247       m -= tmp * pspp_linreg_get_indep_variable_mean (cache, design_matrix_col_to_var (dm, i));
248     }
249   /*
250     Get the covariance matrix of the parameter estimates.
251     Only the upper triangle is necessary.
252   */
253   
254   /*
255     The loops below do not compute the entries related
256     to the estimated intercept.
257   */
258   for (i = 0; i < cache->n_indeps; i++)
259     for (j = i; j < cache->n_indeps; j++)
260       {
261         tmp = -1.0 * cache->mse * gsl_matrix_get (sw, i, j);
262         gsl_matrix_set (cache->cov, i + 1, j + 1, tmp);
263       }
264   /*
265     Get the covariances related to the intercept.
266   */
267   xtx = gsl_matrix_submatrix (sw, 0, 0, cache->n_indeps, cache->n_indeps);
268   xmxtx = gsl_matrix_submatrix (cache->cov, 0, 1, 1, cache->n_indeps);
269   xm = gsl_matrix_calloc (1, cache->n_indeps);
270   for (i = 0; i < xm->size2; i++)
271     {
272       gsl_matrix_set (xm, 0, i, 
273                       pspp_linreg_get_indep_variable_mean (cache, design_matrix_col_to_var (dm, i)));
274     }
275   rc = gsl_blas_dsymm (CblasRight, CblasUpper, cache->mse,
276                        &xtx.matrix, xm, 0.0, &xmxtx.matrix);
277   gsl_matrix_free (xm);
278   if (rc == GSL_SUCCESS)
279     {
280       tmp = cache->mse / cache->n_obs;
281       for (i = 1; i < 1 + cache->n_indeps; i++)
282         {
283           tmp -= gsl_matrix_get (cache->cov, 0, i)
284             * pspp_linreg_get_indep_variable_mean (cache, design_matrix_col_to_var (dm, i - 1));
285         }
286       gsl_matrix_set (cache->cov, 0, 0, tmp);
287       
288       cache->intercept = m;
289     }
290   else
291     {
292       fprintf (stderr, "%s:%d:gsl_blas_dsymm: %s\n",
293                __FILE__, __LINE__, gsl_strerror (rc));
294       exit (rc);
295     }
296 }  
297   
298 /*
299   Fit the linear model via least squares. All pointers passed to pspp_linreg
300   are assumed to be allocated to the correct size and initialized to the
301   values as indicated by opts.
302  */
303 int
304 pspp_linreg (const gsl_vector * Y, const struct design_matrix *dm,
305              const pspp_linreg_opts * opts, pspp_linreg_cache * cache)
306 {
307   int rc;
308   gsl_matrix *design = NULL;
309   gsl_matrix_view xtx;
310   gsl_vector_view xty;
311   gsl_vector_view xi;
312   gsl_vector_view xj;
313   gsl_vector *param_estimates;
314   struct pspp_coeff *coef;
315   const struct variable *v;
316   const union value *val;
317
318   size_t i;
319   size_t j;
320   double tmp;
321   double m;
322   double s;
323   double ss;
324
325   if (cache == NULL)
326     {
327       return GSL_EFAULT;
328     }
329   if (opts->get_depvar_mean_std)
330     {
331       linreg_mean_std (gsl_vector_const_subvector (Y, 0, Y->size),
332                        &m, &s, &ss);
333       cache->depvar_mean = m;
334       cache->depvar_std = s;
335       cache->sst = ss;
336     }
337   cache_init (cache);
338   cache->n_coeffs = dm->m->size2;
339   for (i = 0; i < dm->m->size2; i++)
340     {
341       if (opts->get_indep_mean_std[i])
342         {
343           linreg_mean_std (gsl_matrix_const_column (dm->m, i), &m, &s, &ss);
344           v = design_matrix_col_to_var (dm, i);
345           val = NULL;
346           if (var_is_alpha (v))
347             {
348               j = i - design_matrix_var_to_column (dm, v);
349               val = cat_subscript_to_value (j, v);
350             }
351           coef = pspp_linreg_get_coeff (cache, v, val);
352           pspp_coeff_set_mean (coef, m);
353           pspp_coeff_set_sd (coef, s);
354           gsl_vector_set (cache->ssx, i, ss);
355
356         }
357     }
358
359   if (cache->method == PSPP_LINREG_SWEEP)
360     {
361       gsl_matrix *sw;
362       /*
363          Subtract the means to improve the condition of the design
364          matrix. This requires copying dm->m and Y. We do not divide by the
365          standard deviations of the independent variables here since doing
366          so would cause a miscalculation of the residual sums of
367          squares. Dividing by the standard deviation is done GSL's linear
368          regression functions, so if the design matrix has a poor
369          condition, use QR decomposition.
370
371          The design matrix here does not include a column for the intercept
372          (i.e., a column of 1's). If using PSPP_LINREG_QR, we need that column,
373          so design is allocated here when sweeping, or below if using QR.
374        */
375       design = gsl_matrix_alloc (dm->m->size1, dm->m->size2);
376       for (i = 0; i < dm->m->size2; i++)
377         {
378           v = design_matrix_col_to_var (dm, i);
379           m = pspp_linreg_get_indep_variable_mean (cache, v);
380           for (j = 0; j < dm->m->size1; j++)
381             {
382               tmp = (gsl_matrix_get (dm->m, j, i) - m);
383               gsl_matrix_set (design, j, i, tmp);
384             }
385         }
386       sw = gsl_matrix_calloc (cache->n_indeps + 1, cache->n_indeps + 1);
387       xtx = gsl_matrix_submatrix (sw, 0, 0, cache->n_indeps, cache->n_indeps);
388
389       for (i = 0; i < xtx.matrix.size1; i++)
390         {
391           tmp = gsl_vector_get (cache->ssx, i);
392           gsl_matrix_set (&(xtx.matrix), i, i, tmp);
393           xi = gsl_matrix_column (design, i);
394           for (j = (i + 1); j < xtx.matrix.size2; j++)
395             {
396               xj = gsl_matrix_column (design, j);
397               gsl_blas_ddot (&(xi.vector), &(xj.vector), &tmp);
398               gsl_matrix_set (&(xtx.matrix), i, j, tmp);
399             }
400         }
401
402       gsl_matrix_set (sw, cache->n_indeps, cache->n_indeps, cache->sst);
403       xty = gsl_matrix_column (sw, cache->n_indeps);
404       /*
405          This loop starts at 1, with i=0 outside the loop, so we can get
406          the model sum of squares due to the first independent variable.
407        */
408       xi = gsl_matrix_column (design, 0);
409       gsl_blas_ddot (&(xi.vector), Y, &tmp);
410       gsl_vector_set (&(xty.vector), 0, tmp);
411       tmp *= tmp / gsl_vector_get (cache->ssx, 0);
412       gsl_vector_set (cache->ss_indeps, 0, tmp);
413       for (i = 1; i < cache->n_indeps; i++)
414         {
415           xi = gsl_matrix_column (design, i);
416           gsl_blas_ddot (&(xi.vector), Y, &tmp);
417           gsl_vector_set (&(xty.vector), i, tmp);
418         }
419
420       /*
421          Sweep on the matrix sw, which contains XtX, XtY and YtY.
422        */
423       reg_sweep (sw);
424       post_sweep_computations (cache, dm, sw);
425       gsl_matrix_free (sw);
426     }
427   else if (cache->method == PSPP_LINREG_CONDITIONAL_INVERSE)
428     {
429       /*
430         Use the SVD of X^T X to find a conditional inverse of X^TX. If
431         the SVD is X^T X = U D V^T, then set the conditional inverse
432         to (X^T X)^c = V D^- U^T. D^- is defined as follows: If entry
433         (i, i) has value sigma_i, then entry (i, i) of D^- is 1 /
434         sigma_i if sigma_i > 0, and 0 otherwise. Then solve the normal
435         equations by setting the estimated parameter vector to 
436         (X^TX)^c X^T Y.
437        */
438     }
439   else
440     {
441       gsl_multifit_linear_workspace *wk;
442       /*
443          Use QR decomposition via GSL.
444        */
445
446       param_estimates = gsl_vector_alloc (1 + dm->m->size2);
447       design = gsl_matrix_alloc (dm->m->size1, 1 + dm->m->size2);
448
449       for (j = 0; j < dm->m->size1; j++)
450         {
451           gsl_matrix_set (design, j, 0, 1.0);
452           for (i = 0; i < dm->m->size2; i++)
453             {
454               tmp = gsl_matrix_get (dm->m, j, i);
455               gsl_matrix_set (design, j, i + 1, tmp);
456             }
457         }
458
459       wk = gsl_multifit_linear_alloc (design->size1, design->size2);
460       rc = gsl_multifit_linear (design, Y, param_estimates,
461                                 cache->cov, &(cache->sse), wk);
462       for (i = 0; i < cache->n_coeffs; i++)
463         {
464           cache->coeff[i]->estimate = gsl_vector_get (param_estimates, i + 1);
465         }
466       cache->intercept = gsl_vector_get (param_estimates, 0);
467       if (rc == GSL_SUCCESS)
468         {
469           gsl_multifit_linear_free (wk);
470           gsl_vector_free (param_estimates);
471         }
472       else
473         {
474           fprintf (stderr, "%s:%d: gsl_multifit_linear returned %d\n",
475                    __FILE__, __LINE__, rc);
476         }
477     }
478
479
480   cache->ssm = cache->sst - cache->sse;
481   /*
482      Get the remaining sums of squares for the independent
483      variables.
484    */
485   m = 0;
486   for (i = 1; i < cache->n_indeps; i++)
487     {
488       j = i - 1;
489       m += gsl_vector_get (cache->ss_indeps, j);
490       tmp = cache->ssm - m;
491       gsl_vector_set (cache->ss_indeps, i, tmp);
492     }
493
494   gsl_matrix_free (design);
495   return GSL_SUCCESS;
496 }
497
498 /*
499   Is the coefficient COEF contained in the list of coefficients
500   COEF_LIST?
501  */
502 static int
503 has_coefficient (const struct pspp_coeff **coef_list, const struct pspp_coeff *coef,
504                  size_t n)
505 {
506   size_t i = 0;
507
508   while (i < n)
509     {
510       if (coef_list[i] == coef)
511         {
512           return 1;
513         }
514       i++;
515     }
516   return 0;
517 }
518 /*
519   Predict the value of the dependent variable with the
520   new set of predictors. PREDICTORS must point to a list
521   of variables, each of whose values are stored in VALS,
522   in the same order.
523  */
524 double
525 pspp_linreg_predict (const struct variable **predictors,
526                      const union value **vals, const void *c_, int n_vals)
527 {
528   const pspp_linreg_cache *c = c_;
529   int j;
530   size_t next_coef = 0;
531   const struct pspp_coeff **coef_list;
532   const struct pspp_coeff *coe;
533   double result;
534   double tmp;
535
536   if (predictors == NULL || vals == NULL || c == NULL)
537     {
538       return GSL_NAN;
539     }
540   if (c->coeff == NULL)
541     {
542       /* The stupid model: just guess the mean. */
543       return c->depvar_mean;
544     }
545   coef_list = xnmalloc (c->n_coeffs, sizeof (*coef_list));
546   result = c->intercept;
547
548   /*
549      The loops guard against the possibility that the caller passed us
550      inadequate information, such as too few or too many values, or
551      a redundant list of variable names.
552    */
553   for (j = 0; j < n_vals; j++)
554     {
555       coe = pspp_linreg_get_coeff (c, predictors[j], vals[j]);
556       if (!has_coefficient (coef_list, coe, next_coef))
557         {
558           tmp = pspp_coeff_get_est (coe);
559           if (var_is_numeric (predictors[j]))
560             {
561               tmp *= vals[j]->f;
562             }
563           result += tmp;
564           coef_list[next_coef++] = coe;
565         }
566     }
567   free (coef_list);
568
569   return result;
570 }
571
572 double
573 pspp_linreg_residual (const struct variable **predictors,
574                       const union value **vals,
575                       const union value *obs, const void *c, int n_vals)
576 {
577   double pred;
578   double result;
579
580   if (predictors == NULL || vals == NULL || c == NULL || obs == NULL)
581     {
582       return GSL_NAN;
583     }
584   pred = pspp_linreg_predict (predictors, vals, c, n_vals);
585
586   result = isnan (pred) ? GSL_NAN : (obs->f - pred);
587   return result;
588 }
589
590 /*
591   Which coefficient is associated with V? The VAL argument is relevant
592   only to categorical variables.
593  */
594 struct pspp_coeff *
595 pspp_linreg_get_coeff (const pspp_linreg_cache * c,
596                        const struct variable *v, const union value *val)
597 {
598   if (c == NULL)
599     {
600       return NULL;
601     }
602   if (c->coeff == NULL || c->n_indeps == 0 || v == NULL)
603     {
604       return NULL;
605     }
606   return pspp_coeff_var_to_coeff (v, c->coeff, c->n_coeffs, val);
607 }
608 /*
609   Return the standard deviation of the independent variable.
610  */
611 double pspp_linreg_get_indep_variable_sd (pspp_linreg_cache *c, const struct variable *v)
612 {
613   if (var_is_numeric (v))
614     {
615       const struct pspp_coeff *coef;
616       coef = pspp_linreg_get_coeff (c, v, NULL);
617       return pspp_coeff_get_sd (coef);
618     }
619   return GSL_NAN;
620 }
621
622 void pspp_linreg_set_indep_variable_sd (pspp_linreg_cache *c, const struct variable *v, 
623                                         double s)
624 {
625   if (var_is_numeric (v))
626     {
627       struct pspp_coeff *coef;
628       coef = pspp_linreg_get_coeff (c, v, NULL);
629       pspp_coeff_set_sd (coef, s);
630     }
631 }
632
633 /*
634   Mean of the independent variable.
635  */
636 double pspp_linreg_get_indep_variable_mean (pspp_linreg_cache *c, const struct variable *v)
637 {
638   if (var_is_numeric (v))
639     {
640       struct pspp_coeff *coef;
641       coef = pspp_linreg_get_coeff (c, v, NULL);
642       return pspp_coeff_get_mean (coef);
643     }
644   return GSL_NAN;
645 }
646
647 void pspp_linreg_set_indep_variable_mean (pspp_linreg_cache *c, const struct variable *v, 
648                                           double m)
649 {
650   if (var_is_numeric (v))
651     {
652       struct pspp_coeff *coef;
653       coef = pspp_linreg_get_coeff (c, v, NULL);
654       pspp_coeff_set_mean (coef, m);
655     }
656 }
657
658 /*
659   Make sure the dependent variable is at the last column, and that
660   only variables in the model are in the covariance matrix. 
661  */
662 static struct design_matrix *
663 rearrange_covariance_matrix (const struct covariance_matrix *cm, pspp_linreg_cache *c)
664 {
665   const struct variable **model_vars;
666   struct design_matrix *cov;
667   struct design_matrix *result;
668   size_t *permutation;
669   size_t i;
670   size_t j;
671   size_t k;
672   size_t n_coeffs = 0;
673
674   assert (cm != NULL);
675   cov = covariance_to_design (cm);
676   assert (cov != NULL);
677   assert (c != NULL);
678   assert (cov->m->size1 > 0);
679   assert (cov->m->size2 == cov->m->size1);
680   model_vars = xnmalloc (1 + c->n_indeps, sizeof (*model_vars));
681
682   /*
683     Put the model variables in the right order in MODEL_VARS.
684     Count the number of coefficients.
685    */
686   for (i = 0; i < c->n_indeps; i++)
687     {
688       model_vars[i] = c->indep_vars[i];
689     }
690   model_vars[i] = c->depvar;
691   result = covariance_matrix_create (1 + c->n_indeps, model_vars);
692   permutation = xnmalloc (design_matrix_get_n_cols (result), sizeof (*permutation));
693
694   for (j = 0; j < cov->m->size2; j++)
695     {
696       k = 0;
697       while (k < result->m->size2)
698         {
699           if (design_matrix_col_to_var (cov, j) == design_matrix_col_to_var (result, k)) 
700             {
701               permutation[k] = j;
702             }
703           k++;
704         }
705     }
706   for (i = 0; i < result->m->size1; i++)
707     for (j = 0; j < result->m->size2; j++)
708       {
709         gsl_matrix_set (result->m, i, j, gsl_matrix_get (cov->m, permutation[i], permutation[j]));
710       }
711   free (permutation);
712   free (model_vars);
713   return result;
714 }
715 /*
716   Estimate the model parameters from the covariance matrix only. This
717   method uses less memory than PSPP_LINREG, which requires the entire
718   data set to be stored in memory.
719
720   The function assumes FULL_COV may contain columns corresponding to
721   variables that are not in the model. It fixes this in
722   REARRANG_COVARIANCE_MATRIX. This allows the caller to compute a
723   large covariance matrix once before, then pass it to this without
724   having to alter it. The problem is that this means the caller must
725   set CACHE->N_COEFFS.
726 */
727 void
728 pspp_linreg_with_cov (const struct covariance_matrix *full_cov, 
729                       pspp_linreg_cache * cache)
730 {
731   struct design_matrix *cov;
732
733   assert (full_cov != NULL);
734   assert (cache != NULL);
735
736   cov = rearrange_covariance_matrix (full_cov, cache);
737   cache_init (cache);
738   reg_sweep (cov->m);
739   post_sweep_computations (cache, cov, cov->m);  
740   design_matrix_destroy (cov);
741 }
742
743 double pspp_linreg_mse (const pspp_linreg_cache *c)
744 {
745   assert (c != NULL);
746   return (c->sse / c->dfe);
747 }