|
Packit |
67cb25 |
/* multifit/multiwlinear.c
|
|
Packit |
67cb25 |
*
|
|
Packit |
67cb25 |
* Copyright (C) 2015 Patrick Alken
|
|
Packit |
67cb25 |
*
|
|
Packit |
67cb25 |
* This program is free software; you can redistribute it and/or modify
|
|
Packit |
67cb25 |
* it under the terms of the GNU General Public License as published by
|
|
Packit |
67cb25 |
* the Free Software Foundation; either version 3 of the License, or (at
|
|
Packit |
67cb25 |
* your option) any later version.
|
|
Packit |
67cb25 |
*
|
|
Packit |
67cb25 |
* This program is distributed in the hope that it will be useful, but
|
|
Packit |
67cb25 |
* WITHOUT ANY WARRANTY; without even the implied warranty of
|
|
Packit |
67cb25 |
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
|
|
Packit |
67cb25 |
* General Public License for more details.
|
|
Packit |
67cb25 |
*
|
|
Packit |
67cb25 |
* You should have received a copy of the GNU General Public License
|
|
Packit |
67cb25 |
* along with this program; if not, write to the Free Software
|
|
Packit |
67cb25 |
* Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
|
|
Packit |
67cb25 |
*/
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
#include <config.h>
|
|
Packit |
67cb25 |
#include <gsl/gsl_errno.h>
|
|
Packit |
67cb25 |
#include <gsl/gsl_multifit.h>
|
|
Packit |
67cb25 |
#include <gsl/gsl_blas.h>
|
|
Packit |
67cb25 |
#include <gsl/gsl_vector.h>
|
|
Packit |
67cb25 |
#include <gsl/gsl_matrix.h>
|
|
Packit |
67cb25 |
#include <gsl/gsl_linalg.h>
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
#include "linear_common.c"
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
gsl_multifit_wlinear (const gsl_matrix * X,
|
|
Packit |
67cb25 |
const gsl_vector * w,
|
|
Packit |
67cb25 |
const gsl_vector * y,
|
|
Packit |
67cb25 |
gsl_vector * c,
|
|
Packit |
67cb25 |
gsl_matrix * cov,
|
|
Packit |
67cb25 |
double *chisq, gsl_multifit_linear_workspace * work)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
size_t rank;
|
|
Packit |
67cb25 |
int status = gsl_multifit_wlinear_tsvd(X, w, y, GSL_DBL_EPSILON, c, cov, chisq, &rank, work);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
return status;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
gsl_multifit_wlinear_tsvd (const gsl_matrix * X,
|
|
Packit |
67cb25 |
const gsl_vector * w,
|
|
Packit |
67cb25 |
const gsl_vector * y,
|
|
Packit |
67cb25 |
const double tol,
|
|
Packit |
67cb25 |
gsl_vector * c,
|
|
Packit |
67cb25 |
gsl_matrix * cov,
|
|
Packit |
67cb25 |
double * chisq,
|
|
Packit |
67cb25 |
size_t * rank,
|
|
Packit |
67cb25 |
gsl_multifit_linear_workspace * work)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
const size_t n = X->size1;
|
|
Packit |
67cb25 |
const size_t p = X->size2;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (y->size != n)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
GSL_ERROR("number of observations in y does not match matrix",
|
|
Packit |
67cb25 |
GSL_EBADLEN);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
else if (w->size != n)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
GSL_ERROR("number of weights in w does not match matrix",
|
|
Packit |
67cb25 |
GSL_EBADLEN);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
else if (p != c->size)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
GSL_ERROR ("number of parameters c does not match matrix",
|
|
Packit |
67cb25 |
GSL_EBADLEN);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
else if (tol <= 0)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
GSL_ERROR ("tolerance must be positive", GSL_EINVAL);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
else
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
int status;
|
|
Packit |
67cb25 |
double rnorm, snorm;
|
|
Packit |
67cb25 |
gsl_matrix_view A = gsl_matrix_submatrix(work->A, 0, 0, n, p);
|
|
Packit |
67cb25 |
gsl_vector_view b = gsl_vector_subvector(work->t, 0, n);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* compute A = sqrt(W) X, b = sqrt(W) y */
|
|
Packit |
67cb25 |
status = gsl_multifit_linear_applyW(X, w, y, &A.matrix, &b.vector);
|
|
Packit |
67cb25 |
if (status)
|
|
Packit |
67cb25 |
return status;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* compute SVD of A */
|
|
Packit |
67cb25 |
status = gsl_multifit_linear_bsvd(&A.matrix, work);
|
|
Packit |
67cb25 |
if (status)
|
|
Packit |
67cb25 |
return status;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
status = multifit_linear_solve(X, &b.vector, tol, 0.0, rank,
|
|
Packit |
67cb25 |
c, &rnorm, &snorm, work);
|
|
Packit |
67cb25 |
if (status)
|
|
Packit |
67cb25 |
return status;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
*chisq = rnorm * rnorm;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* variance-covariance matrix cov = s2 * (Q S^-1) (Q S^-1)^T */
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
const size_t p = X->size2;
|
|
Packit |
67cb25 |
size_t i, j;
|
|
Packit |
67cb25 |
gsl_matrix_view QSI = gsl_matrix_submatrix(work->QSI, 0, 0, p, p);
|
|
Packit |
67cb25 |
gsl_vector_view D = gsl_vector_subvector(work->D, 0, p);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
for (i = 0; i < p; i++)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
gsl_vector_view row_i = gsl_matrix_row (&QSI.matrix, i);
|
|
Packit |
67cb25 |
double d_i = gsl_vector_get (&D.vector, i);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
for (j = i; j < p; j++)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
gsl_vector_view row_j = gsl_matrix_row (&QSI.matrix, j);
|
|
Packit |
67cb25 |
double d_j = gsl_vector_get (&D.vector, j);
|
|
Packit |
67cb25 |
double s;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
gsl_blas_ddot (&row_i.vector, &row_j.vector, &s);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
gsl_matrix_set (cov, i, j, s / (d_i * d_j));
|
|
Packit |
67cb25 |
gsl_matrix_set (cov, j, i, s / (d_i * d_j));
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
return GSL_SUCCESS;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
gsl_multifit_wlinear_svd (const gsl_matrix * X,
|
|
Packit |
67cb25 |
const gsl_vector * w,
|
|
Packit |
67cb25 |
const gsl_vector * y,
|
|
Packit |
67cb25 |
double tol,
|
|
Packit |
67cb25 |
size_t * rank,
|
|
Packit |
67cb25 |
gsl_vector * c,
|
|
Packit |
67cb25 |
gsl_matrix * cov,
|
|
Packit |
67cb25 |
double *chisq, gsl_multifit_linear_workspace * work)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
int status = gsl_multifit_wlinear_tsvd(X, w, y, tol, c, cov, chisq, rank, work);
|
|
Packit |
67cb25 |
return status;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
gsl_multifit_wlinear_usvd (const gsl_matrix * X,
|
|
Packit |
67cb25 |
const gsl_vector * w,
|
|
Packit |
67cb25 |
const gsl_vector * y,
|
|
Packit |
67cb25 |
double tol,
|
|
Packit |
67cb25 |
size_t * rank,
|
|
Packit |
67cb25 |
gsl_vector * c,
|
|
Packit |
67cb25 |
gsl_matrix * cov,
|
|
Packit |
67cb25 |
double *chisq, gsl_multifit_linear_workspace * work)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
/* FIXME: this call does actually perform balancing, but this function is deprecated anyway */
|
|
Packit |
67cb25 |
int status = gsl_multifit_wlinear_tsvd(X, w, y, tol, c, cov, chisq, rank, work);
|
|
Packit |
67cb25 |
return status;
|
|
Packit |
67cb25 |
}
|