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/lu.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 lu_base.h 00043 \brief File for LU decomposition and associated solver 00044 */ 00045 00046 #ifdef DOXYGENP 00047 namespace o2scl_linalg { 00048 #endif 00049 00050 /** \brief Return 1 if at least one of the elements in the 00051 diagonal is zero 00052 */ 00053 template<class mat_t> 00054 int diagonal_has_zero(const size_t N, mat_t &A) { 00055 00056 for(size_t i=0;i<N;i++) { 00057 if (O2SCL_IX2(A,i,i)==0.0) return 1; 00058 } 00059 return 0; 00060 } 00061 00062 /** 00063 \brief Compute the LU decomposition of the matrix \c A 00064 00065 On output the diagonal and upper triangular part of the input 00066 matrix A contain the matrix U. The lower triangular part of the 00067 input matrix (excluding the diagonal) contains L. The diagonal 00068 elements of L are unity, and are not stored. 00069 00070 The permutation matrix P is encoded in the permutation p. The 00071 j-th column of the matrix P is given by the k-th column of the 00072 identity matrix, where k = p_j the j-th element of the 00073 permutation vector. The sign of the permutation is given by 00074 signum. It has the value (-1)^n, where n is the number of 00075 interchanges in the permutation. 00076 00077 The algorithm used in the decomposition is Gaussian Elimination 00078 with partial pivoting (Golub & Van Loan, Matrix Computations, 00079 Algorithm 3.4.1). 00080 00081 \future The "swap rows j and i_pivot" section could probably 00082 be made more efficient using a "matrix_row"-like object 00083 as done in GSL. (7/16/09 - I've tried this, and it doesn't 00084 seem to improve the speed significantly.) 00085 */ 00086 template<class mat_t> 00087 int LU_decomp(const size_t N, mat_t &A, o2scl::permutation &p, 00088 int &signum) { 00089 00090 size_t i, j, k; 00091 00092 signum=1; 00093 p.init(); 00094 00095 for (j = 0; j < N - 1; j++) { 00096 00097 /* Find maximum in the j-th column */ 00098 double ajj, max = fabs(O2SCL_IX2(A,j,j)); 00099 size_t i_pivot = j; 00100 00101 for (i = j + 1; i < N; i++) { 00102 double aij = fabs (O2SCL_IX2(A,i,j)); 00103 00104 if (aij > max) { 00105 max = aij; 00106 i_pivot = i; 00107 } 00108 } 00109 00110 if (i_pivot != j) { 00111 00112 // Swap rows j and i_pivot 00113 double temp; 00114 for (k=0;k<N;k++) { 00115 temp=O2SCL_IX2(A,j,k); 00116 O2SCL_IX2(A,j,k)=O2SCL_IX2(A,i_pivot,k); 00117 O2SCL_IX2(A,i_pivot,k)=temp; 00118 } 00119 p.swap(j,i_pivot); 00120 signum=-signum; 00121 } 00122 00123 ajj = O2SCL_IX2(A,j,j); 00124 00125 if (ajj != 0.0) { 00126 for (i = j + 1; i < N; i++) { 00127 double aij = O2SCL_IX2(A,i,j) / ajj; 00128 O2SCL_IX2(A,i,j)=aij; 00129 for (k = j + 1; k < N; k++) { 00130 double aik = O2SCL_IX2(A,i,k); 00131 double ajk = O2SCL_IX2(A,j,k); 00132 O2SCL_IX2(A,i,k)=aik - aij * ajk; 00133 } 00134 } 00135 } 00136 } 00137 00138 return o2scl::gsl_success; 00139 } 00140 00141 /* 00142 Was testing this out, but it doesn't seem much faster 00143 */ 00144 template<class mat_t, class mat_row_t> 00145 int LU_decomp_alt(const size_t N, mat_t &A, o2scl::permutation &p, 00146 int &signum) { 00147 00148 size_t i, j, k; 00149 00150 signum=1; 00151 p.init(); 00152 00153 for (j = 0; j < N - 1; j++) { 00154 00155 /* Find maximum in the j-th column */ 00156 double ajj, max = fabs(O2SCL_IX2(A,j,j)); 00157 size_t i_pivot = j; 00158 00159 for (i = j + 1; i < N; i++) { 00160 double aij = fabs (O2SCL_IX2(A,i,j)); 00161 00162 if (aij > max) { 00163 max = aij; 00164 i_pivot = i; 00165 } 00166 } 00167 00168 if (i_pivot != j) { 00169 00170 // Swap rows j and i_pivot 00171 double temp; 00172 mat_row_t r1(A,j); 00173 mat_row_t r2(A,i_pivot); 00174 for (k=0;k<N;k++) { 00175 temp=O2SCL_IX(r1,k); 00176 O2SCL_IX(r1,k)=O2SCL_IX(r2,k); 00177 O2SCL_IX(r2,k)=temp; 00178 } 00179 p.swap(j,i_pivot); 00180 signum=-signum; 00181 } 00182 00183 ajj = O2SCL_IX2(A,j,j); 00184 00185 if (ajj != 0.0) { 00186 for (i = j + 1; i < N; i++) { 00187 double aij = O2SCL_IX2(A,i,j) / ajj; 00188 O2SCL_IX2(A,i,j)=aij; 00189 for (k = j + 1; k < N; k++) { 00190 double aik = O2SCL_IX2(A,i,k); 00191 double ajk = O2SCL_IX2(A,j,k); 00192 O2SCL_IX2(A,i,k)=aik - aij * ajk; 00193 } 00194 } 00195 } 00196 } 00197 00198 return o2scl::gsl_success; 00199 } 00200 00201 /** 00202 \brief Solve a linear system after LU decomposition 00203 00204 This function solve the square system A x = b using the LU 00205 decomposition of A into (LU, p) given by gsl_linalg_LU_decomp or 00206 gsl_linalg_complex_LU_decomp. 00207 */ 00208 template<class mat_t, class vec_t, class vec2_t> 00209 int LU_solve(const size_t N, const mat_t &LU, const o2scl::permutation &p, 00210 const vec_t &b, vec2_t &x) { 00211 00212 if (diagonal_has_zero(N,LU)) { 00213 O2SCL_ERR_RET("LU matrix is singlar in LU_solve().", 00214 o2scl::gsl_edom); 00215 } 00216 00217 /* Copy x <- b */ 00218 o2scl::vector_copy(N,b,x); 00219 00220 /* Solve for x */ 00221 return LU_svx(N,LU,p,x); 00222 } 00223 00224 /** 00225 \brief Solve a linear system after LU decomposition in place 00226 00227 These functions solve the square system A x = b in-place using 00228 the LU decomposition of A into (LU,p). On input x should contain 00229 the right-hand side b, which is replaced by the solution on 00230 output. 00231 */ 00232 template<class mat_t, class vec_t> 00233 int LU_svx(const size_t N, const mat_t &LU, 00234 const o2scl::permutation &p, vec_t &x) { 00235 00236 if (diagonal_has_zero(N,LU)) { 00237 O2SCL_ERR_RET("LU matrix is singlar in LU_svx().", 00238 o2scl::gsl_edom); 00239 } 00240 00241 /* Apply permutation to RHS */ 00242 p.apply(x); 00243 00244 /* Solve for c using forward-substitution, L c = P b */ 00245 o2scl_cblas::dtrsv(o2scl_cblas::O2cblasRowMajor,o2scl_cblas::O2cblasLower, 00246 o2scl_cblas::O2cblasNoTrans,o2scl_cblas::O2cblasUnit, 00247 N,N,LU,x); 00248 00249 /* Perform back-substitution, U x = c */ 00250 o2scl_cblas::dtrsv(o2scl_cblas::O2cblasRowMajor,o2scl_cblas::O2cblasUpper, 00251 o2scl_cblas::O2cblasNoTrans,o2scl_cblas::O2cblasNonUnit, 00252 N,N,LU,x); 00253 00254 return o2scl::gsl_success; 00255 } 00256 00257 /** 00258 \brief Refine the solution of a linear system 00259 00260 These functions apply an iterative improvement to x, the solution 00261 of A x = b, using the LU decomposition of A into (LU,p). The 00262 initial residual r = A x - b is also computed and stored in 00263 residual. 00264 */ 00265 template<class mat_t, class mat2_t, class vec_t, class vec2_t, class vec3_t> 00266 int LU_refine(const size_t N, const mat_t &A, const mat2_t &LU, 00267 const o2scl::permutation &p, const vec_t &b, vec2_t &x, 00268 vec3_t &residual) { 00269 00270 if (diagonal_has_zero(N,LU)) { 00271 O2SCL_ERR_RET("LU matrix is singlar in LU_refine().", 00272 o2scl::gsl_edom); 00273 } 00274 00275 /* Compute residual, residual = (A * x - b) */ 00276 o2scl::vector_copy(N,b,residual); 00277 o2scl_cblas::dgemv(o2scl_cblas::O2cblasRowMajor, 00278 o2scl_cblas::O2cblasNoTrans, 00279 N,N,1.0,A,x,-1.0,residual); 00280 00281 /* Find correction, delta = - (A^-1) * residual, and apply it */ 00282 00283 int status=LU_svx(N,LU,p,residual); 00284 o2scl_cblas::daxpy(5,-1.0,residual,x); 00285 00286 return status; 00287 } 00288 00289 /** 00290 \brief Compute the inverse of a matrix from its LU decomposition 00291 00292 These functions compute the inverse of a matrix A from its LU 00293 decomposition (LU,p), storing the result in the matrix inverse. 00294 The inverse is computed by solving the system A x = b for each 00295 column of the identity matrix. It is preferable to avoid direct 00296 use of the inverse whenever possible, as the linear solver 00297 functions can obtain the same result more efficiently and 00298 reliably. 00299 00300 \future Could rewrite to avoid mat_col_t, (9/16/09 - However, 00301 the function may be faster if mat_col_t is left in, so it's 00302 unclear what's best.) 00303 */ 00304 template<class mat_t, class mat2_t, class mat_col_t> 00305 int LU_invert(const size_t N, const mat_t &LU, 00306 const o2scl::permutation &p, mat2_t &inverse) { 00307 00308 if (diagonal_has_zero(N,LU)) { 00309 O2SCL_ERR_RET("LU matrix is singlar in LU_invert().", 00310 o2scl::gsl_edom); 00311 } 00312 00313 size_t i; 00314 00315 int status=o2scl::gsl_success; 00316 00317 // Set matrix 'inverse' to the identity 00318 for(i=0;i<N;i++) { 00319 for(size_t j=0;j<N;j++) { 00320 if (i==j) O2SCL_IX2(inverse,i,j)=1.0; 00321 else O2SCL_IX2(inverse,i,j)=0.0; 00322 } 00323 } 00324 00325 for (i = 0; i < N; i++) { 00326 mat_col_t c(inverse,i); 00327 int status_i=LU_svx(N,LU,p,c); 00328 00329 if (status_i) { 00330 status = status_i; 00331 } 00332 } 00333 00334 return status; 00335 } 00336 00337 /** 00338 \brief Compute the determinant of a matrix from its LU decomposition 00339 00340 Compute the determinant of a matrix A from its LU decomposition, 00341 LU. The determinant is computed as the product of the diagonal 00342 elements of U and the sign of the row permutation signum. 00343 */ 00344 template<class mat_t> 00345 double LU_det(const size_t N, const mat_t &LU, int signum) { 00346 00347 size_t i; 00348 00349 double det=(double)signum; 00350 00351 for (i=0;i<N;i++) { 00352 det*=O2SCL_IX2(LU,i,i); 00353 } 00354 00355 return det; 00356 } 00357 00358 /** 00359 \brief Compute the logarithm of the absolute value of the 00360 determinant of a matrix from its LU decomposition 00361 00362 Compute the logarithm of the absolute value of the determinant of 00363 a matrix A, \f$ \ln|\det(A)| \f$, from its LU decomposition, LU. 00364 This function may be useful if the direct computation of the 00365 determinant would overflow or underflow. 00366 */ 00367 template<class mat_t> 00368 double LU_lndet(const size_t N, const mat_t &LU) { 00369 00370 size_t i; 00371 double lndet = 0.0; 00372 00373 for (i = 0; i < N; i++) { 00374 lndet+=log(fabs(O2SCL_IX2(LU,i,i))); 00375 } 00376 00377 return lndet; 00378 } 00379 00380 /** 00381 \brief Compute the sign of the 00382 determinant of a matrix from its LU decomposition 00383 00384 Compute the sign or phase factor of the determinant of a matrix 00385 A, \f$ \det(A)/|\det(A)| \f$, from its LU decomposition, LU. 00386 */ 00387 template<class mat_t> 00388 int LU_sgndet(const size_t N, const mat_t &LU, int signum) { 00389 00390 size_t i; 00391 int s=signum; 00392 00393 for (i=0;i<N;i++) { 00394 double u=O2SCL_IX2(LU,i,i); 00395 if (u<0) { 00396 s*=-1; 00397 } else if (u == 0) { 00398 s=0; 00399 break; 00400 } 00401 } 00402 00403 return s; 00404 } 00405 00406 #ifdef DOXYGENP 00407 } 00408 #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