dense_matrix.C
Go to the documentation of this file.00001 // The libMesh Finite Element Library. 00002 // Copyright (C) 2002-2012 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner 00003 00004 // This library is free software; you can redistribute it and/or 00005 // modify it under the terms of the GNU Lesser General Public 00006 // License as published by the Free Software Foundation; either 00007 // version 2.1 of the License, or (at your option) any later version. 00008 00009 // This library is distributed in the hope that it will be useful, 00010 // but WITHOUT ANY WARRANTY; without even the implied warranty of 00011 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU 00012 // Lesser General Public License for more details. 00013 00014 // You should have received a copy of the GNU Lesser General Public 00015 // License along with this library; if not, write to the Free Software 00016 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA 00017 00018 00019 // C++ Includes 00020 #include <cstdlib> // *must* precede <cmath> for proper std:abs() on PGI, Sun Studio CC 00021 #include <cmath> // for sqrt 00022 00023 // Local Includes 00024 #include "libmesh/dense_matrix.h" 00025 #include "libmesh/dense_vector.h" 00026 #include "libmesh/libmesh.h" 00027 00028 namespace libMesh 00029 { 00030 00031 00032 00033 // ------------------------------------------------------------ 00034 // Dense Matrix member functions 00035 00036 template<typename T> 00037 void DenseMatrix<T>::left_multiply (const DenseMatrixBase<T>& M2) 00038 { 00039 if (this->use_blas_lapack) 00040 this->_multiply_blas(M2, LEFT_MULTIPLY); 00041 else 00042 { 00043 // (*this) <- M2 * (*this) 00044 // Where: 00045 // (*this) = (m x n), 00046 // M2 = (m x p), 00047 // M3 = (p x n) 00048 00049 // M3 is a copy of *this before it gets resize()d 00050 DenseMatrix<T> M3(*this); 00051 00052 // Resize *this so that the result can fit 00053 this->resize (M2.m(), M3.n()); 00054 00055 // Call the multiply function in the base class 00056 this->multiply(*this, M2, M3); 00057 } 00058 } 00059 00060 00061 00062 00063 template<typename T> 00064 void DenseMatrix<T>::left_multiply_transpose(const DenseMatrix<T>& A) 00065 { 00066 if (this->use_blas_lapack) 00067 this->_multiply_blas(A, LEFT_MULTIPLY_TRANSPOSE); 00068 else 00069 { 00070 //Check to see if we are doing (A^T)*A 00071 if (this == &A) 00072 { 00073 //libmesh_here(); 00074 DenseMatrix<T> B(*this); 00075 00076 // Simple but inefficient way 00077 // return this->left_multiply_transpose(B); 00078 00079 // More efficient, but more code way 00080 // If A is mxn, the result will be a square matrix of Size n x n. 00081 const unsigned int n_rows = A.m(); 00082 const unsigned int n_cols = A.n(); 00083 00084 // resize() *this and also zero out all entries. 00085 this->resize(n_cols,n_cols); 00086 00087 // Compute the lower-triangular part 00088 for (unsigned int i=0; i<n_cols; ++i) 00089 for (unsigned int j=0; j<=i; ++j) 00090 for (unsigned int k=0; k<n_rows; ++k) // inner products are over n_rows 00091 (*this)(i,j) += B(k,i)*B(k,j); 00092 00093 // Copy lower-triangular part into upper-triangular part 00094 for (unsigned int i=0; i<n_cols; ++i) 00095 for (unsigned int j=i+1; j<n_cols; ++j) 00096 (*this)(i,j) = (*this)(j,i); 00097 } 00098 00099 else 00100 { 00101 DenseMatrix<T> B(*this); 00102 00103 this->resize (A.n(), B.n()); 00104 00105 libmesh_assert_equal_to (A.m(), B.m()); 00106 libmesh_assert_equal_to (this->m(), A.n()); 00107 libmesh_assert_equal_to (this->n(), B.n()); 00108 00109 const unsigned int m_s = A.n(); 00110 const unsigned int p_s = A.m(); 00111 const unsigned int n_s = this->n(); 00112 00113 // Do it this way because there is a 00114 // decent chance (at least for constraint matrices) 00115 // that A.transpose(i,k) = 0. 00116 for (unsigned int i=0; i<m_s; i++) 00117 for (unsigned int k=0; k<p_s; k++) 00118 if (A.transpose(i,k) != 0.) 00119 for (unsigned int j=0; j<n_s; j++) 00120 (*this)(i,j) += A.transpose(i,k)*B(k,j); 00121 } 00122 } 00123 00124 } 00125 00126 00127 00128 00129 00130 00131 template<typename T> 00132 void DenseMatrix<T>::right_multiply (const DenseMatrixBase<T>& M3) 00133 { 00134 if (this->use_blas_lapack) 00135 this->_multiply_blas(M3, RIGHT_MULTIPLY); 00136 else 00137 { 00138 // (*this) <- M3 * (*this) 00139 // Where: 00140 // (*this) = (m x n), 00141 // M2 = (m x p), 00142 // M3 = (p x n) 00143 00144 // M2 is a copy of *this before it gets resize()d 00145 DenseMatrix<T> M2(*this); 00146 00147 // Resize *this so that the result can fit 00148 this->resize (M2.m(), M3.n()); 00149 00150 this->multiply(*this, M2, M3); 00151 } 00152 } 00153 00154 00155 00156 00157 template<typename T> 00158 void DenseMatrix<T>::right_multiply_transpose (const DenseMatrix<T>& B) 00159 { 00160 if (this->use_blas_lapack) 00161 this->_multiply_blas(B, RIGHT_MULTIPLY_TRANSPOSE); 00162 else 00163 { 00164 //Check to see if we are doing B*(B^T) 00165 if (this == &B) 00166 { 00167 //libmesh_here(); 00168 DenseMatrix<T> A(*this); 00169 00170 // Simple but inefficient way 00171 // return this->right_multiply_transpose(A); 00172 00173 // More efficient, more code 00174 // If B is mxn, the result will be a square matrix of Size m x m. 00175 const unsigned int n_rows = B.m(); 00176 const unsigned int n_cols = B.n(); 00177 00178 // resize() *this and also zero out all entries. 00179 this->resize(n_rows,n_rows); 00180 00181 // Compute the lower-triangular part 00182 for (unsigned int i=0; i<n_rows; ++i) 00183 for (unsigned int j=0; j<=i; ++j) 00184 for (unsigned int k=0; k<n_cols; ++k) // inner products are over n_cols 00185 (*this)(i,j) += A(i,k)*A(j,k); 00186 00187 // Copy lower-triangular part into upper-triangular part 00188 for (unsigned int i=0; i<n_rows; ++i) 00189 for (unsigned int j=i+1; j<n_rows; ++j) 00190 (*this)(i,j) = (*this)(j,i); 00191 } 00192 00193 else 00194 { 00195 DenseMatrix<T> A(*this); 00196 00197 this->resize (A.m(), B.m()); 00198 00199 libmesh_assert_equal_to (A.n(), B.n()); 00200 libmesh_assert_equal_to (this->m(), A.m()); 00201 libmesh_assert_equal_to (this->n(), B.m()); 00202 00203 const unsigned int m_s = A.m(); 00204 const unsigned int p_s = A.n(); 00205 const unsigned int n_s = this->n(); 00206 00207 // Do it this way because there is a 00208 // decent chance (at least for constraint matrices) 00209 // that B.transpose(k,j) = 0. 00210 for (unsigned int j=0; j<n_s; j++) 00211 for (unsigned int k=0; k<p_s; k++) 00212 if (B.transpose(k,j) != 0.) 00213 for (unsigned int i=0; i<m_s; i++) 00214 (*this)(i,j) += A(i,k)*B.transpose(k,j); 00215 } 00216 } 00217 } 00218 00219 00220 00221 template<typename T> 00222 void DenseMatrix<T>::vector_mult (DenseVector<T>& dest, 00223 const DenseVector<T>& arg) const 00224 { 00225 // Make sure the input sizes are compatible 00226 libmesh_assert_equal_to (this->n(), arg.size()); 00227 00228 // Resize and clear dest. 00229 // Note: DenseVector::resize() also zeros the vector. 00230 dest.resize(this->m()); 00231 00232 // Short-circuit if the matrix is empty 00233 if(this->m() == 0 || this->n() == 0) 00234 return; 00235 00236 if (this->use_blas_lapack) 00237 this->_matvec_blas(1., 0., dest, arg); 00238 else 00239 { 00240 const unsigned int n_rows = this->m(); 00241 const unsigned int n_cols = this->n(); 00242 00243 for(unsigned int i=0; i<n_rows; i++) 00244 for(unsigned int j=0; j<n_cols; j++) 00245 dest(i) += (*this)(i,j)*arg(j); 00246 } 00247 } 00248 00249 00250 00251 00252 template<typename T> 00253 void DenseMatrix<T>::vector_mult_transpose (DenseVector<T>& dest, 00254 const DenseVector<T>& arg) const 00255 { 00256 // Make sure the input sizes are compatible 00257 libmesh_assert_equal_to (this->m(), arg.size()); 00258 00259 // Resize and clear dest. 00260 // Note: DenseVector::resize() also zeros the vector. 00261 dest.resize(this->n()); 00262 00263 // Short-circuit if the matrix is empty 00264 if(this->m() == 0) 00265 return; 00266 00267 if (this->use_blas_lapack) 00268 { 00269 this->_matvec_blas(1., 0., dest, arg, /*trans=*/true); 00270 } 00271 else 00272 { 00273 const unsigned int n_rows = this->m(); 00274 const unsigned int n_cols = this->n(); 00275 00276 // WORKS 00277 // for(unsigned int j=0; j<n_cols; j++) 00278 // for(unsigned int i=0; i<n_rows; i++) 00279 // dest(j) += (*this)(i,j)*arg(i); 00280 00281 // ALSO WORKS, (i,j) just swapped 00282 for(unsigned int i=0; i<n_cols; i++) 00283 for(unsigned int j=0; j<n_rows; j++) 00284 dest(i) += (*this)(j,i)*arg(j); 00285 } 00286 } 00287 00288 00289 00290 template<typename T> 00291 void DenseMatrix<T>::vector_mult_add (DenseVector<T>& dest, 00292 const T factor, 00293 const DenseVector<T>& arg) const 00294 { 00295 // Short-circuit if the matrix is empty 00296 if(this->m() == 0) 00297 { 00298 dest.resize(0); 00299 return; 00300 } 00301 00302 if (this->use_blas_lapack) 00303 this->_matvec_blas(factor, 1., dest, arg); 00304 else 00305 { 00306 DenseVector<T> temp(arg.size()); 00307 this->vector_mult(temp, arg); 00308 dest.add(factor, temp); 00309 } 00310 } 00311 00312 00313 00314 template<typename T> 00315 void DenseMatrix<T>::get_principal_submatrix (unsigned int sub_m, 00316 unsigned int sub_n, 00317 DenseMatrix<T>& dest) const 00318 { 00319 libmesh_assert( (sub_m <= this->m()) && (sub_n <= this->n()) ); 00320 00321 dest.resize(sub_m, sub_n); 00322 for(unsigned int i=0; i<sub_m; i++) 00323 for(unsigned int j=0; j<sub_n; j++) 00324 dest(i,j) = (*this)(i,j); 00325 } 00326 00327 00328 00329 template<typename T> 00330 void DenseMatrix<T>::get_principal_submatrix (unsigned int sub_m, DenseMatrix<T>& dest) const 00331 { 00332 get_principal_submatrix(sub_m, sub_m, dest); 00333 } 00334 00335 00336 00337 template<typename T> 00338 void DenseMatrix<T>::get_transpose (DenseMatrix<T>& dest) const 00339 { 00340 dest.resize(this->n(), this->m()); 00341 00342 for (unsigned int i=0; i<dest.m(); i++) 00343 for (unsigned int j=0; j<dest.n(); j++) 00344 dest(i,j) = (*this)(j,i); 00345 } 00346 00347 00348 00349 00350 template<typename T> 00351 void DenseMatrix<T>::lu_solve (const DenseVector<T>& b, 00352 DenseVector<T>& x) 00353 { 00354 // Check to be sure that the matrix is square before attempting 00355 // an LU-solve. In general, one can compute the LU factorization of 00356 // a non-square matrix, but: 00357 // 00358 // Overdetermined systems (m>n) have a solution only if enough of 00359 // the equations are linearly-dependent. 00360 // 00361 // Underdetermined systems (m<n) typically have infinitely many 00362 // solutions. 00363 // 00364 // We don't want to deal with either of these ambiguous cases here... 00365 libmesh_assert_equal_to (this->m(), this->n()); 00366 00367 switch(this->_decomposition_type) 00368 { 00369 case NONE: 00370 { 00371 if (this->use_blas_lapack) 00372 this->_lu_decompose_lapack(); 00373 else 00374 this->_lu_decompose (); 00375 break; 00376 } 00377 00378 case LU_BLAS_LAPACK: 00379 { 00380 // Already factored, just need to call back_substitute. 00381 if (this->use_blas_lapack) 00382 break; 00383 } 00384 00385 case LU: 00386 { 00387 // Already factored, just need to call back_substitute. 00388 if ( !(this->use_blas_lapack) ) 00389 break; 00390 } 00391 00392 default: 00393 { 00394 libMesh::err << "Error! This matrix already has a " 00395 << "different decomposition..." 00396 << std::endl; 00397 libmesh_error(); 00398 } 00399 } 00400 00401 if (this->use_blas_lapack) 00402 this->_lu_back_substitute_lapack (b, x); 00403 else 00404 this->_lu_back_substitute (b, x); 00405 } 00406 00407 00408 00409 00410 00411 00412 template<typename T> 00413 void DenseMatrix<T>::_lu_back_substitute (const DenseVector<T>& b, 00414 DenseVector<T>& x ) const 00415 { 00416 const unsigned int 00417 n_cols = this->n(); 00418 00419 libmesh_assert_equal_to (this->m(), n_cols); 00420 libmesh_assert_equal_to (this->m(), b.size()); 00421 00422 x.resize (n_cols); 00423 00424 // A convenient reference to *this 00425 const DenseMatrix<T>& A = *this; 00426 00427 // Temporary vector storage. We use this instead of 00428 // modifying the RHS. 00429 DenseVector<T> z = b; 00430 00431 // Lower-triangular "top to bottom" solve step, taking into account pivots 00432 for (unsigned int i=0; i<n_cols; ++i) 00433 { 00434 // Swap 00435 if (_pivots[i] != static_cast<int>(i)) 00436 std::swap( z(i), z(_pivots[i]) ); 00437 00438 x(i) = z(i); 00439 00440 for (unsigned int j=0; j<i; ++j) 00441 x(i) -= A(i,j)*x(j); 00442 00443 x(i) /= A(i,i); 00444 } 00445 00446 // Upper-triangular "bottom to top" solve step 00447 const unsigned int last_row = n_cols-1; 00448 00449 for (int i=last_row; i>=0; --i) 00450 { 00451 for (int j=i+1; j<static_cast<int>(n_cols); ++j) 00452 x(i) -= A(i,j)*x(j); 00453 } 00454 } 00455 00456 00457 00458 00459 00460 00461 00462 00463 template<typename T> 00464 void DenseMatrix<T>::_lu_decompose () 00465 { 00466 // If this function was called, there better not be any 00467 // previous decomposition of the matrix. 00468 libmesh_assert_equal_to (this->_decomposition_type, NONE); 00469 00470 // Get the matrix size and make sure it is square 00471 const unsigned int 00472 n_rows = this->m(); 00473 00474 // A convenient reference to *this 00475 DenseMatrix<T>& A = *this; 00476 00477 _pivots.resize(n_rows); 00478 00479 for (unsigned int i=0; i<n_rows; ++i) 00480 { 00481 // Find the pivot row by searching down the i'th column 00482 _pivots[i] = i; 00483 00484 // std::abs(complex) must return a Real! 00485 Real the_max = std::abs( A(i,i) ); 00486 for (unsigned int j=i+1; j<n_rows; ++j) 00487 { 00488 Real candidate_max = std::abs( A(j,i) ); 00489 if (the_max < candidate_max) 00490 { 00491 the_max = candidate_max; 00492 _pivots[i] = j; 00493 } 00494 } 00495 00496 // libMesh::out << "the_max=" << the_max << " found at row " << _pivots[i] << std::endl; 00497 00498 // If the max was found in a different row, interchange rows. 00499 // Here we interchange the *entire* row, in Gaussian elimination 00500 // you would only interchange the subrows A(i,j) and A(p(i),j), for j>i 00501 if (_pivots[i] != static_cast<int>(i)) 00502 { 00503 for (unsigned int j=0; j<n_rows; ++j) 00504 std::swap( A(i,j), A(_pivots[i], j) ); 00505 } 00506 00507 00508 // If the max abs entry found is zero, the matrix is singular 00509 if (A(i,i) == libMesh::zero) 00510 { 00511 libMesh::out << "Matrix A is singular!" << std::endl; 00512 libmesh_error(); 00513 } 00514 00515 // Scale upper triangle entries of row i by the diagonal entry 00516 // Note: don't scale the diagonal entry itself! 00517 const T diag_inv = 1. / A(i,i); 00518 for (unsigned int j=i+1; j<n_rows; ++j) 00519 A(i,j) *= diag_inv; 00520 00521 // Update the remaining sub-matrix A[i+1:m][i+1:m] 00522 // by subtracting off (the diagonal-scaled) 00523 // upper-triangular part of row i, scaled by the 00524 // i'th column entry of each row. In terms of 00525 // row operations, this is: 00526 // for each r > i 00527 // SubRow(r) = SubRow(r) - A(r,i)*SubRow(i) 00528 // 00529 // If we were scaling the i'th column as well, like 00530 // in Gaussian elimination, this would 'zero' the 00531 // entry in the i'th column. 00532 for (unsigned int row=i+1; row<n_rows; ++row) 00533 for (unsigned int col=i+1; col<n_rows; ++col) 00534 A(row,col) -= A(row,i) * A(i,col); 00535 00536 } // end i loop 00537 00538 // Set the flag for LU decomposition 00539 this->_decomposition_type = LU; 00540 } 00541 00542 00543 00544 template<typename T> 00545 void DenseMatrix<T>::svd (DenseVector<T>& sigma) 00546 { 00547 // We use the LAPACK svd implementation 00548 _svd_lapack(sigma); 00549 } 00550 00551 00552 template<typename T> 00553 void DenseMatrix<T>::svd (DenseVector<T>& sigma, DenseMatrix<T>& U, DenseMatrix<T>& VT) 00554 { 00555 // We use the LAPACK svd implementation 00556 _svd_lapack(sigma, U, VT); 00557 } 00558 00559 00560 00561 template<typename T> 00562 T DenseMatrix<T>::det () 00563 { 00564 switch(this->_decomposition_type) 00565 { 00566 case NONE: 00567 { 00568 // First LU decompose the matrix. 00569 // Note that the lu_decompose routine will check to see if the 00570 // matrix is square so we don't worry about it. 00571 if (this->use_blas_lapack) 00572 this->_lu_decompose_lapack(); 00573 else 00574 this->_lu_decompose (); 00575 } 00576 case LU: 00577 case LU_BLAS_LAPACK: 00578 { 00579 // Already decomposed, don't do anything 00580 break; 00581 } 00582 default: 00583 { 00584 libMesh::err << "Error! Can't compute the determinant under " 00585 << "the current decomposition." 00586 << std::endl; 00587 libmesh_error(); 00588 } 00589 } 00590 00591 // A variable to keep track of the running product of diagonal terms. 00592 T determinant = 1.; 00593 00594 // Loop over diagonal terms, computing the product. In practice, 00595 // be careful because this value could easily become too large to 00596 // fit in a double or float. To be safe, one should keep track of 00597 // the power (of 10) of the determinant in a separate variable 00598 // and maintain an order 1 value for the determinant itself. 00599 unsigned int n_interchanges = 0; 00600 for (unsigned int i=0; i<this->m(); i++) 00601 { 00602 if (this->_decomposition_type==LU) 00603 if (_pivots[i] != static_cast<int>(i)) 00604 n_interchanges++; 00605 00606 // Lapack pivots are 1-based! 00607 if (this->_decomposition_type==LU_BLAS_LAPACK) 00608 if (_pivots[i] != static_cast<int>(i+1)) 00609 n_interchanges++; 00610 00611 determinant *= (*this)(i,i); 00612 } 00613 00614 // Compute sign of determinant, depends on number of row interchanges! 00615 // The sign should be (-1)^{n}, where n is the number of interchanges. 00616 Real sign = n_interchanges % 2 == 0 ? 1. : -1.; 00617 00618 return sign*determinant; 00619 } 00620 00621 00622 00623 // The cholesky solve function first decomposes the matrix 00624 // with cholesky_decompose and then uses the cholesky_back_substitute 00625 // routine to find the solution x. 00626 template <typename T> 00627 template <typename T2> 00628 void DenseMatrix<T>::cholesky_solve (const DenseVector<T2>& b, 00629 DenseVector<T2>& x) 00630 { 00631 // Check for a previous decomposition 00632 switch(this->_decomposition_type) 00633 { 00634 case NONE: 00635 { 00636 this->_cholesky_decompose (); 00637 break; 00638 } 00639 00640 case CHOLESKY: 00641 { 00642 // Already factored, just need to call back_substitute. 00643 break; 00644 } 00645 00646 default: 00647 { 00648 libMesh::err << "Error! This matrix already has a " 00649 << "different decomposition..." 00650 << std::endl; 00651 libmesh_error(); 00652 } 00653 } 00654 00655 // Perform back substitution 00656 this->_cholesky_back_substitute (b, x); 00657 } 00658 00659 00660 00661 00662 // This algorithm is based on the Cholesky decomposition in 00663 // the Numerical Recipes in C book. 00664 template<typename T> 00665 void DenseMatrix<T>::_cholesky_decompose () 00666 { 00667 // If we called this function, there better not be any 00668 // previous decomposition of the matrix. 00669 libmesh_assert_equal_to (this->_decomposition_type, NONE); 00670 00671 // Shorthand notation for number of rows and columns. 00672 const unsigned int 00673 n_rows = this->m(), 00674 n_cols = this->n(); 00675 00676 // Just to be really sure... 00677 libmesh_assert_equal_to (n_rows, n_cols); 00678 00679 // A convenient reference to *this 00680 DenseMatrix<T>& A = *this; 00681 00682 for (unsigned int i=0; i<n_rows; ++i) 00683 { 00684 for (unsigned int j=i; j<n_cols; ++j) 00685 { 00686 for (unsigned int k=0; k<i; ++k) 00687 A(i,j) -= A(i,k) * A(j,k); 00688 00689 if (i == j) 00690 { 00691 #ifndef LIBMESH_USE_COMPLEX_NUMBERS 00692 if (A(i,j) <= 0.0) 00693 { 00694 libMesh::err << "Error! Can only use Cholesky decomposition " 00695 << "with symmetric positive definite matrices." 00696 << std::endl; 00697 libmesh_error(); 00698 } 00699 #endif 00700 00701 A(i,i) = std::sqrt(A(i,j)); 00702 } 00703 else 00704 A(j,i) = A(i,j) / A(i,i); 00705 } 00706 } 00707 00708 // Set the flag for CHOLESKY decomposition 00709 this->_decomposition_type = CHOLESKY; 00710 } 00711 00712 00713 00714 template <typename T> 00715 template <typename T2> 00716 void DenseMatrix<T>::_cholesky_back_substitute (const DenseVector<T2>& b, 00717 DenseVector<T2>& x) const 00718 { 00719 // Shorthand notation for number of rows and columns. 00720 const unsigned int 00721 n_rows = this->m(), 00722 n_cols = this->n(); 00723 00724 // Just to be really sure... 00725 libmesh_assert_equal_to (n_rows, n_cols); 00726 00727 // A convenient reference to *this 00728 const DenseMatrix<T>& A = *this; 00729 00730 // Now compute the solution to Ax =b using the factorization. 00731 x.resize(n_rows); 00732 00733 // Solve for Ly=b 00734 for (unsigned int i=0; i<n_cols; ++i) 00735 { 00736 T2 temp = b(i); 00737 00738 for (unsigned int k=0; k<i; ++k) 00739 temp -= A(i,k)*x(k); 00740 00741 x(i) = temp / A(i,i); 00742 } 00743 00744 // Solve for L^T x = y 00745 for (unsigned int i=0; i<n_cols; ++i) 00746 { 00747 const unsigned int ib = (n_cols-1)-i; 00748 00749 for (unsigned int k=(ib+1); k<n_cols; ++k) 00750 x(ib) -= A(k,ib) * x(k); 00751 00752 x(ib) /= A(ib,ib); 00753 } 00754 } 00755 00756 00757 00758 00759 00760 00761 00762 00763 // This routine is commented out since it is not really a memory 00764 // efficient implementation. Also, you don't *need* the inverse 00765 // for anything, instead just use lu_solve to solve Ax=b. 00766 // template<typename T> 00767 // void DenseMatrix<T>::inverse () 00768 // { 00769 // // First LU decompose the matrix 00770 // // Note that the lu_decompose routine will check to see if the 00771 // // matrix is square so we don't worry about it. 00772 // if (!this->_lu_decomposed) 00773 // this->_lu_decompose(); 00774 00775 // // A unit vector which will be used as a rhs 00776 // // to pick off a single value each time. 00777 // DenseVector<T> e; 00778 // e.resize(this->m()); 00779 00780 // // An empty vector which will be used to hold the solution 00781 // // to the back substitutions. 00782 // DenseVector<T> x; 00783 // x.resize(this->m()); 00784 00785 // // An empty dense matrix to store the resulting inverse 00786 // // temporarily until we can overwrite A. 00787 // DenseMatrix<T> inv; 00788 // inv.resize(this->m(), this->n()); 00789 00790 // // Resize the passed in matrix to hold the inverse 00791 // inv.resize(this->m(), this->n()); 00792 00793 // for (unsigned int j=0; j<this->n(); ++j) 00794 // { 00795 // e.zero(); 00796 // e(j) = 1.; 00797 // this->_lu_back_substitute(e, x, false); 00798 // for (unsigned int i=0; i<this->n(); ++i) 00799 // inv(i,j) = x(i); 00800 // } 00801 00802 // // Now overwrite all the entries 00803 // *this = inv; 00804 // } 00805 00806 00807 //-------------------------------------------------------------- 00808 // Explicit instantiations 00809 template class DenseMatrix<Real>; 00810 template void DenseMatrix<Real>::cholesky_solve(const DenseVector<Real>&, DenseVector<Real>&); 00811 template void DenseMatrix<Real>::_cholesky_back_substitute(const DenseVector<Real>&, DenseVector<Real>&) const; 00812 template void DenseMatrix<Real>::cholesky_solve(const DenseVector<Complex>&, DenseVector<Complex>&); 00813 template void DenseMatrix<Real>::_cholesky_back_substitute(const DenseVector<Complex>&, DenseVector<Complex>&) const; 00814 00815 #ifdef LIBMESH_USE_COMPLEX_NUMBERS 00816 template class DenseMatrix<Complex>; 00817 template void DenseMatrix<Complex>::cholesky_solve(const DenseVector<Complex>&,DenseVector<Complex>&); 00818 template void DenseMatrix<Complex>::_cholesky_back_substitute(const DenseVector<Complex>&, DenseVector<Complex>&) const; 00819 #endif 00820 00821 } // namespace libMesh
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:46 UTC
Hosted By: