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