Blame linalg/qrpt.c

Packit 67cb25
/* linalg/qrpt.c
Packit 67cb25
 * 
Packit 67cb25
 * Copyright (C) 1996, 1997, 1998, 1999, 2000, 2007 Gerard Jungman, Brian Gough
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_permute_vector.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
/* Factorise a general M x N matrix A into
Packit 67cb25
 *
Packit 67cb25
 *   A P = Q R
Packit 67cb25
 *
Packit 67cb25
 * where Q is orthogonal (M x M) and R is upper triangular (M x N).
Packit 67cb25
 * When A is rank deficient, r = rank(A) < n, then the permutation is
Packit 67cb25
 * used to ensure that the lower n - r rows of R are zero and the first
Packit 67cb25
 * r columns of Q form an orthonormal basis for A.
Packit 67cb25
 *
Packit 67cb25
 * Q is stored as a packed set of Householder transformations in the
Packit 67cb25
 * strict lower triangular part of the input matrix.
Packit 67cb25
 *
Packit 67cb25
 * R is stored in the diagonal and upper triangle of the input matrix.
Packit 67cb25
 *
Packit 67cb25
 * P: column j of P is column k of the identity matrix, where k =
Packit 67cb25
 * permutation->data[j]
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.  See LAPACK's
Packit 67cb25
 * dgeqpf.f for details.
Packit 67cb25
 * 
Packit 67cb25
 */
