- Home
- About DUNE
- Download
- Documentation
- Community
- Development
00001 // -*- tab-width: 8; indent-tabs-mode: nil; c-basic-offset: 2 -*- 00002 // vi: set et ts=8 sw=2 sts=2: 00003 // $Id: fmatrix.hh 6128 2010-09-08 13:50:00Z christi $ 00004 #ifndef DUNE_DENSEMATRIX_HH 00005 #define DUNE_DENSEMATRIX_HH 00006 00007 #include <cmath> 00008 #include <cstddef> 00009 #include <iostream> 00010 #include <vector> 00011 00012 #include <dune/common/misc.hh> 00013 #include <dune/common/exceptions.hh> 00014 #include <dune/common/fvector.hh> 00015 #include <dune/common/precision.hh> 00016 #include <dune/common/static_assert.hh> 00017 #include <dune/common/classname.hh> 00018 00019 namespace Dune 00020 { 00021 00022 template<typename M> class DenseMatrix; 00023 00024 template<typename M> 00025 struct FieldTraits< DenseMatrix<M> > 00026 { 00027 typedef const typename FieldTraits< typename DenseMatVecTraits<M>::value_type >::field_type field_type; 00028 typedef const typename FieldTraits< typename DenseMatVecTraits<M>::value_type >::real_type real_type; 00029 }; 00030 00031 /* 00032 work around a problem of FieldMatrix/FieldVector, 00033 there is no unique way to obtain the size of a class 00034 */ 00035 template<class K, int N, int M> class FieldMatrix; 00036 template<class K, int N> class FieldVector; 00037 namespace { 00038 template<class V> 00039 struct VectorSize 00040 { 00041 static typename V::size_type size(const V & v) { return v.size(); } 00042 }; 00043 00044 template<class K, int N> 00045 struct VectorSize< const FieldVector<K,N> > 00046 { 00047 typedef FieldVector<K,N> V; 00048 static typename V::size_type size(const V & v) { return N; } 00049 }; 00050 }; 00051 00067 template<typename M, typename T> 00068 void istl_assign_to_fmatrix(DenseMatrix<M>& f, const T& t) 00069 { 00070 DUNE_THROW(NotImplemented, "You need to specialise the method istl_assign_to_fmatrix(DenseMatrix<M>& f, const T& t) " 00071 << "(with M being " << className<M>() << ") " 00072 << "for T == " << className<T>() << "!"); 00073 } 00074 00075 namespace 00076 { 00077 template<bool b> 00078 struct DenseMatrixAssigner 00079 { 00080 template<typename M, typename T> 00081 static void assign(DenseMatrix<M>& fm, const T& t) 00082 { 00083 istl_assign_to_fmatrix(fm, t); 00084 } 00085 00086 }; 00087 00088 00089 template<> 00090 struct DenseMatrixAssigner<true> 00091 { 00092 template<typename M, typename T> 00093 static void assign(DenseMatrix<M>& fm, const T& t) 00094 { 00095 fm = static_cast<const typename DenseMatVecTraits<M>::value_type>(t); 00096 } 00097 }; 00098 } 00099 00101 class FMatrixError : public Exception {}; 00102 00113 template<typename MAT> 00114 class DenseMatrix 00115 { 00116 typedef DenseMatVecTraits<MAT> Traits; 00117 00118 // Curiously recurring template pattern 00119 MAT & asImp() { return static_cast<MAT&>(*this); } 00120 const MAT & asImp() const { return static_cast<const MAT&>(*this); } 00121 00122 public: 00123 //===== type definitions and constants 00124 00126 typedef typename Traits::derived_type derived_type; 00127 00129 typedef typename Traits::value_type value_type; 00130 00132 typedef typename Traits::value_type field_type; 00133 00135 typedef typename Traits::value_type block_type; 00136 00138 typedef typename Traits::size_type size_type; 00139 00141 typedef typename Traits::row_type row_type; 00142 00144 enum { 00146 blocklevel = 1 00147 }; 00148 00149 //===== access to components 00150 00152 row_type & operator[] (size_type i) 00153 { 00154 return asImp().mat_access(i); 00155 } 00156 00157 const row_type & operator[] (size_type i) const 00158 { 00159 return asImp().mat_access(i); 00160 } 00161 00163 size_type size() const 00164 { 00165 return rows(); 00166 } 00167 00168 //===== iterator interface to rows of the matrix 00170 typedef DenseIterator<DenseMatrix,row_type> Iterator; 00172 typedef Iterator iterator; 00174 typedef Iterator RowIterator; 00176 typedef typename row_type::Iterator ColIterator; 00177 00179 Iterator begin () 00180 { 00181 return Iterator(*this,0); 00182 } 00183 00185 Iterator end () 00186 { 00187 return Iterator(*this,rows()); 00188 } 00189 00194 Iterator rbegin() DUNE_DEPRECATED 00195 { 00196 return beforeBegin(); 00197 } 00198 00201 Iterator beforeEnd () 00202 { 00203 return Iterator(*this,rows()-1); 00204 } 00205 00210 Iterator rend () DUNE_DEPRECATED 00211 { 00212 return beforeBegin(); 00213 } 00214 00217 Iterator beforeBegin () 00218 { 00219 return Iterator(*this,-1); 00220 } 00221 00223 typedef DenseIterator<const DenseMatrix,const row_type> ConstIterator; 00225 typedef ConstIterator const_iterator; 00227 typedef ConstIterator ConstRowIterator; 00229 typedef typename row_type::ConstIterator ConstColIterator; 00230 00232 ConstIterator begin () const 00233 { 00234 return ConstIterator(*this,0); 00235 } 00236 00238 ConstIterator end () const 00239 { 00240 return ConstIterator(*this,rows()); 00241 } 00242 00247 ConstIterator rbegin() const DUNE_DEPRECATED 00248 { 00249 return beforeEnd(); 00250 } 00251 00254 ConstIterator beforeEnd () const 00255 { 00256 return ConstIterator(*this,rows()-1); 00257 } 00258 00263 ConstIterator rend () const DUNE_DEPRECATED 00264 { 00265 return beforeBegin(); 00266 } 00267 00270 ConstIterator beforeBegin () const 00271 { 00272 return ConstIterator(*this,-1); 00273 } 00274 00275 //===== assignment from scalar 00276 DenseMatrix& operator= (const field_type& f) 00277 { 00278 for (size_type i=0; i<rows(); i++) 00279 (*this)[i] = f; 00280 return *this; 00281 } 00282 00283 template<typename T> 00284 DenseMatrix& operator= (const T& t) 00285 { 00286 DenseMatrixAssigner<Conversion<T,field_type>::exists>::assign(*this, t); 00287 return *this; 00288 } 00289 //===== vector space arithmetic 00290 00292 DenseMatrix& operator+= (const DenseMatrix& y) 00293 { 00294 for (size_type i=0; i<rows(); i++) 00295 (*this)[i] += y[i]; 00296 return *this; 00297 } 00298 00300 DenseMatrix& operator-= (const DenseMatrix& y) 00301 { 00302 for (size_type i=0; i<rows(); i++) 00303 (*this)[i] -= y[i]; 00304 return *this; 00305 } 00306 00308 DenseMatrix& operator*= (const field_type& k) 00309 { 00310 for (size_type i=0; i<rows(); i++) 00311 (*this)[i] *= k; 00312 return *this; 00313 } 00314 00316 DenseMatrix& operator/= (const field_type& k) 00317 { 00318 for (size_type i=0; i<rows(); i++) 00319 (*this)[i] /= k; 00320 return *this; 00321 } 00322 00324 DenseMatrix &axpy (const field_type &k, const DenseMatrix &y ) 00325 { 00326 for( size_type i = 0; i < rows(); ++i ) 00327 (*this)[ i ].axpy( k, y[ i ] ); 00328 return *this; 00329 } 00330 00332 bool operator== (const DenseMatrix& y) const 00333 { 00334 for (size_type i=0; i<rows(); i++) 00335 if ((*this)[i]!=y[i]) 00336 return false; 00337 return true; 00338 } 00340 bool operator!= (const DenseMatrix& y) const 00341 { 00342 return !operator==(y); 00343 } 00344 00345 00346 //===== linear maps 00347 00349 template<class X, class Y> 00350 void mv (const X& x, Y& y) const 00351 { 00352 #ifdef DUNE_FMatrix_WITH_CHECKING 00353 assert(&x != &y); 00354 if (x.N()!=M()) DUNE_THROW(FMatrixError,"Index out of range"); 00355 if (y.N()!=N()) DUNE_THROW(FMatrixError,"Index out of range"); 00356 #endif 00357 for (size_type i=0; i<rows(); ++i) 00358 { 00359 y[i] = 0; 00360 for (size_type j=0; j<cols(); j++) 00361 y[i] += (*this)[i][j] * x[j]; 00362 } 00363 } 00364 00366 template< class X, class Y > 00367 void mtv ( const X &x, Y &y ) const 00368 { 00369 #ifdef DUNE_FMatrix_WITH_CHECKING 00370 //assert( &x != &y ); 00371 //This assert did not work for me. Compile error: 00372 // comparison between distinct pointer types ‘const 00373 // Dune::FieldVector<double, 3>*’ and ‘Dune::FieldVector<double, 2>*’ lacks a cast 00374 if( x.N() != N() ) 00375 DUNE_THROW( FMatrixError, "Index out of range." ); 00376 if( y.N() != M() ) 00377 DUNE_THROW( FMatrixError, "Index out of range." ); 00378 #endif 00379 for( size_type i = 0; i < cols(); ++i ) 00380 { 00381 y[ i ] = 0; 00382 for( size_type j = 0; j < rows(); ++j ) 00383 y[ i ] += (*this)[ j ][ i ] * x[ j ]; 00384 } 00385 } 00386 00388 template<class X, class Y> 00389 void umv (const X& x, Y& y) const 00390 { 00391 #ifdef DUNE_FMatrix_WITH_CHECKING 00392 if (x.N()!=M()) DUNE_THROW(FMatrixError,"index out of range"); 00393 if (y.N()!=N()) DUNE_THROW(FMatrixError,"index out of range"); 00394 #endif 00395 for (size_type i=0; i<rows(); i++) 00396 for (size_type j=0; j<cols(); j++) 00397 y[i] += (*this)[i][j] * x[j]; 00398 } 00399 00401 template<class X, class Y> 00402 void umtv (const X& x, Y& y) const 00403 { 00404 #ifdef DUNE_FMatrix_WITH_CHECKING 00405 if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range"); 00406 if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range"); 00407 #endif 00408 00409 for (size_type i=0; i<rows(); i++) 00410 for (size_type j=0; j<cols(); j++) 00411 y[j] += (*this)[i][j]*x[i]; 00412 } 00413 00415 template<class X, class Y> 00416 void umhv (const X& x, Y& y) const 00417 { 00418 #ifdef DUNE_FMatrix_WITH_CHECKING 00419 if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range"); 00420 if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range"); 00421 #endif 00422 00423 for (size_type i=0; i<rows(); i++) 00424 for (size_type j=0; j<cols(); j++) 00425 y[j] += conjugateComplex((*this)[i][j])*x[i]; 00426 } 00427 00429 template<class X, class Y> 00430 void mmv (const X& x, Y& y) const 00431 { 00432 #ifdef DUNE_FMatrix_WITH_CHECKING 00433 if (x.N()!=M()) DUNE_THROW(FMatrixError,"index out of range"); 00434 if (y.N()!=N()) DUNE_THROW(FMatrixError,"index out of range"); 00435 #endif 00436 for (size_type i=0; i<rows(); i++) 00437 for (size_type j=0; j<cols(); j++) 00438 y[i] -= (*this)[i][j] * x[j]; 00439 } 00440 00442 template<class X, class Y> 00443 void mmtv (const X& x, Y& y) const 00444 { 00445 #ifdef DUNE_FMatrix_WITH_CHECKING 00446 if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range"); 00447 if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range"); 00448 #endif 00449 00450 for (size_type i=0; i<rows(); i++) 00451 for (size_type j=0; j<cols(); j++) 00452 y[j] -= (*this)[i][j]*x[i]; 00453 } 00454 00456 template<class X, class Y> 00457 void mmhv (const X& x, Y& y) const 00458 { 00459 #ifdef DUNE_FMatrix_WITH_CHECKING 00460 if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range"); 00461 if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range"); 00462 #endif 00463 00464 for (size_type i=0; i<rows(); i++) 00465 for (size_type j=0; j<cols(); j++) 00466 y[j] -= conjugateComplex((*this)[i][j])*x[i]; 00467 } 00468 00470 template<class X, class Y> 00471 void usmv (const field_type& alpha, const X& x, Y& y) const 00472 { 00473 #ifdef DUNE_FMatrix_WITH_CHECKING 00474 if (x.N()!=M()) DUNE_THROW(FMatrixError,"index out of range"); 00475 if (y.N()!=N()) DUNE_THROW(FMatrixError,"index out of range"); 00476 #endif 00477 for (size_type i=0; i<rows(); i++) 00478 for (size_type j=0; j<cols(); j++) 00479 y[i] += alpha * (*this)[i][j] * x[j]; 00480 } 00481 00483 template<class X, class Y> 00484 void usmtv (const field_type& alpha, const X& x, Y& y) const 00485 { 00486 #ifdef DUNE_FMatrix_WITH_CHECKING 00487 if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range"); 00488 if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range"); 00489 #endif 00490 00491 for (size_type i=0; i<rows(); i++) 00492 for (size_type j=0; j<cols(); j++) 00493 y[j] += alpha*(*this)[i][j]*x[i]; 00494 } 00495 00497 template<class X, class Y> 00498 void usmhv (const field_type& alpha, const X& x, Y& y) const 00499 { 00500 #ifdef DUNE_FMatrix_WITH_CHECKING 00501 if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range"); 00502 if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range"); 00503 #endif 00504 00505 for (size_type i=0; i<rows(); i++) 00506 for (size_type j=0; j<cols(); j++) 00507 y[j] += alpha*conjugateComplex((*this)[i][j])*x[i]; 00508 } 00509 00510 //===== norms 00511 00513 typename FieldTraits<value_type>::real_type frobenius_norm () const 00514 { 00515 typename FieldTraits<value_type>::real_type sum=(0.0); 00516 for (size_type i=0; i<rows(); ++i) sum += (*this)[i].two_norm2(); 00517 return fvmeta::sqrt(sum); 00518 } 00519 00521 typename FieldTraits<value_type>::real_type frobenius_norm2 () const 00522 { 00523 typename FieldTraits<value_type>::real_type sum=(0.0); 00524 for (size_type i=0; i<rows(); ++i) sum += (*this)[i].two_norm2(); 00525 return sum; 00526 } 00527 00529 typename FieldTraits<value_type>::real_type infinity_norm () const 00530 { 00531 typename remove_const< typename FieldTraits<value_type>::real_type >::type max=(0.0); 00532 for (size_type i=0; i<rows(); ++i) max = std::max(max,(*this)[i].one_norm()); 00533 return max; 00534 } 00535 00537 typename FieldTraits<value_type>::real_type infinity_norm_real () const 00538 { 00539 typename FieldTraits<value_type>::real_type max(0.0); 00540 for (size_type i=0; i<rows(); ++i) max = std::max(max,(*this)[i].one_norm_real()); 00541 return max; 00542 } 00543 00544 //===== solve 00545 00550 template <class V> 00551 void solve (V& x, const V& b) const; 00552 00557 void invert(); 00558 00560 field_type determinant () const; 00561 00563 template<typename M2> 00564 MAT& leftmultiply (const DenseMatrix<M2>& M) 00565 { 00566 assert(M.rows() == M.cols() && M.rows() == rows()); 00567 MAT C(asImp()); 00568 00569 for (size_type i=0; i<rows(); i++) 00570 for (size_type j=0; j<cols(); j++) { 00571 (*this)[i][j] = 0; 00572 for (size_type k=0; k<rows(); k++) 00573 (*this)[i][j] += M[i][k]*C[k][j]; 00574 } 00575 00576 return asImp(); 00577 } 00578 00580 template<typename M2> 00581 MAT& rightmultiply (const DenseMatrix<M2>& M) 00582 { 00583 assert(M.rows() == M.cols() && M.cols() == cols()); 00584 MAT C(asImp()); 00585 00586 for (size_type i=0; i<rows(); i++) 00587 for (size_type j=0; j<cols(); j++) { 00588 (*this)[i][j] = 0; 00589 for (size_type k=0; k<cols(); k++) 00590 (*this)[i][j] += C[i][k]*M[k][j]; 00591 } 00592 return asImp(); 00593 } 00594 00595 #if 0 00596 00597 template<int l> 00598 DenseMatrix<K,l,cols> leftmultiplyany (const FieldMatrix<K,l,rows>& M) const 00599 { 00600 FieldMatrix<K,l,cols> C; 00601 00602 for (size_type i=0; i<l; i++) { 00603 for (size_type j=0; j<cols(); j++) { 00604 C[i][j] = 0; 00605 for (size_type k=0; k<rows(); k++) 00606 C[i][j] += M[i][k]*(*this)[k][j]; 00607 } 00608 } 00609 return C; 00610 } 00611 00613 template<int l> 00614 FieldMatrix<K,rows,l> rightmultiplyany (const FieldMatrix<K,cols,l>& M) const 00615 { 00616 FieldMatrix<K,rows,l> C; 00617 00618 for (size_type i=0; i<rows(); i++) { 00619 for (size_type j=0; j<l; j++) { 00620 C[i][j] = 0; 00621 for (size_type k=0; k<cols(); k++) 00622 C[i][j] += (*this)[i][k]*M[k][j]; 00623 } 00624 } 00625 return C; 00626 } 00627 #endif 00628 00629 //===== sizes 00630 00632 size_type N () const 00633 { 00634 return rows(); 00635 } 00636 00638 size_type M () const 00639 { 00640 return cols(); 00641 } 00642 00644 size_type rows() const 00645 { 00646 return asImp().mat_rows(); 00647 } 00648 00650 size_type cols() const 00651 { 00652 return asImp().mat_cols(); 00653 } 00654 00655 //===== query 00656 00658 bool exists (size_type i, size_type j) const 00659 { 00660 #ifdef DUNE_FMatrix_WITH_CHECKING 00661 if (i<0 || i>=rows()) DUNE_THROW(FMatrixError,"row index out of range"); 00662 if (j<0 || j>=cols()) DUNE_THROW(FMatrixError,"column index out of range"); 00663 #endif 00664 return true; 00665 } 00666 00667 private: 00668 00669 #ifndef DOXYGEN 00670 struct ElimPivot 00671 { 00672 ElimPivot(std::vector<size_type> & pivot); 00673 00674 void swap(int i, int j); 00675 00676 template<typename T> 00677 void operator()(const T&, int k, int i) 00678 {} 00679 00680 std::vector<size_type> & pivot_; 00681 }; 00682 00683 template<typename V> 00684 struct Elim 00685 { 00686 Elim(V& rhs); 00687 00688 void swap(int i, int j); 00689 00690 void operator()(const typename V::field_type& factor, int k, int i); 00691 00692 V* rhs_; 00693 }; 00694 00695 struct ElimDet 00696 { 00697 ElimDet(field_type& sign) : sign_(sign) 00698 { sign_ = 1; } 00699 00700 void swap(int i, int j) 00701 { sign_ *= -1; } 00702 00703 void operator()(const field_type&, int k, int i) 00704 {} 00705 00706 field_type& sign_; 00707 }; 00708 #endif // DOXYGEN 00709 00710 template<class Func> 00711 void luDecomposition(DenseMatrix<MAT>& A, Func func) const; 00712 }; 00713 00714 #ifndef DOXYGEN 00715 template<typename MAT> 00716 DenseMatrix<MAT>::ElimPivot::ElimPivot(std::vector<size_type> & pivot) 00717 : pivot_(pivot) 00718 { 00719 typedef typename std::vector<size_type>::size_type size_type; 00720 for(size_type i=0; i < pivot_.size(); ++i) pivot_[i]=i; 00721 } 00722 00723 template<typename MAT> 00724 void DenseMatrix<MAT>::ElimPivot::swap(int i, int j) 00725 { 00726 pivot_[i]=j; 00727 } 00728 00729 template<typename MAT> 00730 template<typename V> 00731 DenseMatrix<MAT>::Elim<V>::Elim(V& rhs) 00732 : rhs_(&rhs) 00733 {} 00734 00735 template<typename MAT> 00736 template<typename V> 00737 void DenseMatrix<MAT>::Elim<V>::swap(int i, int j) 00738 { 00739 std::swap((*rhs_)[i], (*rhs_)[j]); 00740 } 00741 00742 template<typename MAT> 00743 template<typename V> 00744 void DenseMatrix<MAT>:: 00745 Elim<V>::operator()(const typename V::field_type& factor, int k, int i) 00746 { 00747 (*rhs_)[k] -= factor*(*rhs_)[i]; 00748 } 00749 template<typename MAT> 00750 template<typename Func> 00751 inline void DenseMatrix<MAT>::luDecomposition(DenseMatrix<MAT>& A, Func func) const 00752 { 00753 typedef typename FieldTraits<value_type>::real_type 00754 real_type; 00755 typename FieldTraits<value_type>::real_type norm = 00756 A.infinity_norm_real(); // for relative thresholds 00757 typename FieldTraits<value_type>::real_type pivthres = 00758 std::max(FMatrixPrecision<real_type>::absolute_limit(),norm*FMatrixPrecision<>::pivoting_limit()); 00759 typename FieldTraits<value_type>::real_type singthres = 00760 std::max(FMatrixPrecision<real_type>::absolute_limit(),norm*FMatrixPrecision<>::singular_limit()); 00761 00762 // LU decomposition of A in A 00763 for (size_type i=0; i<rows(); i++) // loop over all rows 00764 { 00765 typename FieldTraits<value_type>::real_type pivmax=fvmeta::absreal(A[i][i]); 00766 00767 // pivoting ? 00768 if (pivmax<pivthres) 00769 { 00770 // compute maximum of column 00771 size_type imax=i; 00772 typename FieldTraits<value_type>::real_type abs(0.0); 00773 for (size_type k=i+1; k<rows(); k++) 00774 if ((abs=fvmeta::absreal(A[k][i]))>pivmax) 00775 { 00776 pivmax = abs; imax = k; 00777 } 00778 // swap rows 00779 if (imax!=i){ 00780 for (size_type j=0; j<rows(); j++) 00781 std::swap(A[i][j],A[imax][j]); 00782 func.swap(i, imax); // swap the pivot or rhs 00783 } 00784 } 00785 00786 // singular ? 00787 if (pivmax<singthres) 00788 DUNE_THROW(FMatrixError,"matrix is singular"); 00789 00790 // eliminate 00791 for (size_type k=i+1; k<rows(); k++) 00792 { 00793 field_type factor = A[k][i]/A[i][i]; 00794 A[k][i] = factor; 00795 for (size_type j=i+1; j<rows(); j++) 00796 A[k][j] -= factor*A[i][j]; 00797 func(factor, k, i); 00798 } 00799 } 00800 } 00801 00802 template<typename MAT> 00803 template <class V> 00804 inline void DenseMatrix<MAT>::solve(V& x, const V& b) const 00805 { 00806 // never mind those ifs, because they get optimized away 00807 if (rows()!=cols()) 00808 DUNE_THROW(FMatrixError, "Can't solve for a " << rows() << "x" << cols() << " matrix!"); 00809 00810 if (rows()==1) { 00811 00812 #ifdef DUNE_FMatrix_WITH_CHECKING 00813 if (fvmeta::absreal((*this)[0][0])<FMatrixPrecision<>::absolute_limit()) 00814 DUNE_THROW(FMatrixError,"matrix is singular"); 00815 #endif 00816 x[0] = b[0]/(*this)[0][0]; 00817 00818 } 00819 else if (rows()==2) { 00820 00821 field_type detinv = (*this)[0][0]*(*this)[1][1]-(*this)[0][1]*(*this)[1][0]; 00822 #ifdef DUNE_FMatrix_WITH_CHECKING 00823 if (fvmeta::absreal(detinv)<FMatrixPrecision<>::absolute_limit()) 00824 DUNE_THROW(FMatrixError,"matrix is singular"); 00825 #endif 00826 detinv = 1.0/detinv; 00827 00828 x[0] = detinv*((*this)[1][1]*b[0]-(*this)[0][1]*b[1]); 00829 x[1] = detinv*((*this)[0][0]*b[1]-(*this)[1][0]*b[0]); 00830 00831 } 00832 else if (rows()==3) { 00833 00834 field_type d = determinant(); 00835 #ifdef DUNE_FMatrix_WITH_CHECKING 00836 if (fvmeta::absreal(d)<FMatrixPrecision<>::absolute_limit()) 00837 DUNE_THROW(FMatrixError,"matrix is singular"); 00838 #endif 00839 00840 x[0] = (b[0]*(*this)[1][1]*(*this)[2][2] - b[0]*(*this)[2][1]*(*this)[1][2] 00841 - b[1] *(*this)[0][1]*(*this)[2][2] + b[1]*(*this)[2][1]*(*this)[0][2] 00842 + b[2] *(*this)[0][1]*(*this)[1][2] - b[2]*(*this)[1][1]*(*this)[0][2]) / d; 00843 00844 x[1] = ((*this)[0][0]*b[1]*(*this)[2][2] - (*this)[0][0]*b[2]*(*this)[1][2] 00845 - (*this)[1][0] *b[0]*(*this)[2][2] + (*this)[1][0]*b[2]*(*this)[0][2] 00846 + (*this)[2][0] *b[0]*(*this)[1][2] - (*this)[2][0]*b[1]*(*this)[0][2]) / d; 00847 00848 x[2] = ((*this)[0][0]*(*this)[1][1]*b[2] - (*this)[0][0]*(*this)[2][1]*b[1] 00849 - (*this)[1][0] *(*this)[0][1]*b[2] + (*this)[1][0]*(*this)[2][1]*b[0] 00850 + (*this)[2][0] *(*this)[0][1]*b[1] - (*this)[2][0]*(*this)[1][1]*b[0]) / d; 00851 00852 } 00853 else { 00854 00855 V& rhs = x; // use x to store rhs 00856 rhs = b; // copy data 00857 Elim<V> elim(rhs); 00858 MAT A(asImp()); 00859 00860 luDecomposition(A, elim); 00861 00862 // backsolve 00863 for(int i=rows()-1; i>=0; i--){ 00864 for (size_type j=i+1; j<rows(); j++) 00865 rhs[i] -= A[i][j]*x[j]; 00866 x[i] = rhs[i]/A[i][i]; 00867 } 00868 } 00869 } 00870 00871 template<typename MAT> 00872 inline void DenseMatrix<MAT>::invert() 00873 { 00874 // never mind those ifs, because they get optimized away 00875 if (rows()!=cols()) 00876 DUNE_THROW(FMatrixError, "Can't invert a " << rows() << "x" << cols() << " matrix!"); 00877 00878 if (rows()==1) { 00879 00880 #ifdef DUNE_FMatrix_WITH_CHECKING 00881 if (fvmeta::absreal((*this)[0][0])<FMatrixPrecision<>::absolute_limit()) 00882 DUNE_THROW(FMatrixError,"matrix is singular"); 00883 #endif 00884 (*this)[0][0] = 1.0/(*this)[0][0]; 00885 00886 } 00887 else if (rows()==2) { 00888 00889 field_type detinv = (*this)[0][0]*(*this)[1][1]-(*this)[0][1]*(*this)[1][0]; 00890 #ifdef DUNE_FMatrix_WITH_CHECKING 00891 if (fvmeta::absreal(detinv)<FMatrixPrecision<>::absolute_limit()) 00892 DUNE_THROW(FMatrixError,"matrix is singular"); 00893 #endif 00894 detinv = 1.0/detinv; 00895 00896 field_type temp=(*this)[0][0]; 00897 (*this)[0][0] = (*this)[1][1]*detinv; 00898 (*this)[0][1] = -(*this)[0][1]*detinv; 00899 (*this)[1][0] = -(*this)[1][0]*detinv; 00900 (*this)[1][1] = temp*detinv; 00901 00902 } 00903 else { 00904 00905 MAT A(asImp()); 00906 std::vector<size_type> pivot(rows()); 00907 luDecomposition(A, ElimPivot(pivot)); 00908 DenseMatrix<MAT>& L=A; 00909 DenseMatrix<MAT>& U=A; 00910 00911 // initialize inverse 00912 *this=field_type(); 00913 00914 for(size_type i=0; i<rows(); ++i) 00915 (*this)[i][i]=1; 00916 00917 // L Y = I; multiple right hand sides 00918 for (size_type i=0; i<rows(); i++) 00919 for (size_type j=0; j<i; j++) 00920 for (size_type k=0; k<rows(); k++) 00921 (*this)[i][k] -= L[i][j]*(*this)[j][k]; 00922 00923 // U A^{-1} = Y 00924 for (size_type i=rows(); i>0;){ 00925 --i; 00926 for (size_type k=0; k<rows(); k++){ 00927 for (size_type j=i+1; j<rows(); j++) 00928 (*this)[i][k] -= U[i][j]*(*this)[j][k]; 00929 (*this)[i][k] /= U[i][i]; 00930 } 00931 } 00932 00933 for(size_type i=rows(); i>0; ){ 00934 --i; 00935 if(i!=pivot[i]) 00936 for(size_type j=0; j<rows(); ++j) 00937 std::swap((*this)[j][pivot[i]], (*this)[j][i]); 00938 } 00939 } 00940 } 00941 00942 // implementation of the determinant 00943 template<typename MAT> 00944 inline typename DenseMatrix<MAT>::field_type 00945 DenseMatrix<MAT>::determinant() const 00946 { 00947 // never mind those ifs, because they get optimized away 00948 if (rows()!=cols()) 00949 DUNE_THROW(FMatrixError, "There is no determinant for a " << rows() << "x" << cols() << " matrix!"); 00950 00951 if (rows()==1) 00952 return (*this)[0][0]; 00953 00954 if (rows()==2) 00955 return (*this)[0][0]*(*this)[1][1] - (*this)[0][1]*(*this)[1][0]; 00956 00957 if (rows()==3) { 00958 // code generated by maple 00959 field_type t4 = (*this)[0][0] * (*this)[1][1]; 00960 field_type t6 = (*this)[0][0] * (*this)[1][2]; 00961 field_type t8 = (*this)[0][1] * (*this)[1][0]; 00962 field_type t10 = (*this)[0][2] * (*this)[1][0]; 00963 field_type t12 = (*this)[0][1] * (*this)[2][0]; 00964 field_type t14 = (*this)[0][2] * (*this)[2][0]; 00965 00966 return (t4*(*this)[2][2]-t6*(*this)[2][1]-t8*(*this)[2][2]+ 00967 t10*(*this)[2][1]+t12*(*this)[1][2]-t14*(*this)[1][1]); 00968 00969 } 00970 00971 MAT A(asImp()); 00972 field_type det; 00973 try 00974 { 00975 luDecomposition(A, ElimDet(det)); 00976 } 00977 catch (FMatrixError&) 00978 { 00979 return 0; 00980 } 00981 for (size_type i = 0; i < rows(); ++i) 00982 det *= A[i][i]; 00983 return det; 00984 } 00985 00986 #endif // DOXYGEN 00987 00988 namespace DenseMatrixHelp { 00989 #if 0 00990 00991 template <typename K> 00992 static inline K invertMatrix (const FieldMatrix<K,1,1> &matrix, FieldMatrix<K,1,1> &inverse) 00993 { 00994 inverse[0][0] = 1.0/matrix[0][0]; 00995 return matrix[0][0]; 00996 } 00997 00999 template <typename K> 01000 static inline K invertMatrix_retTransposed (const FieldMatrix<K,1,1> &matrix, FieldMatrix<K,1,1> &inverse) 01001 { 01002 return invertMatrix(matrix,inverse); 01003 } 01004 01005 01007 template <typename K> 01008 static inline K invertMatrix (const FieldMatrix<K,2,2> &matrix, FieldMatrix<K,2,2> &inverse) 01009 { 01010 // code generated by maple 01011 field_type det = (matrix[0][0]*matrix[1][1] - matrix[0][1]*matrix[1][0]); 01012 field_type det_1 = 1.0/det; 01013 inverse[0][0] = matrix[1][1] * det_1; 01014 inverse[0][1] = - matrix[0][1] * det_1; 01015 inverse[1][0] = - matrix[1][0] * det_1; 01016 inverse[1][1] = matrix[0][0] * det_1; 01017 return det; 01018 } 01019 01022 template <typename K> 01023 static inline K invertMatrix_retTransposed (const FieldMatrix<K,2,2> &matrix, FieldMatrix<K,2,2> &inverse) 01024 { 01025 // code generated by maple 01026 field_type det = (matrix[0][0]*matrix[1][1] - matrix[0][1]*matrix[1][0]); 01027 field_type det_1 = 1.0/det; 01028 inverse[0][0] = matrix[1][1] * det_1; 01029 inverse[1][0] = - matrix[0][1] * det_1; 01030 inverse[0][1] = - matrix[1][0] * det_1; 01031 inverse[1][1] = matrix[0][0] * det_1; 01032 return det; 01033 } 01034 01036 template <typename K> 01037 static inline K invertMatrix (const FieldMatrix<K,3,3> &matrix, FieldMatrix<K,3,3> &inverse) 01038 { 01039 // code generated by maple 01040 field_type t4 = matrix[0][0] * matrix[1][1]; 01041 field_type t6 = matrix[0][0] * matrix[1][2]; 01042 field_type t8 = matrix[0][1] * matrix[1][0]; 01043 field_type t10 = matrix[0][2] * matrix[1][0]; 01044 field_type t12 = matrix[0][1] * matrix[2][0]; 01045 field_type t14 = matrix[0][2] * matrix[2][0]; 01046 01047 field_type det = (t4*matrix[2][2]-t6*matrix[2][1]-t8*matrix[2][2]+ 01048 t10*matrix[2][1]+t12*matrix[1][2]-t14*matrix[1][1]); 01049 field_type t17 = 1.0/det; 01050 01051 inverse[0][0] = (matrix[1][1] * matrix[2][2] - matrix[1][2] * matrix[2][1])*t17; 01052 inverse[0][1] = -(matrix[0][1] * matrix[2][2] - matrix[0][2] * matrix[2][1])*t17; 01053 inverse[0][2] = (matrix[0][1] * matrix[1][2] - matrix[0][2] * matrix[1][1])*t17; 01054 inverse[1][0] = -(matrix[1][0] * matrix[2][2] - matrix[1][2] * matrix[2][0])*t17; 01055 inverse[1][1] = (matrix[0][0] * matrix[2][2] - t14) * t17; 01056 inverse[1][2] = -(t6-t10) * t17; 01057 inverse[2][0] = (matrix[1][0] * matrix[2][1] - matrix[1][1] * matrix[2][0]) * t17; 01058 inverse[2][1] = -(matrix[0][0] * matrix[2][1] - t12) * t17; 01059 inverse[2][2] = (t4-t8) * t17; 01060 01061 return det; 01062 } 01063 01065 template <typename K> 01066 static inline K invertMatrix_retTransposed (const FieldMatrix<K,3,3> &matrix, FieldMatrix<K,3,3> &inverse) 01067 { 01068 // code generated by maple 01069 field_type t4 = matrix[0][0] * matrix[1][1]; 01070 field_type t6 = matrix[0][0] * matrix[1][2]; 01071 field_type t8 = matrix[0][1] * matrix[1][0]; 01072 field_type t10 = matrix[0][2] * matrix[1][0]; 01073 field_type t12 = matrix[0][1] * matrix[2][0]; 01074 field_type t14 = matrix[0][2] * matrix[2][0]; 01075 01076 field_type det = (t4*matrix[2][2]-t6*matrix[2][1]-t8*matrix[2][2]+ 01077 t10*matrix[2][1]+t12*matrix[1][2]-t14*matrix[1][1]); 01078 field_type t17 = 1.0/det; 01079 01080 inverse[0][0] = (matrix[1][1] * matrix[2][2] - matrix[1][2] * matrix[2][1])*t17; 01081 inverse[1][0] = -(matrix[0][1] * matrix[2][2] - matrix[0][2] * matrix[2][1])*t17; 01082 inverse[2][0] = (matrix[0][1] * matrix[1][2] - matrix[0][2] * matrix[1][1])*t17; 01083 inverse[0][1] = -(matrix[1][0] * matrix[2][2] - matrix[1][2] * matrix[2][0])*t17; 01084 inverse[1][1] = (matrix[0][0] * matrix[2][2] - t14) * t17; 01085 inverse[2][1] = -(t6-t10) * t17; 01086 inverse[0][2] = (matrix[1][0] * matrix[2][1] - matrix[1][1] * matrix[2][0]) * t17; 01087 inverse[1][2] = -(matrix[0][0] * matrix[2][1] - t12) * t17; 01088 inverse[2][2] = (t4-t8) * t17; 01089 01090 return det; 01091 } 01092 01094 template< class K, int m, int n, int p > 01095 static inline void multMatrix ( const FieldMatrix< K, m, n > &A, 01096 const FieldMatrix< K, n, p > &B, 01097 FieldMatrix< K, m, p > &ret ) 01098 { 01099 typedef typename FieldMatrix< K, m, p > :: size_type size_type; 01100 01101 for( size_type i = 0; i < m; ++i ) 01102 { 01103 for( size_type j = 0; j < p; ++j ) 01104 { 01105 ret[ i ][ j ] = K( 0 ); 01106 for( size_type k = 0; k < n; ++k ) 01107 ret[ i ][ j ] += A[ i ][ k ] * B[ k ][ j ]; 01108 } 01109 } 01110 } 01111 01113 template <typename K, int rows, int cols> 01114 static inline void multTransposedMatrix(const FieldMatrix<K,rows,cols> &matrix, FieldMatrix<K,cols,cols>& ret) 01115 { 01116 typedef typename FieldMatrix<K,rows,cols>::size_type size_type; 01117 01118 for(size_type i=0; i<cols(); i++) 01119 for(size_type j=0; j<cols(); j++) 01120 { 01121 ret[i][j]=0.0; 01122 for(size_type k=0; k<rows(); k++) 01123 ret[i][j]+=matrix[k][i]*matrix[k][j]; 01124 } 01125 } 01126 #endif 01127 01129 template <typename MAT, typename V1, typename V2> 01130 static inline void multAssign(const DenseMatrix<MAT> &matrix, const DenseVector<V1> & x, DenseVector<V2> & ret) 01131 { 01132 assert(x.size() == matrix.cols()); 01133 assert(ret.size() == matrix.rows()); 01134 typedef typename DenseMatrix<MAT>::size_type size_type; 01135 01136 for(size_type i=0; i<matrix.rows(); ++i) 01137 { 01138 ret[i] = 0.0; 01139 for(size_type j=0; j<matrix.cols(); ++j) 01140 { 01141 ret[i] += matrix[i][j]*x[j]; 01142 } 01143 } 01144 } 01145 01146 #if 0 01147 01148 template <typename K, int rows, int cols> 01149 static inline void multAssignTransposed( const FieldMatrix<K,rows,cols> &matrix, const FieldVector<K,rows> & x, FieldVector<K,cols> & ret) 01150 { 01151 typedef typename FieldMatrix<K,rows,cols>::size_type size_type; 01152 01153 for(size_type i=0; i<cols(); ++i) 01154 { 01155 ret[i] = 0.0; 01156 for(size_type j=0; j<rows(); ++j) 01157 ret[i] += matrix[j][i]*x[j]; 01158 } 01159 } 01160 01162 template <typename K, int rows, int cols> 01163 static inline FieldVector<K,rows> mult(const FieldMatrix<K,rows,cols> &matrix, const FieldVector<K,cols> & x) 01164 { 01165 FieldVector<K,rows> ret; 01166 multAssign(matrix,x,ret); 01167 return ret; 01168 } 01169 01171 template <typename K, int rows, int cols> 01172 static inline FieldVector<K,cols> multTransposed(const FieldMatrix<K,rows,cols> &matrix, const FieldVector<K,rows> & x) 01173 { 01174 FieldVector<K,cols> ret; 01175 multAssignTransposed( matrix, x, ret ); 01176 return ret; 01177 } 01178 #endif 01179 01180 } // end namespace DenseMatrixHelp 01181 01183 template<typename MAT> 01184 std::ostream& operator<< (std::ostream& s, const DenseMatrix<MAT>& a) 01185 { 01186 for (typename DenseMatrix<MAT>::size_type i=0; i<a.rows(); i++) 01187 s << a[i] << std::endl; 01188 return s; 01189 } 01190 01193 } // end namespace Dune 01194 01195 #endif
Generated on Fri Apr 29 2011 with Doxygen (ver 1.7.1) [doxygen-log,error-log].