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