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