Object-oriented Scientific Computing Library: Version 0.910
qr_base.h
Go to the documentation of this file.
00001 /*
00002   -------------------------------------------------------------------
00003 
00004   Copyright (C) 2006-2012, Andrew W. Steiner
00005 
00006   This file is part of O2scl.
00007 
00008   O2scl is free software; you can redistribute it and/or modify
00009   it under the terms of the GNU General Public License as published by
00010   the Free Software Foundation; either version 3 of the License, or
00011   (at your option) any later version.
00012 
00013   O2scl is distributed in the hope that it will be useful,
00014   but WITHOUT ANY WARRANTY; without even the implied warranty of
00015   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00016   GNU General Public License for more details.
00017 
00018   You should have received a copy of the GNU General Public License
00019   along with O2scl. If not, see <http://www.gnu.org/licenses/>.
00020 
00021   -------------------------------------------------------------------
00022 */
00023 /* linalg/qr.c
00024  * 
00025  * Copyright (C) 1996, 1997, 1998, 1999, 2000, 2007 Gerard Jungman, 
00026  * Brian Gough
00027  * 
00028  * This program is free software; you can redistribute it and/or modify
00029  * it under the terms of the GNU General Public License as published by
00030  * the Free Software Foundation; either version 3 of the License, or (at
00031  * your option) any later version.
00032  * 
00033  * This program is distributed in the hope that it will be useful, but
00034  * WITHOUT ANY WARRANTY; without even the implied warranty of
00035  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00036  * General Public License for more details.
00037  * 
00038  * You should have received a copy of the GNU General Public License
00039  * along with this program; if not, write to the Free Software
00040  * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 
00041  * 02110-1301, USA.
00042  */
00043 /** \file qr_base.h
00044     \brief File for QR decomposition and associated solver
00045 */
00046 
00047 #include <o2scl/householder.h>
00048 #include <o2scl/givens.h>
00049 
00050 #ifdef DOXYGENP
00051 namespace o2scl_linalg {
00052 #endif
00053 
00054   /** \brief Compute the QR decomposition of matrix \c A
00055   */
00056   template<class mat_t, class vec_t>
00057     void QR_decomp(size_t M, size_t N, mat_t &A, vec_t &tau) {
00058     
00059     size_t i;
00060     for(i=0;i<GSL_MIN(M,N);i++) {
00061       O2SCL_IX(tau,i)=householder_transform_subcol(A,i,i,M);
00062       if (i+1<N) {
00063         householder_hm_sub(A,i,i+1,M,N,A,i,i,O2SCL_IX(tau,i));
00064       }
00065     }
00066     return;
00067   }
00068   
00069   /** \brief Solve the system A x = b using the QR factorization
00070   */
00071   template<class mat_t, class vec_t, class vec2_t, class vec3_t>
00072     void QR_solve(size_t N, const mat_t &QR, const vec_t &tau, 
00073                  const vec2_t &b, vec3_t &x) {
00074     o2scl::vector_copy(N,b,x);
00075     QR_svx(N,N,QR,tau,x);
00076     return;
00077   }
00078   
00079   /** \brief Solve the system A x = b in place using the QR factorization
00080   */
00081   template<class mat_t, class vec_t, class vec2_t>
00082     void QR_svx(size_t M, size_t N, const mat_t &QR, const vec_t &tau, 
00083                vec2_t &x) {
00084     QR_QTvec(M,N,QR,tau,x);
00085     o2scl_cblas::dtrsv(o2scl_cblas::o2cblas_RowMajor,
00086                        o2scl_cblas::o2cblas_Upper,
00087                        o2scl_cblas::o2cblas_NoTrans,
00088                        o2scl_cblas::o2cblas_NonUnit,
00089                        M,N,QR,x);
00090     return;
00091   }
00092 
00093   /** \brief Form the product Q^T v from a QR factorized matrix
00094   */
00095   template<class mat_t, class vec_t, class vec2_t>
00096     void QR_QTvec(const size_t M, const size_t N, 
00097                  const mat_t &QR, const vec_t &tau, vec2_t &v) {
00098     
00099     size_t i;
00100     
00101     // compute Q^T v 
00102     for (i = 0; i < GSL_MIN(M,N); i++) {
00103       householder_hv_subcol(QR,v,O2SCL_IX(tau,i),i,M);
00104     }
00105     return;
00106   }
00107 
00108   /** \brief Unpack the QR matrix to the individual Q and R components
00109   */
00110   template<class mat1_t, class mat2_t, class mat3_t, class vec_t>
00111     void QR_unpack(const size_t M, const size_t N,
00112                   const mat1_t &QR, const vec_t &tau, mat2_t &Q, 
00113                   mat3_t &R) {
00114     
00115     size_t i, j;
00116     
00117     /* Initialize Q to the identity */
00118     
00119     for(i=0;i<M;i++) {
00120       for(j=0;j<M;j++) {
00121         if (i==j) O2SCL_IX2(Q,i,j)=1.0;
00122         else O2SCL_IX2(Q,i,j)=0.0;
00123       }
00124     }
00125     
00126     for (i = GSL_MIN (M, N); i-- > 0;) {
00127       householder_hm_sub2(M,i,O2SCL_IX(tau,i),QR,Q);
00128     }
00129     
00130     /*  Form the right triangular matrix R from a packed QR matrix */
00131     
00132     for (i = 0; i < M; i++) {
00133       for (j = 0; j < i && j < N; j++) {
00134         O2SCL_IX2(R,i,j)=0.0;
00135       }
00136       for (j = i; j < N; j++) {
00137         O2SCL_IX2(R,i,j)=O2SCL_IX2(QR,i,j);
00138       }
00139     }
00140     
00141     return;
00142   }
00143 
00144   /** \brief Update a QR factorisation for A= Q R,  A' = A + u v^T,
00145 
00146       The parameters \c M and \c N are the number of rows and columns 
00147       of the matrix \c R.
00148 
00149       \verbatim
00150       * Q' R' = QR + u v^T
00151       *       = Q (R + Q^T u v^T)
00152       *       = Q (R + w v^T)
00153       *
00154       * where w = Q^T u.
00155       *
00156       * Algorithm from Golub and Van Loan, "Matrix Computations", Section
00157       * 12.5 (Updating Matrix Factorizations, Rank-One Changes)
00158       \endverbatim
00159   */
00160   template<class mat1_t, class mat2_t, class vec1_t, class vec2_t> 
00161     void QR_update(size_t M, size_t N, mat1_t &Q, mat2_t &R, 
00162                   vec1_t &w, vec2_t &v) {
00163     
00164     size_t j;
00165     // Integer to avoid problems with decreasing loop below
00166     int k;
00167     double w0;
00168          
00169     /* 
00170        Apply Given's rotations to reduce w to (|w|, 0, 0, ... , 0)
00171        
00172        J_1^T .... J_(n-1)^T w = +/- |w| e_1
00173        
00174        simultaneously applied to R,  H = J_1^T ... J^T_(n-1) R
00175        so that H is upper Hessenberg.  (12.5.2) 
00176     */
00177          
00178     // Loop from k = M-1 to 1
00179     for (k = M - 1; k > 0; k--) {
00180 
00181       double c, s;
00182       double wk = O2SCL_IX(w,k);
00183       double wkm1 = O2SCL_IX(w,k-1);
00184              
00185       create_givens (wkm1, wk, c, s);
00186       apply_givens_vec (w, k - 1, k, c, s);
00187       apply_givens_qr (M, N, Q, R, k - 1, k, c, s);
00188     }
00189     
00190     w0 = O2SCL_IX(w,0);
00191          
00192     /* Add in w v^T  (Equation 12.5.3) */
00193          
00194     for (j = 0; j < N; j++) {
00195       double r0j = O2SCL_IX2(R,0,j);
00196       double vj = O2SCL_IX(v,j);
00197       O2SCL_IX2(R,0,j)=r0j+w0*vj;
00198     }
00199          
00200     /* Apply Givens transformations R' = G_(n-1)^T ... G_1^T H
00201        Equation 12.5.4 */
00202     for (k = 1; k < ((int)(GSL_MIN(M,N+1))); k++) {
00203 
00204       double c, s;
00205       double diag = O2SCL_IX2(R,k-1,k-1);
00206       double offdiag = O2SCL_IX2(R,k,k-1);
00207              
00208       create_givens(diag, offdiag, c, s);
00209       apply_givens_qr(M, N, Q, R, k - 1, k, c, s);
00210     
00211       O2SCL_IX2(R,k,k-1)=0.0;
00212     }
00213          
00214     return;
00215   }
00216   
00217 #ifdef DOXYGENP
00218 }
00219 #endif
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines

Documentation generated with Doxygen. Provided under the GNU Free Documentation License (see License Information).

Get Object-oriented Scientific Computing
Lib at SourceForge.net. Fast, secure and Free Open Source software
downloads.