Blame linalg/lq.c

Packit 67cb25
/* linalg/lq.c
Packit 67cb25
 * 
Packit 67cb25
 * Copyright (C) 1996, 1997, 1998, 1999, 2000, 2007 Gerard Jungman, Brian Gough
Packit 67cb25
 * Copyright (C) 2004 Joerg Wensch, modifications for LQ.
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 <stdlib.h>
Packit 67cb25
#include <string.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_blas.h>
Packit 67cb25
Packit 67cb25
#include <gsl/gsl_linalg.h>
Packit 67cb25
Packit 67cb25
#include "apply_givens.c"
Packit 67cb25
Packit 67cb25
/* Note: The standard in numerical linear algebra is to solve A x = b
Packit 67cb25
 * resp. ||A x - b||_2 -> min by QR-decompositions where x, b are
Packit 67cb25
 * column vectors.
Packit 67cb25
 *
Packit 67cb25
 * When the matrix A has a large number of rows it is much more
Packit 67cb25
 * efficient to work with the transposed matrix A^T and to solve the
Packit 67cb25
 * system x^T A = b^T resp. ||x^T A - b^T||_2 -> min.  This is caused
Packit 67cb25
 * by the row-oriented format in which GSL stores matrices.  Therefore
Packit 67cb25
 * the QR-decomposition of A has to be replaced by a LQ decomposition
Packit 67cb25
 * of A^T
Packit 67cb25
 *
Packit 67cb25
 * The purpose of this package is to provide the algorithms to compute
Packit 67cb25
 * the LQ-decomposition and to solve the linear equations resp. least
Packit 67cb25
 * squares problems.  The dimensions N, M of the matrix are switched
Packit 67cb25
 * because here A will probably be a transposed matrix.  We write x^T,
Packit 67cb25
 * b^T,... for vectors the comments to emphasize that they are row
Packit 67cb25
 * vectors.
Packit 67cb25
 *
Packit 67cb25
 * It may even be useful to transpose your matrix explicitly (assumed
Packit 67cb25
 * that there are no memory restrictions) because this takes O(M x N)
Packit 67cb25
 * computing time where the decompostion takes O(M x N^2) computing
Packit 67cb25
 * time.  */
Packit 67cb25
Packit 67cb25
/* Factorise a general N x M matrix A into
Packit 67cb25
 *  
Packit 67cb25
 *   A = L Q
Packit 67cb25
 *
Packit 67cb25
 * where Q is orthogonal (M x M) and L is lower triangular (N x M).
Packit 67cb25
 *
Packit 67cb25
 * Q is stored as a packed set of Householder transformations in the
Packit 67cb25
 * strict upper triangular part of the input matrix.
Packit 67cb25
 *
Packit 67cb25
 * R is stored in the diagonal and lower triangle of the input matrix.
Packit 67cb25
 *
Packit 67cb25
 * The full matrix for Q can be obtained as the product
Packit 67cb25
 *
Packit 67cb25
 *       Q = Q_k .. Q_2 Q_1
Packit 67cb25
 *
Packit 67cb25
 * where k = MIN(M,N) and
Packit 67cb25
 *
Packit 67cb25
 *       Q_i = (I - tau_i * v_i * v_i')
Packit 67cb25
 *
Packit 67cb25
 * and where v_i is a Householder vector
Packit 67cb25
 *
Packit 67cb25
 *       v_i = [1, m(i+1,i), m(i+2,i), ... , m(M,i)]
Packit 67cb25
 *
Packit 67cb25
 * This storage scheme is the same as in LAPACK.  */
