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