Blame multilarge/multilarge.c

Packit 67cb25
/* multilarge.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_math.h>
Packit 67cb25
#include <gsl/gsl_vector.h>
Packit 67cb25
#include <gsl/gsl_matrix.h>
Packit 67cb25
#include <gsl/gsl_errno.h>
Packit 67cb25
#include <gsl/gsl_multifit.h>
Packit 67cb25
#include <gsl/gsl_multilarge.h>
Packit 67cb25
#include <gsl/gsl_blas.h>
Packit 67cb25
Packit 67cb25
gsl_multilarge_linear_workspace *
Packit 67cb25
gsl_multilarge_linear_alloc(const gsl_multilarge_linear_type *T,
Packit 67cb25
                            const size_t p)
Packit 67cb25
{
Packit 67cb25
  gsl_multilarge_linear_workspace *w;
Packit 67cb25
Packit 67cb25
  w = calloc(1, sizeof(gsl_multilarge_linear_workspace));
Packit 67cb25
  if (w == NULL)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR_NULL("failed to allocate space for workspace",
Packit 67cb25
                     GSL_ENOMEM);
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
  w->type = T;
Packit 67cb25
Packit 67cb25
  w->state = w->type->alloc(p);
Packit 67cb25
  if (w->state == NULL)
Packit 67cb25
    {
Packit 67cb25
      gsl_multilarge_linear_free(w);
Packit 67cb25
      GSL_ERROR_NULL("failed to allocate space for multilarge state",
Packit 67cb25
                     GSL_ENOMEM);
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
  w->p = p;
Packit 67cb25
Packit 67cb25
  /* initialize newly allocated state */
Packit 67cb25
  gsl_multilarge_linear_reset(w);