Packit 67cb25
Packit 67cb25
int
Packit 67cb25
gsl_linalg_LQ_decomp (gsl_matrix * A, gsl_vector * tau)
Packit 67cb25
{
Packit 67cb25
  const size_t N = A->size1;
Packit 67cb25
  const size_t M = A->size2;
Packit 67cb25
Packit 67cb25
  if (tau->size != GSL_MIN (M, N))
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("size of tau must be MIN(M,N)", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else
Packit 67cb25
    {
Packit 67cb25
      size_t i;
Packit 67cb25
Packit 67cb25
      for (i = 0; i < GSL_MIN (M, N); i++)
Packit 67cb25
        {
Packit 67cb25
          /* Compute the Householder transformation to reduce the j-th
Packit 67cb25
             column of the matrix to a multiple of the j-th unit vector */
Packit 67cb25
Packit 67cb25
          gsl_vector_view c_full = gsl_matrix_row (A, i);
Packit 67cb25
          gsl_vector_view c = gsl_vector_subvector (&(c_full.vector), i, M-i);
Packit 67cb25
Packit 67cb25
          double tau_i = gsl_linalg_householder_transform (&(c.vector));
Packit 67cb25
Packit 67cb25
          gsl_vector_set (tau, i, tau_i);
Packit 67cb25
Packit 67cb25
          /* Apply the transformation to the remaining columns and
Packit 67cb25
             update the norms */
Packit 67cb25
Packit 67cb25
          if (i + 1 < N)
Packit 67cb25
            {
Packit 67cb25
              gsl_matrix_view m = gsl_matrix_submatrix (A, i + 1, i, N - (i + 1), M - i );
Packit 67cb25
              gsl_linalg_householder_mh (tau_i, &(c.vector), &(m.matrix));
Packit 67cb25
            }
Packit 67cb25
        }
Packit 67cb25
Packit 67cb25
      return GSL_SUCCESS;
Packit 67cb25
    }
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
/* Solves the system x^T A = b^T using the LQ factorisation,
Packit 67cb25
Packit 67cb25
 *  x^T L = b^T Q^T
Packit 67cb25
 *
Packit 67cb25
 * to obtain x. Based on SLATEC code. 
Packit 67cb25
 */
Packit 67cb25
Packit 67cb25
Packit 67cb25
int
Packit 67cb25
gsl_linalg_LQ_solve_T (const gsl_matrix * LQ, const gsl_vector * tau, const gsl_vector * b, gsl_vector * x)
Packit 67cb25
{
Packit 67cb25
  if (LQ->size1 != LQ->size2)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("LQ matrix must be square", GSL_ENOTSQR);
Packit 67cb25
    }
Packit 67cb25
  else if (LQ->size2 != b->size)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("matrix size must match b size", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else if (LQ->size1 != x->size)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("matrix size must match solution size", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else
Packit 67cb25
    {
Packit 67cb25
      /* Copy x <- b */
Packit 67cb25
Packit 67cb25
      gsl_vector_memcpy (x, b);
Packit 67cb25
Packit 67cb25
      /* Solve for x */
Packit 67cb25
Packit 67cb25
      gsl_linalg_LQ_svx_T (LQ, tau, x);
Packit 67cb25
Packit 67cb25
      return GSL_SUCCESS;
Packit 67cb25
    }
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
/* Solves the system x^T A = b^T in place using the LQ factorisation,
Packit 67cb25
 *
Packit 67cb25
 *  x^T L = b^T Q^T
Packit 67cb25
 *
Packit 67cb25
 * to obtain x. Based on SLATEC code.
Packit 67cb25
 */
Packit 67cb25
Packit 67cb25
int
Packit 67cb25
gsl_linalg_LQ_svx_T (const gsl_matrix * LQ, const gsl_vector * tau, gsl_vector * x)
Packit 67cb25
{
Packit 67cb25
Packit 67cb25
  if (LQ->size1 != LQ->size2)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("LQ matrix must be square", GSL_ENOTSQR);
Packit 67cb25
    }
Packit 67cb25
  else if (LQ->size1 != x->size)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("matrix size must match x/rhs size", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else
Packit 67cb25
    {
Packit 67cb25
      /* compute rhs = Q^T b */
Packit 67cb25
Packit 67cb25
      gsl_linalg_LQ_vecQT (LQ, tau, x);
Packit 67cb25
Packit 67cb25
      /* Solve R x = rhs, storing x in-place */
Packit 67cb25
Packit 67cb25
      gsl_blas_dtrsv (CblasLower, CblasTrans, CblasNonUnit, LQ, x);
Packit 67cb25
Packit 67cb25
      return GSL_SUCCESS;
Packit 67cb25
    }
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
Packit 67cb25
/* Find the least squares solution to the overdetermined system
Packit 67cb25
 *
Packit 67cb25
 *   x^T A = b^T
Packit 67cb25
 *
Packit 67cb25
 * for M >= N using the LQ factorization A = L Q.
Packit 67cb25
 */
Packit 67cb25
Packit 67cb25
int
Packit 67cb25
gsl_linalg_LQ_lssolve_T (const gsl_matrix * LQ, const gsl_vector * tau, const gsl_vector * b, gsl_vector * x, gsl_vector * residual)
Packit 67cb25
{
Packit 67cb25
  const size_t N = LQ->size1;
Packit 67cb25
  const size_t M = LQ->size2;
Packit 67cb25
Packit 67cb25
  if (M < N)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("LQ matrix must have M>=N", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else if (M != b->size)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("matrix size must match b size", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else if (N != x->size)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("matrix size must match solution size", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else if (M != residual->size)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("matrix size must match residual size", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else
Packit 67cb25
    {
Packit 67cb25
      gsl_matrix_const_view L = gsl_matrix_const_submatrix (LQ, 0, 0, N, N);
Packit 67cb25
      gsl_vector_view c = gsl_vector_subvector(residual, 0, N);
Packit 67cb25
Packit 67cb25
      gsl_vector_memcpy(residual, b);
Packit 67cb25
Packit 67cb25
      /* compute rhs = b^T Q^T */
Packit 67cb25
Packit 67cb25
      gsl_linalg_LQ_vecQT (LQ, tau, residual);
Packit 67cb25
Packit 67cb25
      /* Solve x^T L = rhs */
Packit 67cb25
Packit 67cb25
      gsl_vector_memcpy(x, &(c.vector));
Packit 67cb25
Packit 67cb25
      gsl_blas_dtrsv (CblasLower, CblasTrans, CblasNonUnit, &(L.matrix), x);
Packit 67cb25
Packit 67cb25
      /* Compute residual = b^T - x^T A = (b^T Q^T - x^T L) Q */
Packit 67cb25
      
Packit 67cb25
      gsl_vector_set_zero(&(c.vector));
Packit 67cb25
Packit 67cb25
      gsl_linalg_LQ_vecQ(LQ, tau, residual);
Packit 67cb25
Packit 67cb25
      return GSL_SUCCESS;
Packit 67cb25
    }
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
Packit 67cb25
int
Packit 67cb25
gsl_linalg_LQ_Lsolve_T (const gsl_matrix * LQ, const gsl_vector * b, gsl_vector * x)
Packit 67cb25
{
Packit 67cb25
  if (LQ->size1 != LQ->size2)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("LQ matrix must be square", GSL_ENOTSQR);
Packit 67cb25
    }
Packit 67cb25
  else if (LQ->size1 != b->size)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("matrix size must match b size", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else if (LQ->size1 != x->size)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("matrix size must match x size", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else
Packit 67cb25
    {
Packit 67cb25
      /* Copy x <- b */
Packit 67cb25
Packit 67cb25
      gsl_vector_memcpy (x, b);
Packit 67cb25
Packit 67cb25
      /* Solve R x = b, storing x in-place */
Packit 67cb25
Packit 67cb25
      gsl_blas_dtrsv (CblasLower, CblasTrans, CblasNonUnit, LQ, x);
Packit 67cb25
Packit 67cb25
      return GSL_SUCCESS;
Packit 67cb25
    }
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
Packit 67cb25
int
Packit 67cb25
gsl_linalg_LQ_Lsvx_T (const gsl_matrix * LQ, gsl_vector * x)
Packit 67cb25
{
Packit 67cb25
  if (LQ->size1 != LQ->size2)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("LQ matrix must be square", GSL_ENOTSQR);
Packit 67cb25
    }
Packit 67cb25
  else if (LQ->size2 != x->size)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("matrix size must match rhs size", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else
Packit 67cb25
    {
Packit 67cb25
      /* Solve x^T L = b^T, storing x in-place */
Packit 67cb25
Packit 67cb25
      gsl_blas_dtrsv (CblasLower, CblasTrans, CblasNonUnit, LQ, x);
Packit 67cb25
Packit 67cb25
      return GSL_SUCCESS;
Packit 67cb25
    }
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
int
Packit 67cb25
gsl_linalg_L_solve_T (const gsl_matrix * L, const gsl_vector * b, gsl_vector * x)
Packit 67cb25
{
Packit 67cb25
  if (L->size1 != L->size2)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("R matrix must be square", GSL_ENOTSQR);
Packit 67cb25
    }
Packit 67cb25
  else if (L->size2 != b->size)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("matrix size must match b size", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else if (L->size1 != x->size)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("matrix size must match solution size", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else
Packit 67cb25
    {
Packit 67cb25
      /* Copy x <- b */
Packit 67cb25
Packit 67cb25
      gsl_vector_memcpy (x, b);
Packit 67cb25
Packit 67cb25
      /* Solve R x = b, storing x inplace in b */
Packit 67cb25
Packit 67cb25
      gsl_blas_dtrsv (CblasLower, CblasTrans, CblasNonUnit, L, x);
Packit 67cb25
Packit 67cb25
      return GSL_SUCCESS;
Packit 67cb25
    }
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
Packit 67cb25
Packit 67cb25
Packit 67cb25
int
Packit 67cb25
gsl_linalg_LQ_vecQT (const gsl_matrix * LQ, const gsl_vector * tau, gsl_vector * v)
Packit 67cb25
{
Packit 67cb25
  const size_t N = LQ->size1;
Packit 67cb25
  const size_t M = LQ->size2;
Packit 67cb25
Packit 67cb25
  if (tau->size != GSL_MIN (M, N))
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("size of tau must be MIN(M,N)", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else if (v->size != M)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("vector size must be M", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else
Packit 67cb25
    {
Packit 67cb25
      size_t i;
Packit 67cb25
Packit 67cb25
      /* compute v Q^T  */
Packit 67cb25
Packit 67cb25
      for (i = 0; i < GSL_MIN (M, N); i++)
Packit 67cb25
        {
Packit 67cb25
          gsl_vector_const_view c = gsl_matrix_const_row (LQ, i);
Packit 67cb25
          gsl_vector_const_view h = gsl_vector_const_subvector (&(c.vector),
Packit 67cb25
                                                                i, M - i);
Packit 67cb25
          gsl_vector_view w = gsl_vector_subvector (v, i, M - i);
Packit 67cb25
          double ti = gsl_vector_get (tau, i);
Packit 67cb25
          gsl_linalg_householder_hv (ti, &(h.vector), &(w.vector));
Packit 67cb25
        }
Packit 67cb25
      return GSL_SUCCESS;
Packit 67cb25
    }
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
int
Packit 67cb25
gsl_linalg_LQ_vecQ (const gsl_matrix * LQ, const gsl_vector * tau, gsl_vector * v)
Packit 67cb25
{
Packit 67cb25
  const size_t N = LQ->size1;
Packit 67cb25
  const size_t M = LQ->size2;
Packit 67cb25
Packit 67cb25
  if (tau->size != GSL_MIN (M, N))
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("size of tau must be MIN(M,N)", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else if (v->size != M)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("vector size must be M", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else
Packit 67cb25
    {
Packit 67cb25
      size_t i;
Packit 67cb25
Packit 67cb25
      /* compute v Q^T  */
Packit 67cb25
      
Packit 67cb25
      for (i =  GSL_MIN (M, N); i-- > 0;) 
Packit 67cb25
        {
Packit 67cb25
          gsl_vector_const_view c = gsl_matrix_const_row (LQ, i);
Packit 67cb25
          gsl_vector_const_view h = gsl_vector_const_subvector (&(c.vector),
Packit 67cb25
                                                                i, M - i);
Packit 67cb25
          gsl_vector_view w = gsl_vector_subvector (v, i, M - i);
Packit 67cb25
          double ti = gsl_vector_get (tau, i);
Packit 67cb25
          gsl_linalg_householder_hv (ti, &(h.vector), &(w.vector));
Packit 67cb25
        }
Packit 67cb25
      return GSL_SUCCESS;
Packit 67cb25
    }
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
Packit 67cb25
/*  Form the orthogonal matrix Q from the packed LQ matrix */
Packit 67cb25
Packit 67cb25
int
Packit 67cb25
gsl_linalg_LQ_unpack (const gsl_matrix * LQ, const gsl_vector * tau, gsl_matrix * Q, gsl_matrix * L)
Packit 67cb25
{
Packit 67cb25
  const size_t N = LQ->size1;
Packit 67cb25
  const size_t M = LQ->size2;
Packit 67cb25
Packit 67cb25
  if (Q->size1 != M || Q->size2 != M)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("Q matrix must be M x M", GSL_ENOTSQR);
Packit 67cb25
    }
Packit 67cb25
  else if (L->size1 != N || L->size2 != M)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("R matrix must be N x M", GSL_ENOTSQR);
Packit 67cb25
    }
Packit 67cb25
  else if (tau->size != GSL_MIN (M, N))
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("size of tau must be MIN(M,N)", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else
Packit 67cb25
    {
Packit 67cb25
      size_t i, j, l_border;
Packit 67cb25
Packit 67cb25
      /* Initialize Q to the identity */
Packit 67cb25
Packit 67cb25
      gsl_matrix_set_identity (Q);
Packit 67cb25
Packit 67cb25
      for (i = GSL_MIN (M, N); i-- > 0;)
Packit 67cb25
        {
Packit 67cb25
          gsl_vector_const_view c = gsl_matrix_const_row (LQ, i);
Packit 67cb25
          gsl_vector_const_view h = gsl_vector_const_subvector (&c.vector,
Packit 67cb25
                                                                i, M - i);
Packit 67cb25
          gsl_matrix_view m = gsl_matrix_submatrix (Q, i, i, M - i, M - i);
Packit 67cb25
          double ti = gsl_vector_get (tau, i);
Packit 67cb25
          gsl_linalg_householder_mh (ti, &h.vector, &m.matrix);
Packit 67cb25
        }
Packit 67cb25
Packit 67cb25
      /*  Form the lower triangular matrix L from a packed LQ matrix */
Packit 67cb25
Packit 67cb25
      for (i = 0; i < N; i++)
Packit 67cb25
        {
Packit 67cb25
	    l_border=GSL_MIN(i,M-1);
Packit 67cb25
		for (j = 0; j <= l_border ; j++)
Packit 67cb25
		    gsl_matrix_set (L, i, j, gsl_matrix_get (LQ, i, j));
Packit 67cb25
Packit 67cb25
	    for (j = l_border+1; j < M; j++)
Packit 67cb25
		gsl_matrix_set (L, i, j, 0.0);
Packit 67cb25
        }
Packit 67cb25
Packit 67cb25
      return GSL_SUCCESS;
Packit 67cb25
    }
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
Packit 67cb25
/* Update a LQ factorisation for A= L Q ,  A' = A + v u^T,
Packit 67cb25
Packit 67cb25
 * L' Q' = LQ + v u^T
Packit 67cb25
 *       = (L + v u^T Q^T) Q
Packit 67cb25
 *       = (L + v w^T) Q
Packit 67cb25
 *
Packit 67cb25
 * where w = Q u.
Packit 67cb25
 *
Packit 67cb25
 * Algorithm from Golub and Van Loan, "Matrix Computations", Section
Packit 67cb25
 * 12.5 (Updating Matrix Factorizations, Rank-One Changes)
Packit 67cb25
 */
Packit 67cb25
Packit 67cb25
int
Packit 67cb25
gsl_linalg_LQ_update (gsl_matrix * Q, gsl_matrix * L,
Packit 67cb25
                      const gsl_vector * v, gsl_vector * w)
Packit 67cb25
{
Packit 67cb25
  const size_t N = L->size1;
Packit 67cb25
  const size_t M = L->size2;
Packit 67cb25
Packit 67cb25
  if (Q->size1 != M || Q->size2 != M)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("Q matrix must be N x N if L is M x N", GSL_ENOTSQR);
Packit 67cb25
    }
Packit 67cb25
  else if (w->size != M)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("w must be length N if L is M x N", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else if (v->size != N)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("v must be length M if L is M x N", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else
Packit 67cb25
    {
Packit 67cb25
      size_t j, k;
Packit 67cb25
      double w0;
Packit 67cb25
Packit 67cb25
      /* Apply Given's rotations to reduce w to (|w|, 0, 0, ... , 0)
Packit 67cb25
Packit 67cb25
         J_1^T .... J_(n-1)^T w = +/- |w| e_1
Packit 67cb25
Packit 67cb25
         simultaneously applied to L,  H = J_1^T ... J^T_(n-1) L
Packit 67cb25
         so that H is upper Hessenberg.  (12.5.2) */
Packit 67cb25
      
Packit 67cb25
      for (k = M - 1; k > 0; k--)  /* loop from k = M-1 to 1 */
Packit 67cb25
        {
Packit 67cb25
          double c, s;
Packit 67cb25
          double wk = gsl_vector_get (w, k);
Packit 67cb25
          double wkm1 = gsl_vector_get (w, k - 1);
Packit 67cb25
Packit 67cb25
          gsl_linalg_givens (wkm1, wk, &c, &s);
Packit 67cb25
          gsl_linalg_givens_gv (w, k - 1, k, c, s);
Packit 67cb25
          apply_givens_lq (M, N, Q, L, k - 1, k, c, s);
Packit 67cb25
       }
Packit 67cb25
Packit 67cb25
      w0 = gsl_vector_get (w, 0);
Packit 67cb25
Packit 67cb25
      /* Add in v w^T  (Equation 12.5.3) */
Packit 67cb25
Packit 67cb25
      for (j = 0; j < N; j++)
Packit 67cb25
        {
Packit 67cb25
          double lj0 = gsl_matrix_get (L, j, 0);
Packit 67cb25
          double vj = gsl_vector_get (v, j);
Packit 67cb25
          gsl_matrix_set (L, j, 0, lj0 + w0 * vj);
Packit 67cb25
        }
Packit 67cb25
Packit 67cb25
      /* Apply Givens transformations L' = G_(n-1)^T ... G_1^T H
Packit 67cb25
         Equation 12.5.4 */
Packit 67cb25
Packit 67cb25
      for (k = 1; k < GSL_MIN(M,N+1); k++)
Packit 67cb25
        {
Packit 67cb25
          double c, s;
Packit 67cb25
          double diag = gsl_matrix_get (L, k - 1, k - 1);
Packit 67cb25
          double offdiag = gsl_matrix_get (L, k - 1 , k);
Packit 67cb25
Packit 67cb25
          gsl_linalg_givens (diag, offdiag, &c, &s);
Packit 67cb25
          apply_givens_lq (M, N, Q, L, k - 1, k, c, s);
Packit 67cb25
Packit 67cb25
          gsl_matrix_set (L, k - 1, k, 0.0);    /* exact zero of G^T */
Packit 67cb25
        }
Packit 67cb25
Packit 67cb25
      return GSL_SUCCESS;
Packit 67cb25
    }
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
int
Packit 67cb25
gsl_linalg_LQ_LQsolve (gsl_matrix * Q, gsl_matrix * L, const gsl_vector * b, gsl_vector * x)
Packit 67cb25
{
Packit 67cb25
  const size_t N = L->size1;
Packit 67cb25
  const size_t M = L->size2;
Packit 67cb25
Packit 67cb25
  if (M != N)
Packit 67cb25
    {
Packit 67cb25
      return GSL_ENOTSQR;
Packit 67cb25
    }
Packit 67cb25
  else if (Q->size1 != M || b->size != M || x->size != M)
Packit 67cb25
    {
Packit 67cb25
      return GSL_EBADLEN;
Packit 67cb25
    }
Packit 67cb25
  else
Packit 67cb25
    {
Packit 67cb25
      /* compute sol = b^T Q^T */
Packit 67cb25
Packit 67cb25
      gsl_blas_dgemv (CblasNoTrans, 1.0, Q, b, 0.0, x);
Packit 67cb25
Packit 67cb25
      /* Solve x^T L = sol, storing x in-place */
Packit 67cb25
Packit 67cb25
      gsl_blas_dtrsv (CblasLower, CblasTrans, CblasNonUnit, L, x);
Packit 67cb25
Packit 67cb25
      return GSL_SUCCESS;
Packit 67cb25
    }
Packit 67cb25
}