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 lu_base.h 00024 \brief File for LU decomposition and associated solver 00025 */ 00026 00027 #ifdef DOXYGENP 00028 namespace o2scl_linalg { 00029 #endif 00030 00031 /** 00032 \brief Compute the LU decomposition of the matrix \c A 00033 00034 On output the diagonal and upper triangular part of the input 00035 matrix A contain the matrix U. The lower triangular part of the 00036 input matrix (excluding the diagonal) contains L. The diagonal 00037 elements of L are unity, and are not stored. 00038 00039 The permutation matrix P is encoded in the permutation p. The 00040 j-th column of the matrix P is given by the k-th column of the 00041 identity matrix, where k = p_j the j-th element of the 00042 permutation vector. The sign of the permutation is given by 00043 signum. It has the value (-1)^n, where n is the number of 00044 interchanges in the permutation. 00045 00046 The algorithm used in the decomposition is Gaussian Elimination 00047 with partial pivoting (Golub & Van Loan, Matrix Computations, 00048 Algorithm 3.4.1). 00049 */ 00050 template<class mat_t> 00051 int LU_decomp(const size_t N, mat_t &A, o2scl::permutation &p, 00052 int &signum) { 00053 00054 size_t i, j, k; 00055 00056 signum=1; 00057 p.init(); 00058 00059 for (j = 0; j < N - 1; j++) { 00060 00061 /* Find maximum in the j-th column */ 00062 double ajj, max = fabs(O2SCL_IX2(A,j,j)); 00063 size_t i_pivot = j; 00064 00065 for (i = j + 1; i < N; i++) { 00066 double aij = fabs (O2SCL_IX2(A,i,j)); 00067 00068 if (aij > max) { 00069 max = aij; 00070 i_pivot = i; 00071 } 00072 } 00073 00074 if (i_pivot != j) { 00075 00076 // Swap rows j and i_pivot 00077 double temp; 00078 for (k=0;k<N;k++) { 00079 temp=O2SCL_IX2(A,j,k); 00080 O2SCL_IX2(A,j,k)=O2SCL_IX2(A,i_pivot,k); 00081 O2SCL_IX2(A,i_pivot,k)=temp; 00082 } 00083 p.swap(j,i_pivot); 00084 signum=-signum; 00085 } 00086 00087 ajj = O2SCL_IX2(A,j,j); 00088 00089 if (ajj != 0.0) { 00090 for (i = j + 1; i < N; i++) { 00091 double aij = O2SCL_IX2(A,i,j) / ajj; 00092 O2SCL_IX2(A,i,j)=aij; 00093 for (k = j + 1; k < N; k++) { 00094 double aik = O2SCL_IX2(A,i,k); 00095 double ajk = O2SCL_IX2(A,j,k); 00096 O2SCL_IX2(A,i,k)=aik - aij * ajk; 00097 } 00098 } 00099 } 00100 } 00101 00102 return o2scl::gsl_success; 00103 } 00104 00105 /** 00106 \brief Solve a linear system after LU decomposition 00107 00108 This function solve the square system A x = b using the LU 00109 decomposition of A into (LU, p) given by gsl_linalg_LU_decomp or 00110 gsl_linalg_complex_LU_decomp. 00111 */ 00112 template<class mat_t, class vec_t> 00113 int LU_solve(const size_t N, const mat_t &LU, const o2scl::permutation &p, 00114 const vec_t &b, vec_t &x) { 00115 00116 /* Copy x <- b */ 00117 o2scl::vector_copy(N,b,x); 00118 00119 /* Solve for x */ 00120 LU_svx(N,LU,p,x); 00121 00122 return o2scl::gsl_success; 00123 } 00124 00125 /** 00126 \brief Solve a linear system after LU decomposition in place 00127 00128 These functions solve the square system A x = b in-place using 00129 the LU decomposition of A into (LU,p). On input x should contain 00130 the right-hand side b, which is replaced by the solution on 00131 output. 00132 */ 00133 template<class mat_t, class vec_t> 00134 int LU_svx(const size_t N, const mat_t &LU, 00135 const o2scl::permutation &p, vec_t &x) { 00136 00137 /* Apply permutation to RHS */ 00138 p.apply(x); 00139 00140 /* Solve for c using forward-substitution, L c = P b */ 00141 o2scl_cblas::dtrsv(o2scl_cblas::O2cblasRowMajor,o2scl_cblas::O2cblasLower, 00142 o2scl_cblas::O2cblasNoTrans,o2scl_cblas::O2cblasUnit, 00143 N,N,LU,x); 00144 00145 /* Perform back-substitution, U x = c */ 00146 o2scl_cblas::dtrsv(o2scl_cblas::O2cblasRowMajor,o2scl_cblas::O2cblasUpper, 00147 o2scl_cblas::O2cblasNoTrans,o2scl_cblas::O2cblasNonUnit, 00148 N,N,LU,x); 00149 00150 return o2scl::gsl_success; 00151 } 00152 00153 /** 00154 \brief Refine the solution of a linear system 00155 00156 These functions apply an iterative improvement to x, the solution 00157 of A x = b, using the LU decomposition of A into (LU,p). The 00158 initial residual r = A x - b is also computed and stored in 00159 residual. 00160 */ 00161 template<class mat_t, class vec_t> 00162 int LU_refine(const size_t N, const mat_t &A, const mat_t &LU, 00163 const o2scl::permutation &p, const vec_t &b, vec_t &x, 00164 vec_t &residual) { 00165 00166 /* Compute residual, residual = (A * x - b) */ 00167 o2scl::vector_copy(N,b,residual); 00168 o2scl_cblas::dgemv(o2scl_cblas::O2cblasRowMajor, 00169 o2scl_cblas::O2cblasNoTrans, 00170 N,N,1.0,A,x,-1.0,residual); 00171 00172 /* Find correction, delta = - (A^-1) * residual, and apply it */ 00173 00174 LU_svx(N,LU,p,residual); 00175 o2scl_cblas::daxpy(5,-1.0,residual,x); 00176 00177 return o2scl::gsl_success; 00178 } 00179 00180 /** 00181 \brief Compute the inverse of a matrix from its LU decomposition 00182 00183 These functions compute the inverse of a matrix A from its LU 00184 decomposition (LU,p), storing the result in the matrix inverse. 00185 The inverse is computed by solving the system A x = b for each 00186 column of the identity matrix. It is preferable to avoid direct 00187 use of the inverse whenever possible, as the linear solver 00188 functions can obtain the same result more efficiently and 00189 reliably (consult any introductory textbook on numerical linear 00190 algebra for details). 00191 00192 \future could rewrite to avoid mat_col_t 00193 */ 00194 template<class mat_t, class mat_col_t> 00195 int LU_invert(const size_t N, const mat_t &LU, 00196 const o2scl::permutation &p, mat_t &inverse) { 00197 size_t i; 00198 00199 int status=o2scl::gsl_success; 00200 00201 // Set matrix 'inverse' to the identity 00202 for(i=0;i<N;i++) { 00203 for(size_t j=0;j<N;j++) { 00204 if (i==j) O2SCL_IX2(inverse,i,j)=1.0; 00205 else O2SCL_IX2(inverse,i,j)=0.0; 00206 } 00207 } 00208 00209 for (i = 0; i < N; i++) { 00210 mat_col_t c(inverse,i); 00211 int status_i=LU_svx(N,LU,p,c); 00212 00213 if (status_i) { 00214 status = status_i; 00215 } 00216 } 00217 00218 return status; 00219 } 00220 00221 /** 00222 \brief Compute the determinant of a matrix from its LU decomposition 00223 00224 These functions compute the determinant of a matrix A from its LU 00225 decomposition, LU. The determinant is computed as the product of 00226 the diagonal elements of U and the sign of the row permutation 00227 signum. 00228 */ 00229 template<class mat_t> 00230 double LU_det(const size_t N, const mat_t &LU, int signum) { 00231 00232 size_t i; 00233 00234 double det=(double)signum; 00235 00236 for (i=0;i<N;i++) { 00237 det*=O2SCL_IX2(LU,i,i); 00238 } 00239 00240 return det; 00241 } 00242 00243 /** 00244 \brief Compute the logarithm of the absolute value of the 00245 determinant of a matrix from its LU decomposition 00246 00247 These functions compute the logarithm of the absolute value of 00248 the determinant of a matrix A, \f$ \ln|\det(A)| \f$, from its LU 00249 decomposition, LU. This function may be useful if the direct 00250 computation of the determinant would overflow or underflow. 00251 */ 00252 template<class mat_t> 00253 double LU_lndet(const size_t N, const mat_t &LU) { 00254 00255 size_t i; 00256 double lndet = 0.0; 00257 00258 for (i = 0; i < N; i++) { 00259 lndet+=log(fabs(O2SCL_IX2(LU,i,i))); 00260 } 00261 00262 return lndet; 00263 } 00264 00265 /** 00266 \brief Compute the sign of the 00267 determinant of a matrix from its LU decomposition 00268 00269 These functions compute the sign or phase factor of the 00270 determinant of a matrix A, \f$ \det(A)/|\det(A)| \f$, from its 00271 LU decomposition, LU. 00272 */ 00273 template<class mat_t> 00274 int LU_sgndet(const size_t N, const mat_t &LU, int signum) { 00275 00276 size_t i; 00277 int s=signum; 00278 00279 for (i=0;i<N;i++) { 00280 double u=O2SCL_IX2(LU,i,i); 00281 if (u<0) { 00282 s*=-1; 00283 } else if (u == 0) { 00284 s=0; 00285 break; 00286 } 00287 } 00288 00289 return s; 00290 } 00291 00292 #ifdef DOXYGENP 00293 } 00294 #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