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