X-Git-Url: https://pintos-os.org/cgi-bin/gitweb.cgi?a=blobdiff_plain;f=lib%2Flinreg%2Flinreg.c;h=17833e57db0dd0bd877a395de00bfd5d8d6a3c0d;hb=000e3e8c5818476c3afbc75fad9347aefb6e902a;hp=2a8aff3f40262ed441499988a6866f87b4aecc24;hpb=d7b5d9144738a5a8989d45a01f4e458a78b68c0b;p=pspp diff --git a/lib/linreg/linreg.c b/lib/linreg/linreg.c index 2a8aff3f40..17833e57db 100644 --- a/lib/linreg/linreg.c +++ b/lib/linreg/linreg.c @@ -90,7 +90,6 @@ pspp_linreg_cache_alloc (size_t n, size_t p) pspp_linreg_cache *c; c = (pspp_linreg_cache *) malloc (sizeof (pspp_linreg_cache)); - c->param_estimates = gsl_vector_alloc (p + 1); c->indep_means = gsl_vector_alloc (p); c->indep_std = gsl_vector_alloc (p); c->ssx = gsl_vector_alloc (p); /* Sums of squares for the independent @@ -113,12 +112,13 @@ pspp_linreg_cache_alloc (size_t n, size_t p) void pspp_linreg_cache_free (pspp_linreg_cache * c) { - gsl_vector_free (c->param_estimates); + int i; + gsl_vector_free (c->indep_means); gsl_vector_free (c->indep_std); gsl_vector_free (c->ss_indeps); gsl_matrix_free (c->cov); - free (c->coeff); + pspp_linreg_coeff_free (c->coeff); free (c); } @@ -139,6 +139,7 @@ pspp_linreg (const gsl_vector * Y, const gsl_matrix * X, gsl_vector_view xty; gsl_vector_view xi; gsl_vector_view xj; + gsl_vector *param_estimates; size_t i; size_t j; @@ -249,7 +250,6 @@ pspp_linreg (const gsl_vector * Y, const gsl_matrix * X, { tmp = gsl_matrix_get (sw, i, cache->n_indeps); cache->coeff[i + 1].estimate = tmp; - gsl_vector_set (cache->param_estimates, i + 1, tmp); m -= tmp * gsl_vector_get (cache->indep_means, i); } /* @@ -285,7 +285,7 @@ pspp_linreg (const gsl_vector * Y, const gsl_matrix * X, } gsl_matrix_set (cache->cov, 0, 0, tmp); - gsl_vector_set (cache->param_estimates, 0, m); + cache->coeff[0].estimate = m; } else { @@ -300,6 +300,8 @@ pspp_linreg (const gsl_vector * Y, const gsl_matrix * X, /* Use QR decomposition via GSL. */ + + param_estimates = gsl_vector_alloc (1 + X->size2); design = gsl_matrix_alloc (X->size1, 1 + X->size2); for (j = 0; j < X->size1; j++) @@ -313,11 +315,16 @@ pspp_linreg (const gsl_vector * Y, const gsl_matrix * X, } gsl_multifit_linear_workspace *wk = gsl_multifit_linear_alloc (design->size1, design->size2); - rc = gsl_multifit_linear (design, Y, cache->param_estimates, + rc = gsl_multifit_linear (design, Y, param_estimates, cache->cov, &(cache->sse), wk); + for (i = 0; i < cache->n_coeffs; i++) + { + cache->coeff[i].estimate = gsl_vector_get (param_estimates, i); + } if (rc == GSL_SUCCESS) { gsl_multifit_linear_free (wk); + gsl_vector_free (param_estimates); } else {