Packit 67cb25
Packit 67cb25
  return w;
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
void
Packit 67cb25
gsl_multilarge_linear_free(gsl_multilarge_linear_workspace *w)
Packit 67cb25
{
Packit 67cb25
  RETURN_IF_NULL(w);
Packit 67cb25
Packit 67cb25
  if (w->state)
Packit 67cb25
    w->type->free(w->state);
Packit 67cb25
Packit 67cb25
  free(w);
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
const char *
Packit 67cb25
gsl_multilarge_linear_name(const gsl_multilarge_linear_workspace *w)
Packit 67cb25
{
Packit 67cb25
  return w->type->name;
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
int
Packit 67cb25
gsl_multilarge_linear_reset(gsl_multilarge_linear_workspace *w)
Packit 67cb25
{
Packit 67cb25
  int status = w->type->reset(w->state);
Packit 67cb25
  return status;
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
int
Packit 67cb25
gsl_multilarge_linear_accumulate(gsl_matrix * X, gsl_vector * y,
Packit 67cb25
                                 gsl_multilarge_linear_workspace * w)
Packit 67cb25
{
Packit 67cb25
  int status = w->type->accumulate(X, y, w->state);
Packit 67cb25
  return status;
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
int
Packit 67cb25
gsl_multilarge_linear_solve(const double lambda, gsl_vector * c,
Packit 67cb25
                            double * rnorm, double * snorm,
Packit 67cb25
                            gsl_multilarge_linear_workspace * w)
Packit 67cb25
{
Packit 67cb25
  int status = w->type->solve(lambda, c, rnorm, snorm, w->state);
Packit 67cb25
  return status;
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
int
Packit 67cb25
gsl_multilarge_linear_rcond(double *rcond, gsl_multilarge_linear_workspace * w)
Packit 67cb25
{
Packit 67cb25
  int status = w->type->rcond(rcond, w->state);
Packit 67cb25
  return status;
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
int
Packit 67cb25
gsl_multilarge_linear_lcurve(gsl_vector * reg_param, gsl_vector * rho,
Packit 67cb25
                             gsl_vector * eta,
Packit 67cb25
                             gsl_multilarge_linear_workspace * w)
Packit 67cb25
{
Packit 67cb25
  const size_t len = reg_param->size;
Packit 67cb25
Packit 67cb25
  if (len != rho->size)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("reg_param and rho have different sizes", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else if (len != eta->size)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("reg_param and eta have different sizes", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else
Packit 67cb25
    {
Packit 67cb25
      int status = w->type->lcurve(reg_param, rho, eta, w->state);
Packit 67cb25
      return status;
Packit 67cb25
    }
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
/*
Packit 67cb25
gsl_multilarge_linear_wstdform1()
Packit 67cb25
  Using regularization matrix
Packit 67cb25
L = diag(l_1,l_2,...,l_p), transform to Tikhonov standard form:
Packit 67cb25
Packit 67cb25
X~ = sqrt(W) X L^{-1}
Packit 67cb25
y~ = sqrt(W) y
Packit 67cb25
c~ = L c
Packit 67cb25
Packit 67cb25
Inputs: L    - Tikhonov matrix as a vector of diagonal elements p-by-1;
Packit 67cb25
               or NULL for L = I
Packit 67cb25
        X    - least squares matrix n-by-p
Packit 67cb25
        y    - right hand side vector n-by-1
Packit 67cb25
        w    - weight vector n-by-1; or NULL for W = I
Packit 67cb25
        Xs   - least squares matrix in standard form X~ n-by-p
Packit 67cb25
        ys   - right hand side vector in standard form y~ n-by-1
Packit 67cb25
        work - workspace
Packit 67cb25
Packit 67cb25
Return: success/error
Packit 67cb25
Packit 67cb25
Notes:
Packit 67cb25
1) It is allowed for X = Xs and y = ys
Packit 67cb25
*/
Packit 67cb25
Packit 67cb25
int
Packit 67cb25
gsl_multilarge_linear_wstdform1 (const gsl_vector * L,
Packit 67cb25
                                 const gsl_matrix * X,
Packit 67cb25
                                 const gsl_vector * w,
Packit 67cb25
                                 const gsl_vector * y,
Packit 67cb25
                                 gsl_matrix * Xs,
Packit 67cb25
                                 gsl_vector * ys,
Packit 67cb25
                                 gsl_multilarge_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 (L != NULL && p != L->size)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR("L vector does not match X", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else if (n != y->size)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR("y vector does not match X", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else if (w != NULL && n != w->size)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR("weight vector does not match X", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else if (n != Xs->size1 || p != Xs->size2)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR("Xs matrix dimensions do not match X", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else if (n != ys->size)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR("ys vector must be length n", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else
Packit 67cb25
    {
Packit 67cb25
      int status = GSL_SUCCESS;
Packit 67cb25
Packit 67cb25
      /* compute Xs = sqrt(W) X and ys = sqrt(W) y */
Packit 67cb25
      status = gsl_multifit_linear_applyW(X, w, y, Xs, ys);
Packit 67cb25
      if (status)
Packit 67cb25
        return status;
Packit 67cb25
Packit 67cb25
      if (L != NULL)
Packit 67cb25
        {
Packit 67cb25
          size_t j;
Packit 67cb25
Packit 67cb25
          /* construct X~ = sqrt(W) X * L^{-1} matrix */
Packit 67cb25
          for (j = 0; j < p; ++j)
Packit 67cb25
            {
Packit 67cb25
              gsl_vector_view Xj = gsl_matrix_column(Xs, j);
Packit 67cb25
              double lj = gsl_vector_get(L, j);
Packit 67cb25
Packit 67cb25
              if (lj == 0.0)
Packit 67cb25
                {
Packit 67cb25
                  GSL_ERROR("L matrix is singular", GSL_EDOM);
Packit 67cb25
                }
Packit 67cb25
Packit 67cb25
              gsl_vector_scale(&Xj.vector, 1.0 / lj);
Packit 67cb25
            }
Packit 67cb25
        }
Packit 67cb25
Packit 67cb25
      return status;
Packit 67cb25
    }
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
int
Packit 67cb25
gsl_multilarge_linear_stdform1 (const gsl_vector * L,
Packit 67cb25
                                const gsl_matrix * X,
Packit 67cb25
                                const gsl_vector * y,
Packit 67cb25
                                gsl_matrix * Xs,
Packit 67cb25
                                gsl_vector * ys,
Packit 67cb25
                                gsl_multilarge_linear_workspace * work)
Packit 67cb25
{
Packit 67cb25
  int status;
Packit 67cb25
Packit 67cb25
  status = gsl_multilarge_linear_wstdform1(L, X, NULL, y, Xs, ys, work);
Packit 67cb25
Packit 67cb25
  return status;
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
int
Packit 67cb25
gsl_multilarge_linear_L_decomp (gsl_matrix * L, gsl_vector * tau)
Packit 67cb25
{
Packit 67cb25
  const size_t m = L->size1;
Packit 67cb25
  const size_t p = L->size2;
Packit 67cb25
Packit 67cb25
  if (m < p)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR("m < p not yet supported", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else
Packit 67cb25
    {
Packit 67cb25
      int status;
Packit 67cb25
Packit 67cb25
      status = gsl_multifit_linear_L_decomp(L, tau);
Packit 67cb25
Packit 67cb25
      return status;
Packit 67cb25
    }
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
int
Packit 67cb25
gsl_multilarge_linear_wstdform2 (const gsl_matrix * LQR,
Packit 67cb25
                                 const gsl_vector * Ltau,
Packit 67cb25
                                 const gsl_matrix * X,
Packit 67cb25
                                 const gsl_vector * w,
Packit 67cb25
                                 const gsl_vector * y,
Packit 67cb25
                                 gsl_matrix * Xs,
Packit 67cb25
                                 gsl_vector * ys,
Packit 67cb25
                                 gsl_multilarge_linear_workspace * work)
Packit 67cb25
{
Packit 67cb25
  const size_t m = LQR->size1;
Packit 67cb25
  const size_t n = X->size1;
Packit 67cb25
  const size_t p = X->size2;
Packit 67cb25
Packit 67cb25
  if (p != work->p)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR("X has wrong number of columns", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else if (p != LQR->size2)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR("LQR and X matrices have different numbers of columns", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else if (n != y->size)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR("y vector does not match X", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else if (w != NULL && n != w->size)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR("weights vector must be length n", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else if (m < p)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR("m < p not yet supported", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else if (n != Xs->size1 || p != Xs->size2)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR("Xs matrix must be n-by-p", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else if (n != ys->size)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR("ys vector must have length n", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else
Packit 67cb25
    {
Packit 67cb25
      int status;
Packit 67cb25
      size_t i;
Packit 67cb25
      gsl_matrix_const_view R = gsl_matrix_const_submatrix(LQR, 0, 0, p, p);
Packit 67cb25
Packit 67cb25
      /* compute Xs = sqrt(W) X and ys = sqrt(W) y */
Packit 67cb25
      status = gsl_multifit_linear_applyW(X, w, y, Xs, ys);
Packit 67cb25
      if (status)
Packit 67cb25
        return status;
Packit 67cb25
Packit 67cb25
      /* compute X~ = X R^{-1} using QR decomposition of L */
Packit 67cb25
      for (i = 0; i < n; ++i)
Packit 67cb25
        {
Packit 67cb25
          gsl_vector_view v = gsl_matrix_row(Xs, i);
Packit 67cb25
Packit 67cb25
          /* solve: R^T y = X_i */
Packit 67cb25
          gsl_blas_dtrsv(CblasUpper, CblasTrans, CblasNonUnit, &R.matrix, &v.vector);
Packit 67cb25
        }
Packit 67cb25
Packit 67cb25
      return GSL_SUCCESS;
Packit 67cb25
    }
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
int
Packit 67cb25
gsl_multilarge_linear_stdform2 (const gsl_matrix * LQR,
Packit 67cb25
                                const gsl_vector * Ltau,
Packit 67cb25
                                const gsl_matrix * X,
Packit 67cb25
                                const gsl_vector * y,
Packit 67cb25
                                gsl_matrix * Xs,
Packit 67cb25
                                gsl_vector * ys,
Packit 67cb25
                                gsl_multilarge_linear_workspace * work)
Packit 67cb25
{
Packit 67cb25
  int status;
Packit 67cb25
Packit 67cb25
  status = gsl_multilarge_linear_wstdform2(LQR, Ltau, X, NULL, y, Xs, ys, work);
Packit 67cb25
Packit 67cb25
  return status;
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
/*
Packit 67cb25
gsl_multilarge_linear_genform1()
Packit 67cb25
  Backtransform regularized solution vector using matrix
Packit 67cb25
L = diag(L)
Packit 67cb25
*/
Packit 67cb25
Packit 67cb25
int
Packit 67cb25
gsl_multilarge_linear_genform1 (const gsl_vector * L,
Packit 67cb25
                                const gsl_vector * cs,
Packit 67cb25
                                gsl_vector * c,
Packit 67cb25
                                gsl_multilarge_linear_workspace * work)
Packit 67cb25
{
Packit 67cb25
  if (L->size != work->p)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR("L vector does not match workspace", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else if (L->size != cs->size)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR("cs vector does not match L", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else if (L->size != c->size)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR("c vector does not match L", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else
Packit 67cb25
    {
Packit 67cb25
      /* compute true solution vector c = L^{-1} c~ */
Packit 67cb25
      gsl_vector_memcpy(c, cs);
Packit 67cb25
      gsl_vector_div(c, L);
Packit 67cb25
Packit 67cb25
      return GSL_SUCCESS;
Packit 67cb25
    }
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
int
Packit 67cb25
gsl_multilarge_linear_genform2 (const gsl_matrix * LQR,
Packit 67cb25
                                const gsl_vector * Ltau,
Packit 67cb25
                                const gsl_vector * cs,
Packit 67cb25
                                gsl_vector * c,
Packit 67cb25
                                gsl_multilarge_linear_workspace * work)
Packit 67cb25
{
Packit 67cb25
  const size_t m = LQR->size1;
Packit 67cb25
  const size_t p = LQR->size2;
Packit 67cb25
Packit 67cb25
  if (p != c->size)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR("c vector does not match LQR", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else if (m < p)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR("m < p not yet supported", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else if (p != cs->size)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR("cs vector size does not match c", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else
Packit 67cb25
    {
Packit 67cb25
      int s;
Packit 67cb25
      gsl_matrix_const_view R = gsl_matrix_const_submatrix(LQR, 0, 0, p, p); /* R factor of L */
Packit 67cb25
Packit 67cb25
      /* solve R c = cs for true solution c, using QR decomposition of L */
Packit 67cb25
      gsl_vector_memcpy(c, cs);
Packit 67cb25
      s = gsl_blas_dtrsv(CblasUpper, CblasNoTrans, CblasNonUnit, &R.matrix, c);
Packit 67cb25
Packit 67cb25
      return s;
Packit 67cb25
    }
Packit 67cb25
}