00001 /* 00002 ------------------------------------------------------------------- 00003 00004 Copyright (C) 2006, 2007, 2008, 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 /** \file qr_base.h 00024 \brief File for QR decomposition and associated solver 00025 */ 00026 00027 #include <o2scl/householder.h> 00028 #include <o2scl/givens.h> 00029 00030 #ifdef DOXYGENP 00031 namespace o2scl_linalg { 00032 #endif 00033 00034 /** 00035 \brief Compute the QR decomposition of matrix \c A 00036 */ 00037 template<class mat_t, class vec_t> 00038 int QR_decomp(size_t M, size_t N, mat_t &A, vec_t &tau) { 00039 00040 size_t i; 00041 for(i=0;i<GSL_MIN(M,N);i++) { 00042 O2SCL_IX(tau,i)=householder_transform_subcol(A,i,i,M-i); 00043 if (i+1<N) { 00044 householder_hm_sub(A,i,i+1,M-i,N-(i+1),A,i,i,O2SCL_IX(tau,i)); 00045 } 00046 } 00047 return 0; 00048 } 00049 00050 /** 00051 \brief Solve the system A x = b using the QR factorization 00052 */ 00053 template<class mat_t, class vec_t> 00054 int QR_solve(size_t N, const mat_t &QR, const vec_t &tau, 00055 const vec_t &b, vec_t &x) { 00056 o2scl::vector_copy(N,b,x); 00057 QR_svx(N,N,QR,tau,x); 00058 return 0; 00059 } 00060 00061 /** 00062 \brief Solve the system A x = b in place using the QR factorization 00063 */ 00064 template<class mat_t, class vec_t> 00065 int QR_svx(size_t M, size_t N, const mat_t &QR, const vec_t &tau, 00066 vec_t &x) { 00067 QR_QTvec(M,N,QR,tau,x); 00068 o2scl_cblas::dtrsv(o2scl_cblas::O2cblasRowMajor, 00069 o2scl_cblas::O2cblasUpper, 00070 o2scl_cblas::O2cblasNoTrans, 00071 o2scl_cblas::O2cblasNonUnit, 00072 M,N,QR,x); 00073 return 0; 00074 } 00075 00076 /** 00077 \brief Form the product Q^T v from a QR factorized matrix 00078 */ 00079 template<class mat_t, class vec_t> 00080 int QR_QTvec(const size_t M, const size_t N, 00081 const mat_t &QR, const vec_t &tau, vec_t &v) { 00082 00083 size_t i; 00084 00085 // compute Q^T v 00086 for (i = 0; i < GSL_MIN(M,N); i++) { 00087 householder_hv_sub(QR,v,O2SCL_IX(tau,i),i,M); 00088 } 00089 return o2scl::gsl_success; 00090 } 00091 00092 /** 00093 \brief Unpack the QR matrix to the individual Q and R components 00094 */ 00095 template<class mat1_t, class mat2_t, class mat3_t, class vec_t> 00096 int QR_unpack(const size_t M, const size_t N, 00097 const mat1_t &QR, const vec_t &tau, mat2_t &Q, 00098 mat3_t &R) { 00099 00100 size_t i, j; 00101 00102 /* Initialize Q to the identity */ 00103 00104 for(i=0;i<M;i++) { 00105 for(j=0;j<M;j++) { 00106 if (i==j) O2SCL_IX2(Q,i,j)=1.0; 00107 else O2SCL_IX2(Q,i,j)=0.0; 00108 } 00109 } 00110 00111 for (i = GSL_MIN (M, N); i > 0 && i--;) { 00112 householder_hm_sub2(M,i,O2SCL_IX(tau,i),QR,Q); 00113 } 00114 00115 /* Form the right triangular matrix R from a packed QR matrix */ 00116 00117 for (i = 0; i < M; i++) { 00118 for (j = 0; j < i && j < N; j++) { 00119 O2SCL_IX2(R,i,j)=0.0; 00120 } 00121 for (j = i; j < N; j++) { 00122 O2SCL_IX2(R,i,j)=O2SCL_IX2(QR,i,j); 00123 } 00124 } 00125 00126 return o2scl::gsl_success; 00127 } 00128 00129 /** 00130 \brief Update a QR factorisation for A= Q R , A' = A + u v^T, 00131 00132 \c M and \c N are the number of rows and columns of \c R. 00133 00134 \verbatim 00135 * Q' R' = QR + u v^T 00136 * = Q (R + Q^T u v^T) 00137 * = Q (R + w v^T) 00138 * 00139 * where w = Q^T u. 00140 * 00141 * Algorithm from Golub and Van Loan, "Matrix Computations", Section 00142 * 12.5 (Updating Matrix Factorizations, Rank-One Changes) 00143 \endverbatim 00144 */ 00145 00146 template<class mat1_t, class mat2_t, class vec1_t, class vec2_t> 00147 int QR_update(size_t M, size_t N, mat1_t &Q, mat2_t &R, 00148 vec1_t &w, vec2_t &v) { 00149 00150 size_t j, k; 00151 double w0; 00152 00153 /* 00154 Apply Given's rotations to reduce w to (|w|, 0, 0, ... , 0) 00155 00156 J_1^T .... J_(n-1)^T w = +/- |w| e_1 00157 00158 simultaneously applied to R, H = J_1^T ... J^T_(n-1) R 00159 so that H is upper Hessenberg. (12.5.2) 00160 */ 00161 00162 for (k = M - 1; k > 0; k--) { 00163 00164 double c, s; 00165 double wk = O2SCL_IX(w,k); 00166 double wkm1 = O2SCL_IX(w,k-1); 00167 00168 create_givens (wkm1, wk, c, s); 00169 apply_givens_vec (w, k - 1, k, c, s); 00170 apply_givens_qr (M, N, Q, R, k - 1, k, c, s); 00171 } 00172 00173 w0 = O2SCL_IX(w,0); 00174 00175 /* Add in w v^T (Equation 12.5.3) */ 00176 00177 for (j = 0; j < N; j++) { 00178 double r0j = O2SCL_IX2(R,0,j); 00179 double vj = O2SCL_IX(v,j); 00180 O2SCL_IX2(R,0,j)=r0j+w0*vj; 00181 } 00182 00183 /* Apply Givens transformations R' = G_(n-1)^T ... G_1^T H 00184 Equation 12.5.4 */ 00185 for (k = 1; k < GSL_MIN(M,N+1); k++) { 00186 00187 double c, s; 00188 double diag = O2SCL_IX2(R,k-1,k-1); 00189 double offdiag = O2SCL_IX2(R,k,k-1); 00190 00191 create_givens (diag, offdiag, c, s); 00192 apply_givens_qr (M, N, Q, R, k - 1, k, c, s); 00193 00194 O2SCL_IX2(R,k,j-1)=0.0; 00195 } 00196 00197 return o2scl::gsl_success; 00198 00199 } 00200 00201 #ifdef DOXYGENP 00202 } 00203 #endif
Documentation generated with Doxygen and provided under the GNU Free Documentation License. See License Information for details.
Project hosting provided by
,
O2scl Sourceforge Project Page