Packit 67cb25
Packit 67cb25
int
Packit 67cb25
gsl_linalg_QRPT_decomp (gsl_matrix * A, gsl_vector * tau, gsl_permutation * p, int *signum, gsl_vector * norm)
Packit 67cb25
{
Packit 67cb25
  const size_t M = A->size1;
Packit 67cb25
  const size_t N = 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 if (p->size != N)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("permutation size must be N", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else if (norm->size != N)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("norm size must be N", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else
Packit 67cb25
    {
Packit 67cb25
      size_t i;
Packit 67cb25
Packit 67cb25
      *signum = 1;
Packit 67cb25
Packit 67cb25
      gsl_permutation_init (p); /* set to identity */
Packit 67cb25
Packit 67cb25
      /* Compute column norms and store in workspace */
Packit 67cb25
Packit 67cb25
      for (i = 0; i < N; i++)
Packit 67cb25
        {
Packit 67cb25
          gsl_vector_view c = gsl_matrix_column (A, i);
Packit 67cb25
          double x = gsl_blas_dnrm2 (&c.vector);
Packit 67cb25
          gsl_vector_set (norm, i, x);
Packit 67cb25
        }
Packit 67cb25
Packit 67cb25
      for (i = 0; i < GSL_MIN (M, N); i++)
Packit 67cb25
        {
Packit 67cb25
          /* Bring the column of largest norm into the pivot position */
Packit 67cb25
Packit 67cb25
          double max_norm = gsl_vector_get(norm, i);
Packit 67cb25
          size_t j, kmax = i;
Packit 67cb25
Packit 67cb25
          for (j = i + 1; j < N; j++)
Packit 67cb25
            {
Packit 67cb25
              double x = gsl_vector_get (norm, j);
Packit 67cb25
Packit 67cb25
              if (x > max_norm)
Packit 67cb25
                {
Packit 67cb25
                  max_norm = x;
Packit 67cb25
                  kmax = j;
Packit 67cb25
                }
Packit 67cb25
            }
Packit 67cb25
Packit 67cb25
          if (kmax != i)
Packit 67cb25
            {
Packit 67cb25
              gsl_matrix_swap_columns (A, i, kmax);
Packit 67cb25
              gsl_permutation_swap (p, i, kmax);
Packit 67cb25
              gsl_vector_swap_elements(norm,i,kmax);
Packit 67cb25
Packit 67cb25
              (*signum) = -(*signum);
Packit 67cb25
            }
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
          {
Packit 67cb25
            gsl_vector_view c_full = gsl_matrix_column (A, i);
Packit 67cb25
            gsl_vector_view c = gsl_vector_subvector (&c_full.vector, 
Packit 67cb25
                                                      i, M - i);
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 */
Packit 67cb25
Packit 67cb25
            if (i + 1 < N)
Packit 67cb25
              {
Packit 67cb25
                gsl_matrix_view m = gsl_matrix_submatrix (A, i, i + 1, M - i, N - (i+1));
Packit 67cb25
Packit 67cb25
                gsl_linalg_householder_hm (tau_i, &c.vector, &m.matrix);
Packit 67cb25
              }
Packit 67cb25
          }
Packit 67cb25
Packit 67cb25
          /* Update the norms of the remaining columns too */
Packit 67cb25
Packit 67cb25
          if (i + 1 < M) 
Packit 67cb25
            {
Packit 67cb25
              for (j = i + 1; j < N; j++)
Packit 67cb25
                {
Packit 67cb25
                  double x = gsl_vector_get (norm, j);
Packit 67cb25
Packit 67cb25
                  if (x > 0.0)
Packit 67cb25
                    {
Packit 67cb25
                      double y = 0;
Packit 67cb25
                      double temp= gsl_matrix_get (A, i, j) / x;
Packit 67cb25
                  
Packit 67cb25
                      if (fabs (temp) >= 1)
Packit 67cb25
                        y = 0.0;
Packit 67cb25
                      else
Packit 67cb25
                        y = x * sqrt (1 - temp * temp);
Packit 67cb25
                      
Packit 67cb25
                      /* recompute norm to prevent loss of accuracy */
Packit 67cb25
Packit 67cb25
                      if (fabs (y / x) < sqrt (20.0) * GSL_SQRT_DBL_EPSILON)
Packit 67cb25
                        {
Packit 67cb25
                          gsl_vector_view c_full = gsl_matrix_column (A, j);
Packit 67cb25
                          gsl_vector_view c = 
Packit 67cb25
                            gsl_vector_subvector(&c_full.vector,
Packit 67cb25
                                                 i+1, M - (i+1));
Packit 67cb25
                          y = gsl_blas_dnrm2 (&c.vector);
Packit 67cb25
                        }
Packit 67cb25
                  
Packit 67cb25
                      gsl_vector_set (norm, j, y);
Packit 67cb25
                    }
Packit 67cb25
                }
Packit 67cb25
            }
Packit 67cb25
        }
Packit 67cb25
Packit 67cb25
      return GSL_SUCCESS;
Packit 67cb25
    }
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
int
Packit 67cb25
gsl_linalg_QRPT_decomp2 (const gsl_matrix * A, gsl_matrix * q, gsl_matrix * r, gsl_vector * tau, gsl_permutation * p, int *signum, gsl_vector * norm)
Packit 67cb25
{
Packit 67cb25
  const size_t M = A->size1;
Packit 67cb25
  const size_t N = A->size2;
Packit 67cb25
Packit 67cb25
  if (q->size1 != M || q->size2 !=M) 
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("q must be M x M", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else if (r->size1 != M || r->size2 !=N)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("r must be M x N", GSL_EBADLEN);
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 if (p->size != N)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("permutation size must be N", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else if (norm->size != N)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("norm size must be N", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
  gsl_matrix_memcpy (r, A);
Packit 67cb25
Packit 67cb25
  gsl_linalg_QRPT_decomp (r, tau, p, signum, norm);
Packit 67cb25
Packit 67cb25
  /* FIXME:  aliased arguments depends on behavior of unpack routine! */
Packit 67cb25
Packit 67cb25
  gsl_linalg_QR_unpack (r, tau, q, r);
Packit 67cb25
Packit 67cb25
  return GSL_SUCCESS;
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
Packit 67cb25
/* Solves the system A x = b using the Q R P^T factorisation,
Packit 67cb25
Packit 67cb25
   R z = Q^T b
Packit 67cb25
Packit 67cb25
   x = P z;
Packit 67cb25
Packit 67cb25
   to obtain x. Based on SLATEC code. */
Packit 67cb25
Packit 67cb25
int
Packit 67cb25
gsl_linalg_QRPT_solve (const gsl_matrix * QR,
Packit 67cb25
                       const gsl_vector * tau,
Packit 67cb25
                       const gsl_permutation * p,
Packit 67cb25
                       const gsl_vector * b,
Packit 67cb25
                       gsl_vector * x)
Packit 67cb25
{
Packit 67cb25
  if (QR->size1 != QR->size2)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("QR matrix must be square", GSL_ENOTSQR);
Packit 67cb25
    }
Packit 67cb25
  else if (QR->size1 != p->size)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("matrix size must match permutation size", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else if (QR->size1 != b->size)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("matrix size must match b size", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else if (QR->size2 != 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
      gsl_vector_memcpy (x, b);
Packit 67cb25
Packit 67cb25
      gsl_linalg_QRPT_svx (QR, tau, p, x);
Packit 67cb25
      
Packit 67cb25
      return GSL_SUCCESS;
Packit 67cb25
    }
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
int
Packit 67cb25
gsl_linalg_QRPT_svx (const gsl_matrix * QR,
Packit 67cb25
                     const gsl_vector * tau,
Packit 67cb25
                     const gsl_permutation * p,
Packit 67cb25
                     gsl_vector * x)
Packit 67cb25
{
Packit 67cb25
  if (QR->size1 != QR->size2)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("QR matrix must be square", GSL_ENOTSQR);
Packit 67cb25
    }
Packit 67cb25
  else if (QR->size1 != p->size)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("matrix size must match permutation size", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else if (QR->size2 != 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
      /* compute sol = Q^T b */
Packit 67cb25
Packit 67cb25
      gsl_linalg_QR_QTvec (QR, tau, x);
Packit 67cb25
Packit 67cb25
      /* Solve R x = sol, storing x inplace in sol */
Packit 67cb25
Packit 67cb25
      gsl_blas_dtrsv (CblasUpper, CblasNoTrans, CblasNonUnit, QR, x);
Packit 67cb25
Packit 67cb25
      gsl_permute_vector_inverse (p, x);
Packit 67cb25
Packit 67cb25
      return GSL_SUCCESS;
Packit 67cb25
    }
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
/* Find the least squares solution to the overdetermined system 
Packit 67cb25
 *
Packit 67cb25
 *   A x = b 
Packit 67cb25
 *  
Packit 67cb25
 * for M >= N using the QRPT factorization A P = Q R. Assumes
Packit 67cb25
 * A has full column rank.
Packit 67cb25
 */
Packit 67cb25
Packit 67cb25
int
Packit 67cb25
gsl_linalg_QRPT_lssolve (const gsl_matrix * QR, const gsl_vector * tau, const gsl_permutation * p,
Packit 67cb25
                         const gsl_vector * b, gsl_vector * x, gsl_vector * residual)
Packit 67cb25
{
Packit 67cb25
  const size_t N = QR->size2;
Packit 67cb25
  int status = gsl_linalg_QRPT_lssolve2(QR, tau, p, b, N, x, residual);
Packit 67cb25
  return status;
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
/* Find the least squares solution to the overdetermined system 
Packit 67cb25
 *
Packit 67cb25
 *   A x = b 
Packit 67cb25
 *  
Packit 67cb25
 * for M >= N using the QRPT factorization A P = Q R, where A
Packit 67cb25
 * is assumed rank deficient with a given rank.
Packit 67cb25
 */
Packit 67cb25
Packit 67cb25
int
Packit 67cb25
gsl_linalg_QRPT_lssolve2 (const gsl_matrix * QR, const gsl_vector * tau, const gsl_permutation * p,
Packit 67cb25
                          const gsl_vector * b, const size_t rank, gsl_vector * x, gsl_vector * residual)
Packit 67cb25
{
Packit 67cb25
  const size_t M = QR->size1;
Packit 67cb25
  const size_t N = QR->size2;
Packit 67cb25
Packit 67cb25
  if (M < N)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("QR 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 (rank == 0 || rank > N)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("rank must have 0 < rank <= N", 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 R11 = gsl_matrix_const_submatrix (QR, 0, 0, rank, rank);
Packit 67cb25
      gsl_vector_view QTb1 = gsl_vector_subvector(residual, 0, rank);
Packit 67cb25
      gsl_vector_view x1 = gsl_vector_subvector(x, 0, rank);
Packit 67cb25
      size_t i;
Packit 67cb25
Packit 67cb25
      /* compute work = Q^T b */
Packit 67cb25
      gsl_vector_memcpy(residual, b);
Packit 67cb25
      gsl_linalg_QR_QTvec (QR, tau, residual);
Packit 67cb25
Packit 67cb25
      /* solve R_{11} x(1:r) = [Q^T b](1:r) */
Packit 67cb25
      gsl_vector_memcpy(&(x1.vector), &(QTb1.vector));
Packit 67cb25
      gsl_blas_dtrsv (CblasUpper, CblasNoTrans, CblasNonUnit, &(R11.matrix), &(x1.vector));
Packit 67cb25
Packit 67cb25
      /* x(r+1:N) = 0 */
Packit 67cb25
      for (i = rank; i < N; ++i)
Packit 67cb25
        gsl_vector_set(x, i, 0.0);
Packit 67cb25
Packit 67cb25
      /* compute x = P y */
Packit 67cb25
      gsl_permute_vector_inverse (p, x);
Packit 67cb25
Packit 67cb25
      /* compute residual = b - A x = Q (Q^T b - R x) */
Packit 67cb25
      gsl_vector_set_zero(&(QTb1.vector));
Packit 67cb25
      gsl_linalg_QR_Qvec(QR, tau, residual);
Packit 67cb25
Packit 67cb25
      return GSL_SUCCESS;
Packit 67cb25
    }
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
int
Packit 67cb25
gsl_linalg_QRPT_QRsolve (const gsl_matrix * Q, const gsl_matrix * R,
Packit 67cb25
                         const gsl_permutation * p,
Packit 67cb25
                         const gsl_vector * b,
Packit 67cb25
                         gsl_vector * x)
Packit 67cb25
{
Packit 67cb25
  if (Q->size1 != Q->size2 || R->size1 != R->size2)
Packit 67cb25
    {
Packit 67cb25
      return GSL_ENOTSQR;
Packit 67cb25
    }
Packit 67cb25
  else if (Q->size1 != p->size || Q->size1 != R->size1
Packit 67cb25
           || Q->size1 != b->size)
Packit 67cb25
    {
Packit 67cb25
      return GSL_EBADLEN;
Packit 67cb25
    }
Packit 67cb25
  else
Packit 67cb25
    {
Packit 67cb25
      /* compute b' = Q^T b */
Packit 67cb25
Packit 67cb25
      gsl_blas_dgemv (CblasTrans, 1.0, Q, b, 0.0, x);
Packit 67cb25
Packit 67cb25
      /* Solve R x = b', storing x inplace */
Packit 67cb25
Packit 67cb25
      gsl_blas_dtrsv (CblasUpper, CblasNoTrans, CblasNonUnit, R, x);
Packit 67cb25
Packit 67cb25
      /* Apply permutation to solution in place */
Packit 67cb25
Packit 67cb25
      gsl_permute_vector_inverse (p, x);
Packit 67cb25
Packit 67cb25
      return GSL_SUCCESS;
Packit 67cb25
    }
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
int
Packit 67cb25
gsl_linalg_QRPT_Rsolve (const gsl_matrix * QR,
Packit 67cb25
                        const gsl_permutation * p,
Packit 67cb25
                        const gsl_vector * b,
Packit 67cb25
                        gsl_vector * x)
Packit 67cb25
{
Packit 67cb25
  if (QR->size1 != QR->size2)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("QR matrix must be square", GSL_ENOTSQR);
Packit 67cb25
    }
Packit 67cb25
  else if (QR->size1 != b->size)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("matrix size must match b size", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else if (QR->size2 != x->size)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("matrix size must match x size", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else if (p->size != x->size)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("permutation 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 inplace */
Packit 67cb25
Packit 67cb25
      gsl_blas_dtrsv (CblasUpper, CblasNoTrans, CblasNonUnit, QR, x);
Packit 67cb25
Packit 67cb25
      gsl_permute_vector_inverse (p, x);
Packit 67cb25
Packit 67cb25
      return GSL_SUCCESS;
Packit 67cb25
    }
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
Packit 67cb25
int
Packit 67cb25
gsl_linalg_QRPT_Rsvx (const gsl_matrix * QR,
Packit 67cb25
                      const gsl_permutation * p,
Packit 67cb25
                      gsl_vector * x)
Packit 67cb25
{
Packit 67cb25
  if (QR->size1 != QR->size2)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("QR matrix must be square", GSL_ENOTSQR);
Packit 67cb25
    }
Packit 67cb25
  else if (QR->size2 != x->size)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("matrix size must match x size", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else if (p->size != x->size)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("permutation size must match x size", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else
Packit 67cb25
    {
Packit 67cb25
      /* Solve R x = b, storing x inplace */
Packit 67cb25
Packit 67cb25
      gsl_blas_dtrsv (CblasUpper, CblasNoTrans, CblasNonUnit, QR, x);
Packit 67cb25
Packit 67cb25
      gsl_permute_vector_inverse (p, x);
Packit 67cb25
Packit 67cb25
      return GSL_SUCCESS;
Packit 67cb25
    }
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
Packit 67cb25
Packit 67cb25
/* Update a Q R P^T factorisation for A P= Q R ,  A' = A + u v^T,
Packit 67cb25
Packit 67cb25
   Q' R' P^-1 = QR P^-1 + u v^T
Packit 67cb25
   = Q (R + Q^T u v^T P ) P^-1
Packit 67cb25
   = Q (R + w v^T P) P^-1
Packit 67cb25
Packit 67cb25
   where w = Q^T 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
int
Packit 67cb25
gsl_linalg_QRPT_update (gsl_matrix * Q, gsl_matrix * R,
Packit 67cb25
                        const gsl_permutation * p,
Packit 67cb25
                        gsl_vector * w, const gsl_vector * v)
Packit 67cb25
{
Packit 67cb25
  const size_t M = R->size1;
Packit 67cb25
  const size_t N = R->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 if R 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 M if R 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 N if R 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 R,  H = J_1^T ... J^T_(n-1) R
Packit 67cb25
         so that H is upper Hessenberg.  (12.5.2) */
Packit 67cb25
Packit 67cb25
      for (k = M - 1; k > 0; k--)
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_qr (M, N, Q, R, k - 1, k, c, s);
Packit 67cb25
        }
Packit 67cb25
Packit 67cb25
      w0 = gsl_vector_get (w, 0);
Packit 67cb25
Packit 67cb25
      /* Add in w v^T  (Equation 12.5.3) */
Packit 67cb25
Packit 67cb25
      for (j = 0; j < N; j++)
Packit 67cb25
        {
Packit 67cb25
          double r0j = gsl_matrix_get (R, 0, j);
Packit 67cb25
          size_t p_j = gsl_permutation_get (p, j);
Packit 67cb25
          double vj = gsl_vector_get (v, p_j);
Packit 67cb25
          gsl_matrix_set (R, 0, j, r0j + w0 * vj);
Packit 67cb25
        }
Packit 67cb25
Packit 67cb25
      /* Apply Givens transformations R' = 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 (R, k - 1, k - 1);
Packit 67cb25
          double offdiag = gsl_matrix_get (R, k, k - 1);
Packit 67cb25
Packit 67cb25
          gsl_linalg_givens (diag, offdiag, &c, &s);
Packit 67cb25
          apply_givens_qr (M, N, Q, R, k - 1, k, c, s);
Packit 67cb25
Packit 67cb25
          gsl_matrix_set (R, k, k - 1, 0.0);    /* exact zero of G^T */
Packit 67cb25
        }
Packit 67cb25
Packit 67cb25
      return GSL_SUCCESS;
Packit 67cb25
    }
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
/*
Packit 67cb25
gsl_linalg_QRPT_rank()
Packit 67cb25
  Estimate rank of triangular matrix from QRPT decomposition
Packit 67cb25
Packit 67cb25
Inputs: QR  - QRPT decomposed matrix
Packit 67cb25
        tol - tolerance for rank determination; if < 0,
Packit 67cb25
              a default value is used:
Packit 67cb25
              20 * (M + N) * eps(max(|diag(R)|))
Packit 67cb25
Packit 67cb25
Return: rank estimate
Packit 67cb25
*/
Packit 67cb25
Packit 67cb25
size_t
Packit 67cb25
gsl_linalg_QRPT_rank (const gsl_matrix * QR, const double tol)
Packit 67cb25
{
Packit 67cb25
  const size_t M = QR->size1;
Packit 67cb25
  const size_t N = QR->size2;
Packit 67cb25
  gsl_vector_const_view diag = gsl_matrix_const_diagonal(QR);
Packit 67cb25
  double eps;
Packit 67cb25
  size_t i;
Packit 67cb25
  size_t r = 0;
Packit 67cb25
Packit 67cb25
  if (tol < 0.0)
Packit 67cb25
    {
Packit 67cb25
      double min, max, absmax;
Packit 67cb25
      int ee;
Packit 67cb25
Packit 67cb25
      gsl_vector_minmax(&diag.vector, &min, &max;;
Packit 67cb25
      absmax = GSL_MAX(fabs(min), fabs(max));
Packit 67cb25
      ee = (int) (log(absmax) / M_LN2); /* can't use log2 since its not ANSI */
Packit 67cb25
Packit 67cb25
      eps = 20.0 * (M + N) * pow(2.0, (double) ee) * GSL_DBL_EPSILON;
Packit 67cb25
    }
Packit 67cb25
  else
Packit 67cb25
    eps = tol;
Packit 67cb25
Packit 67cb25
  /* count number of diagonal elements with |di| > eps */
Packit 67cb25
  for (i = 0; i < GSL_MIN(M, N); ++i)
Packit 67cb25
    {
Packit 67cb25
      double di = gsl_vector_get(&diag.vector, i);
Packit 67cb25
      if (fabs(di) > eps)
Packit 67cb25
        ++r;
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
  return r;
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
int
Packit 67cb25
gsl_linalg_QRPT_rcond(const gsl_matrix * QR, double * rcond, gsl_vector * work)
Packit 67cb25
{
Packit 67cb25
  const size_t M = QR->size1;
Packit 67cb25
  const size_t N = QR->size2;
Packit 67cb25
Packit 67cb25
  if (M < N)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("M must be >= N", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else if (work->size != 3 * N)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("work vector must have length 3*N", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else
Packit 67cb25
    {
Packit 67cb25
      gsl_matrix_const_view R = gsl_matrix_const_submatrix (QR, 0, 0, N, N);
Packit 67cb25
      int status;
Packit 67cb25
Packit 67cb25
      status = gsl_linalg_tri_upper_rcond(&R.matrix, rcond, work);
Packit 67cb25
Packit 67cb25
      return status;
Packit 67cb25
    }
Packit 67cb25
}