dense_matrix.h
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 00020 #ifndef LIBMESH_DENSE_MATRIX_H 00021 #define LIBMESH_DENSE_MATRIX_H 00022 00023 // Local Includes 00024 #include "libmesh/libmesh_common.h" 00025 #include "libmesh/dense_matrix_base.h" 00026 00027 // C++ includes 00028 #include <vector> 00029 #include <algorithm> 00030 #include <cstring> // std::memset 00031 00032 namespace libMesh 00033 { 00034 00035 // Forward Declarations 00036 template <typename T> class DenseVector; 00037 00038 00039 00048 // ------------------------------------------------------------ 00049 // Dense Matrix class definition 00050 template<typename T> 00051 class DenseMatrix : public DenseMatrixBase<T> 00052 { 00053 public: 00054 00058 DenseMatrix(const unsigned int new_m=0, 00059 const unsigned int new_n=0); 00060 00064 //DenseMatrix (const DenseMatrix<T>& other_matrix); 00065 00069 virtual ~DenseMatrix() {} 00070 00071 00075 virtual void zero(); 00076 00080 T operator() (const unsigned int i, 00081 const unsigned int j) const; 00082 00086 T & operator() (const unsigned int i, 00087 const unsigned int j); 00088 00092 virtual T el(const unsigned int i, 00093 const unsigned int j) const { return (*this)(i,j); } 00094 00098 virtual T & el(const unsigned int i, 00099 const unsigned int j) { return (*this)(i,j); } 00100 00104 virtual void left_multiply (const DenseMatrixBase<T>& M2); 00105 00109 virtual void right_multiply (const DenseMatrixBase<T>& M3); 00110 00115 void vector_mult(DenseVector<T>& dest, const DenseVector<T>& arg) const; 00116 00121 void vector_mult_transpose(DenseVector<T>& dest, const DenseVector<T>& arg) const; 00122 00127 void vector_mult_add (DenseVector<T>& dest, 00128 const T factor, 00129 const DenseVector<T>& arg) const; 00130 00134 void get_principal_submatrix (unsigned int sub_m, unsigned int sub_n, DenseMatrix<T>& dest) const; 00135 00139 void get_principal_submatrix (unsigned int sub_m, DenseMatrix<T>& dest) const; 00140 00144 DenseMatrix<T>& operator = (const DenseMatrix<T>& other_matrix); 00145 00149 void swap(DenseMatrix<T>& other_matrix); 00150 00155 void resize(const unsigned int new_m, 00156 const unsigned int new_n); 00157 00161 void scale (const T factor); 00162 00166 DenseMatrix<T>& operator *= (const T factor); 00167 00171 void add (const T factor, 00172 const DenseMatrix<T>& mat); 00173 00177 bool operator== (const DenseMatrix<T> &mat) const; 00178 00182 bool operator!= (const DenseMatrix<T> &mat) const; 00183 00187 DenseMatrix<T>& operator+= (const DenseMatrix<T> &mat); 00188 00192 DenseMatrix<T>& operator-= (const DenseMatrix<T> &mat); 00193 00199 Real min () const; 00200 00206 Real max () const; 00207 00218 Real l1_norm () const; 00219 00231 Real linfty_norm () const; 00232 00236 void left_multiply_transpose (const DenseMatrix<T>& A); 00237 00238 00242 void right_multiply_transpose (const DenseMatrix<T>& A); 00243 00247 T transpose (const unsigned int i, 00248 const unsigned int j) const; 00249 00253 void get_transpose(DenseMatrix<T>& dest) const; 00254 00260 std::vector<T>& get_values() { return _val; } 00261 00265 const std::vector<T>& get_values() const { return _val; } 00266 00273 void condense(const unsigned int i, 00274 const unsigned int j, 00275 const T val, 00276 DenseVector<T>& rhs) 00277 { DenseMatrixBase<T>::condense (i, j, val, rhs); } 00278 00284 void lu_solve (const DenseVector<T>& b, 00285 DenseVector<T>& x); 00286 00287 00288 00298 template <typename T2> 00299 void cholesky_solve(const DenseVector<T2>& b, 00300 DenseVector<T2>& x); 00301 00302 00311 void svd(DenseVector<T>& sigma); 00312 00313 00325 void svd(DenseVector<T>& sigma, DenseMatrix<T>& U, DenseMatrix<T>& VT); 00326 00327 00333 T det(); 00334 00343 // void inverse(); 00344 00350 bool use_blas_lapack; 00351 00352 private: 00353 00357 std::vector<T> _val; 00358 00364 void _lu_decompose (); 00365 00371 void _lu_back_substitute (const DenseVector<T>& b, 00372 DenseVector<T>& x) const; 00373 00380 void _cholesky_decompose(); 00381 00388 template <typename T2> 00389 void _cholesky_back_substitute(const DenseVector<T2>& b, 00390 DenseVector<T2>& x) const; 00391 00399 enum DecompositionType {LU=0, CHOLESKY=1, LU_BLAS_LAPACK, NONE}; 00400 00405 DecompositionType _decomposition_type; 00406 00411 enum _BLAS_Multiply_Flag { 00412 LEFT_MULTIPLY = 0, 00413 RIGHT_MULTIPLY, 00414 LEFT_MULTIPLY_TRANSPOSE, 00415 RIGHT_MULTIPLY_TRANSPOSE 00416 }; 00417 00425 void _multiply_blas(const DenseMatrixBase<T>& other, 00426 _BLAS_Multiply_Flag flag); 00427 00437 void _lu_decompose_lapack(); 00438 00444 void _svd_lapack(DenseVector<T>& sigma); 00445 00451 void _svd_lapack(DenseVector<T>& sigma, 00452 DenseMatrix<T>& U, 00453 DenseMatrix<T>& VT); 00454 00459 void _svd_helper (char JOBU, 00460 char JOBVT, 00461 std::vector<T>& sigma_val, 00462 std::vector<T>& U_val, 00463 std::vector<T>& VT_val); 00464 00470 std::vector<int> _pivots; 00471 00481 void _lu_back_substitute_lapack (const DenseVector<T>& b, 00482 DenseVector<T>& x); 00483 00496 void _matvec_blas(T alpha, T beta, 00497 DenseVector<T>& dest, 00498 const DenseVector<T>& arg, 00499 bool trans=false) const; 00500 }; 00501 00502 00503 00504 00505 00506 // ------------------------------------------------------------ 00510 namespace DenseMatrices 00511 { 00512 00517 typedef DenseMatrix<Real> RealDenseMatrix; 00518 00528 typedef DenseMatrix<Complex> ComplexDenseMatrix; 00529 00530 } 00531 00532 00533 00534 using namespace DenseMatrices; 00535 00536 00537 00538 00539 00540 // ------------------------------------------------------------ 00541 // Dense Matrix member functions 00542 template<typename T> 00543 inline 00544 DenseMatrix<T>::DenseMatrix(const unsigned int new_m, 00545 const unsigned int new_n) 00546 : DenseMatrixBase<T>(new_m,new_n), 00547 #if defined(LIBMESH_HAVE_PETSC) && defined(LIBMESH_USE_REAL_NUMBERS) && defined(LIBMESH_DEFAULT_DOUBLE_PRECISION) 00548 use_blas_lapack(true), 00549 #else 00550 use_blas_lapack(false), 00551 #endif 00552 _val(), 00553 _decomposition_type(NONE), 00554 _pivots() 00555 { 00556 this->resize(new_m,new_n); 00557 } 00558 00559 00560 00561 // FIXME[JWP]: This copy ctor has not been maintained along with 00562 // the rest of the class... 00563 // Can we just use the compiler-generated copy ctor here? 00564 // template<typename T> 00565 // inline 00566 // DenseMatrix<T>::DenseMatrix (const DenseMatrix<T>& other_matrix) 00567 // : DenseMatrixBase<T>(other_matrix._m, other_matrix._n) 00568 // { 00569 // _val = other_matrix._val; 00570 // } 00571 00572 00573 00574 template<typename T> 00575 inline 00576 void DenseMatrix<T>::swap(DenseMatrix<T>& other_matrix) 00577 { 00578 std::swap(this->_m, other_matrix._m); 00579 std::swap(this->_n, other_matrix._n); 00580 _val.swap(other_matrix._val); 00581 DecompositionType _temp = _decomposition_type; 00582 _decomposition_type = other_matrix._decomposition_type; 00583 other_matrix._decomposition_type = _temp; 00584 } 00585 00586 00587 00588 template<typename T> 00589 inline 00590 void DenseMatrix<T>::resize(const unsigned int new_m, 00591 const unsigned int new_n) 00592 { 00593 _val.resize(new_m*new_n); 00594 00595 this->_m = new_m; 00596 this->_n = new_n; 00597 00598 _decomposition_type = NONE; 00599 this->zero(); 00600 } 00601 00602 00603 00604 template<typename T> 00605 inline 00606 void DenseMatrix<T>::zero() 00607 { 00608 _decomposition_type = NONE; 00609 00610 // Just doing this ifdef to be completely safe 00611 #ifndef LIBMESH_USE_COMPLEX_NUMBERS 00612 if(_val.size()) 00613 std::memset(&_val[0], 0, sizeof(T) * _val.size()); 00614 #else 00615 std::fill (_val.begin(), _val.end(), 0.); 00616 #endif 00617 } 00618 00619 00620 00621 template<typename T> 00622 inline 00623 DenseMatrix<T>& DenseMatrix<T>::operator = (const DenseMatrix<T>& other_matrix) 00624 { 00625 this->_m = other_matrix._m; 00626 this->_n = other_matrix._n; 00627 00628 _val = other_matrix._val; 00629 _decomposition_type = other_matrix._decomposition_type; 00630 00631 return *this; 00632 } 00633 00634 00635 00636 template<typename T> 00637 inline 00638 T DenseMatrix<T>::operator () (const unsigned int i, 00639 const unsigned int j) const 00640 { 00641 libmesh_assert_less (i*j, _val.size()); 00642 libmesh_assert_less (i, this->_m); 00643 libmesh_assert_less (j, this->_n); 00644 00645 00646 // return _val[(i) + (this->_m)*(j)]; // col-major 00647 return _val[(i)*(this->_n) + (j)]; // row-major 00648 } 00649 00650 00651 00652 template<typename T> 00653 inline 00654 T & DenseMatrix<T>::operator () (const unsigned int i, 00655 const unsigned int j) 00656 { 00657 libmesh_assert_less (i*j, _val.size()); 00658 libmesh_assert_less (i, this->_m); 00659 libmesh_assert_less (j, this->_n); 00660 00661 //return _val[(i) + (this->_m)*(j)]; // col-major 00662 return _val[(i)*(this->_n) + (j)]; // row-major 00663 } 00664 00665 00666 00667 00668 00669 template<typename T> 00670 inline 00671 void DenseMatrix<T>::scale (const T factor) 00672 { 00673 for (unsigned int i=0; i<_val.size(); i++) 00674 _val[i] *= factor; 00675 } 00676 00677 00678 00679 template<typename T> 00680 inline 00681 DenseMatrix<T>& DenseMatrix<T>::operator *= (const T factor) 00682 { 00683 this->scale(factor); 00684 return *this; 00685 } 00686 00687 00688 00689 template<typename T> 00690 inline 00691 void DenseMatrix<T>::add (const T factor, const DenseMatrix<T>& mat) 00692 { 00693 for (unsigned int i=0; i<_val.size(); i++) 00694 _val[i] += factor * mat._val[i]; 00695 } 00696 00697 00698 00699 template<typename T> 00700 inline 00701 bool DenseMatrix<T>::operator == (const DenseMatrix<T> &mat) const 00702 { 00703 for (unsigned int i=0; i<_val.size(); i++) 00704 if (_val[i] != mat._val[i]) 00705 return false; 00706 00707 return true; 00708 } 00709 00710 00711 00712 template<typename T> 00713 inline 00714 bool DenseMatrix<T>::operator != (const DenseMatrix<T> &mat) const 00715 { 00716 for (unsigned int i=0; i<_val.size(); i++) 00717 if (_val[i] != mat._val[i]) 00718 return true; 00719 00720 return false; 00721 } 00722 00723 00724 00725 template<typename T> 00726 inline 00727 DenseMatrix<T>& DenseMatrix<T>::operator += (const DenseMatrix<T> &mat) 00728 { 00729 for (unsigned int i=0; i<_val.size(); i++) 00730 _val[i] += mat._val[i]; 00731 00732 return *this; 00733 } 00734 00735 00736 00737 template<typename T> 00738 inline 00739 DenseMatrix<T>& DenseMatrix<T>::operator -= (const DenseMatrix<T> &mat) 00740 { 00741 for (unsigned int i=0; i<_val.size(); i++) 00742 _val[i] -= mat._val[i]; 00743 00744 return *this; 00745 } 00746 00747 00748 00749 template<typename T> 00750 inline 00751 Real DenseMatrix<T>::min () const 00752 { 00753 libmesh_assert (this->_m); 00754 libmesh_assert (this->_n); 00755 Real my_min = libmesh_real((*this)(0,0)); 00756 00757 for (unsigned int i=0; i!=this->_m; i++) 00758 { 00759 for (unsigned int j=0; j!=this->_n; j++) 00760 { 00761 Real current = libmesh_real((*this)(i,j)); 00762 my_min = (my_min < current? my_min : current); 00763 } 00764 } 00765 return my_min; 00766 } 00767 00768 00769 00770 template<typename T> 00771 inline 00772 Real DenseMatrix<T>::max () const 00773 { 00774 libmesh_assert (this->_m); 00775 libmesh_assert (this->_n); 00776 Real my_max = libmesh_real((*this)(0,0)); 00777 00778 for (unsigned int i=0; i!=this->_m; i++) 00779 { 00780 for (unsigned int j=0; j!=this->_n; j++) 00781 { 00782 Real current = libmesh_real((*this)(i,j)); 00783 my_max = (my_max > current? my_max : current); 00784 } 00785 } 00786 return my_max; 00787 } 00788 00789 00790 00791 template<typename T> 00792 inline 00793 Real DenseMatrix<T>::l1_norm () const 00794 { 00795 libmesh_assert (this->_m); 00796 libmesh_assert (this->_n); 00797 00798 Real columnsum = 0.; 00799 for (unsigned int i=0; i!=this->_m; i++) 00800 { 00801 columnsum += std::abs((*this)(i,0)); 00802 } 00803 Real my_max = columnsum; 00804 for (unsigned int j=1; j!=this->_n; j++) 00805 { 00806 columnsum = 0.; 00807 for (unsigned int i=0; i!=this->_m; i++) 00808 { 00809 columnsum += std::abs((*this)(i,j)); 00810 } 00811 my_max = (my_max > columnsum? my_max : columnsum); 00812 } 00813 return my_max; 00814 } 00815 00816 00817 00818 template<typename T> 00819 inline 00820 Real DenseMatrix<T>::linfty_norm () const 00821 { 00822 libmesh_assert (this->_m); 00823 libmesh_assert (this->_n); 00824 00825 Real rowsum = 0.; 00826 for (unsigned int j=0; j!=this->_n; j++) 00827 { 00828 rowsum += std::abs((*this)(0,j)); 00829 } 00830 Real my_max = rowsum; 00831 for (unsigned int i=1; i!=this->_m; i++) 00832 { 00833 rowsum = 0.; 00834 for (unsigned int j=0; j!=this->_n; j++) 00835 { 00836 rowsum += std::abs((*this)(i,j)); 00837 } 00838 my_max = (my_max > rowsum? my_max : rowsum); 00839 } 00840 return my_max; 00841 } 00842 00843 00844 00845 template<typename T> 00846 inline 00847 T DenseMatrix<T>::transpose (const unsigned int i, 00848 const unsigned int j) const 00849 { 00850 // Implement in terms of operator() 00851 return (*this)(j,i); 00852 } 00853 00854 00855 00856 00857 00858 // template<typename T> 00859 // inline 00860 // void DenseMatrix<T>::condense(const unsigned int iv, 00861 // const unsigned int jv, 00862 // const T val, 00863 // DenseVector<T>& rhs) 00864 // { 00865 // libmesh_assert_equal_to (this->_m, rhs.size()); 00866 // libmesh_assert_equal_to (iv, jv); 00867 00868 00869 // // move the known value into the RHS 00870 // // and zero the column 00871 // for (unsigned int i=0; i<this->m(); i++) 00872 // { 00873 // rhs(i) -= ((*this)(i,jv))*val; 00874 // (*this)(i,jv) = 0.; 00875 // } 00876 00877 // // zero the row 00878 // for (unsigned int j=0; j<this->n(); j++) 00879 // (*this)(iv,j) = 0.; 00880 00881 // (*this)(iv,jv) = 1.; 00882 // rhs(iv) = val; 00883 00884 // } 00885 00886 00887 00888 00889 } // namespace libMesh 00890 00891 00892 00893 00894 #endif // LIBMESH_DENSE_MATRIX_H 00895
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:46 UTC
Hosted By: