Dune Core Modules (2.3.1)

densematrix.hh
Go to the documentation of this file.
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3// $Id: fmatrix.hh 6128 2010-09-08 13:50:00Z christi $
4#ifndef DUNE_DENSEMATRIX_HH
5#define DUNE_DENSEMATRIX_HH
6
7#include <cmath>
8#include <cstddef>
9#include <iostream>
10#include <vector>
11
17#include <dune/common/math.hh>
18#include <dune/common/unused.hh>
19
20
21namespace Dune
22{
23
24 template<typename M> class DenseMatrix;
25
26 template<typename M>
27 struct FieldTraits< DenseMatrix<M> >
28 {
29 typedef const typename FieldTraits< typename DenseMatVecTraits<M>::value_type >::field_type field_type;
30 typedef const typename FieldTraits< typename DenseMatVecTraits<M>::value_type >::real_type real_type;
31 };
32
33 /*
34 work around a problem of FieldMatrix/FieldVector,
35 there is no unique way to obtain the size of a class
36 */
37 template<class K, int N, int M> class FieldMatrix;
38 template<class K, int N> class FieldVector;
39 namespace {
40 template<class V>
41 struct VectorSize
42 {
43 static typename V::size_type size(const V & v) { return v.size(); }
44 };
45
46 template<class K, int N>
47 struct VectorSize< const FieldVector<K,N> >
48 {
49 typedef FieldVector<K,N> V;
50 static typename V::size_type size(const V & v) { return N; }
51 };
52 }
53
72 template< class DenseMatrix, class RHS >
74
75
76
77 template< class DenseMatrix, class K, int N, int M >
78 void istl_assign_to_fmatrix ( DenseMatrix &denseMatrix, const K (&values)[ M ][ N ] )
79 {
80 for( int i = 0; i < N; ++i )
81 for( int j = 0; j < M; ++j )
82 denseMatrix[ i ][ j ] = values[ i ][ j ];
83 }
84
85
86
87#ifndef DOXYGEN
88 namespace
89 {
90 template< class DenseMatrix, class RHS,
92 class DenseMatrixAssignerImplementation;
93
94 template< class DenseMatrix, class RHS >
95 class DenseMatrixAssignerImplementation< DenseMatrix, RHS, true >
96 {
97 public:
98 static void apply ( DenseMatrix &denseMatrix, const RHS &rhs )
99 {
100 typedef typename DenseMatrix::field_type field_type;
101
102 const field_type value = static_cast< field_type >( rhs );
103
104 const std::size_t size = denseMatrix.size();
105 for( std::size_t i = 0; i < size; ++i )
106 denseMatrix[ i ] = value;
107 }
108 };
109
110 template< class DenseMatrix, class RHS >
111 class DenseMatrixAssignerImplementation< DenseMatrix, RHS, false >
112 {
113 template< class M, class T>
114 struct have_istl_assign_to_fmatrix
115 {
116 struct yes { char dummy[ 1 ]; };
117 struct no { char dummy[ 2 ]; };
118
119 template< class C>
120 static C &get_ref();
121
122 template< class C>
123 static yes test( decltype( istl_assign_to_fmatrix( get_ref< M >(), get_ref< C >() ) ) * );
124 template< class C >
125 static no test(...);
126
127 public:
128 static const bool v = sizeof( test< const T >( 0 ) ) == sizeof( yes );
129 };
130
131 template< class M, class T, bool = have_istl_assign_to_fmatrix< M, T >::v >
132 struct DefaultImplementation;
133
134 // forward to istl_assign_to_fmatrix()
135 template< class M, class T >
136 struct DefaultImplementation< M, T, true >
137 {
138 static void apply ( M &m, const T &t )
139 {
140 istl_assign_to_fmatrix( m, t );
141 }
142 };
143
144 // static_cast
145 template< class M, class T >
146 struct DefaultImplementation< M, T, false >
147 {
148 static void apply ( M &m, const T &t )
149 {
150 dune_static_assert( (Conversion< const T, const M >::exists), "No template specialization of DenseMatrixAssigner found" );
151 m = static_cast< const M & >( t );
152 }
153 };
154
155 public:
156 static void apply ( DenseMatrix &denseMatrix, const RHS &rhs )
157 {
158 DefaultImplementation< DenseMatrix, RHS >::apply( denseMatrix, rhs );
159 }
160 };
161 }
162
163
164
165 template< class DenseMatrix, class RHS >
166 struct DenseMatrixAssigner
167 {
168 static void apply ( DenseMatrix &denseMatrix, const RHS &rhs )
169 {
170 DenseMatrixAssignerImplementation< DenseMatrix, RHS >::apply( denseMatrix, rhs );
171 }
172 };
173#endif // #ifndef DOXYGEN
174
175
176
178 class FMatrixError : public Exception {};
179
190 template<typename MAT>
192 {
193 typedef DenseMatVecTraits<MAT> Traits;
194
195 // Curiously recurring template pattern
196 MAT & asImp() { return static_cast<MAT&>(*this); }
197 const MAT & asImp() const { return static_cast<const MAT&>(*this); }
198
199 public:
200 //===== type definitions and constants
201
203 typedef typename Traits::derived_type derived_type;
204
206 typedef typename Traits::value_type value_type;
207
209 typedef typename Traits::value_type field_type;
210
212 typedef typename Traits::value_type block_type;
213
215 typedef typename Traits::size_type size_type;
216
218 typedef typename Traits::row_type row_type;
219
221 typedef typename Traits::row_reference row_reference;
222
224 typedef typename Traits::const_row_reference const_row_reference;
225
227 enum {
229 blocklevel = 1
230 };
231
232 //===== access to components
233
236 {
237 return asImp().mat_access(i);
238 }
239
241 {
242 return asImp().mat_access(i);
243 }
244
247 {
248 return rows();
249 }
250
251 //===== iterator interface to rows of the matrix
259 typedef typename row_type::Iterator ColIterator;
260
263 {
264 return Iterator(*this,0);
265 }
266
269 {
270 return Iterator(*this,rows());
271 }
272
276 {
277 return Iterator(*this,rows()-1);
278 }
279
283 {
284 return Iterator(*this,-1);
285 }
286
294 typedef typename row_type::ConstIterator ConstColIterator;
295
298 {
299 return ConstIterator(*this,0);
300 }
301
304 {
305 return ConstIterator(*this,rows());
306 }
307
311 {
312 return ConstIterator(*this,rows()-1);
313 }
314
318 {
319 return ConstIterator(*this,-1);
320 }
321
322 //===== assignment
323
324 template< class RHS >
325 DenseMatrix &operator= ( const RHS &rhs )
326 {
328 return *this;
329 }
330
331 //===== vector space arithmetic
332
334 template <class Other>
336 {
337 for (size_type i=0; i<rows(); i++)
338 (*this)[i] += y[i];
339 return *this;
340 }
341
343 template <class Other>
345 {
346 for (size_type i=0; i<rows(); i++)
347 (*this)[i] -= y[i];
348 return *this;
349 }
350
353 {
354 for (size_type i=0; i<rows(); i++)
355 (*this)[i] *= k;
356 return *this;
357 }
358
361 {
362 for (size_type i=0; i<rows(); i++)
363 (*this)[i] /= k;
364 return *this;
365 }
366
368 template <class Other>
370 {
371 for( size_type i = 0; i < rows(); ++i )
372 (*this)[ i ].axpy( k, y[ i ] );
373 return *this;
374 }
375
377 template <class Other>
378 bool operator== (const DenseMatrix<Other>& y) const
379 {
380 for (size_type i=0; i<rows(); i++)
381 if ((*this)[i]!=y[i])
382 return false;
383 return true;
384 }
386 template <class Other>
387 bool operator!= (const DenseMatrix<Other>& y) const
388 {
389 return !operator==(y);
390 }
391
392
393 //===== linear maps
394
396 template<class X, class Y>
397 void mv (const X& x, Y& y) const
398 {
399#ifdef DUNE_FMatrix_WITH_CHECKING
400 if (x.N()!=M()) DUNE_THROW(FMatrixError,"Index out of range");
401 if (y.N()!=N()) DUNE_THROW(FMatrixError,"Index out of range");
402#endif
403 for (size_type i=0; i<rows(); ++i)
404 {
405 y[i] = 0;
406 for (size_type j=0; j<cols(); j++)
407 y[i] += (*this)[i][j] * x[j];
408 }
409 }
410
412 template< class X, class Y >
413 void mtv ( const X &x, Y &y ) const
414 {
415#ifdef DUNE_FMatrix_WITH_CHECKING
416 //assert( &x != &y );
417 //This assert did not work for me. Compile error:
418 // comparison between distinct pointer types ‘const
419 // Dune::FieldVector<double, 3>*’ and ‘Dune::FieldVector<double, 2>*’ lacks a cast
420 if( x.N() != N() )
421 DUNE_THROW( FMatrixError, "Index out of range." );
422 if( y.N() != M() )
423 DUNE_THROW( FMatrixError, "Index out of range." );
424#endif
425 for( size_type i = 0; i < cols(); ++i )
426 {
427 y[ i ] = 0;
428 for( size_type j = 0; j < rows(); ++j )
429 y[ i ] += (*this)[ j ][ i ] * x[ j ];
430 }
431 }
432
434 template<class X, class Y>
435 void umv (const X& x, Y& y) const
436 {
437#ifdef DUNE_FMatrix_WITH_CHECKING
438 if (x.N()!=M())
439 DUNE_THROW(FMatrixError,"y += A x -- index out of range (sizes: x: " << x.N() << ", y: " << y.N() << ", A: " << this->N() << " x " << this->M() << ")" << std::endl);
440 if (y.N()!=N())
441 DUNE_THROW(FMatrixError,"y += A x -- index out of range (sizes: x: " << x.N() << ", y: " << y.N() << ", A: " << this->N() << " x " << this->M() << ")" << std::endl);
442#endif
443 for (size_type i=0; i<rows(); i++)
444 for (size_type j=0; j<cols(); j++)
445 y[i] += (*this)[i][j] * x[j];
446 }
447
449 template<class X, class Y>
450 void umtv (const X& x, Y& y) const
451 {
452#ifdef DUNE_FMatrix_WITH_CHECKING
453 if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
454 if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
455#endif
456
457 for (size_type i=0; i<rows(); i++)
458 for (size_type j=0; j<cols(); j++)
459 y[j] += (*this)[i][j]*x[i];
460 }
461
463 template<class X, class Y>
464 void umhv (const X& x, Y& y) const
465 {
466#ifdef DUNE_FMatrix_WITH_CHECKING
467 if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
468 if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
469#endif
470
471 for (size_type i=0; i<rows(); i++)
472 for (size_type j=0; j<cols(); j++)
473 y[j] += conjugateComplex((*this)[i][j])*x[i];
474 }
475
477 template<class X, class Y>
478 void mmv (const X& x, Y& y) const
479 {
480#ifdef DUNE_FMatrix_WITH_CHECKING
481 if (x.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
482 if (y.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
483#endif
484 for (size_type i=0; i<rows(); i++)
485 for (size_type j=0; j<cols(); j++)
486 y[i] -= (*this)[i][j] * x[j];
487 }
488
490 template<class X, class Y>
491 void mmtv (const X& x, Y& y) const
492 {
493#ifdef DUNE_FMatrix_WITH_CHECKING
494 if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
495 if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
496#endif
497
498 for (size_type i=0; i<rows(); i++)
499 for (size_type j=0; j<cols(); j++)
500 y[j] -= (*this)[i][j]*x[i];
501 }
502
504 template<class X, class Y>
505 void mmhv (const X& x, Y& y) const
506 {
507#ifdef DUNE_FMatrix_WITH_CHECKING
508 if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
509 if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
510#endif
511
512 for (size_type i=0; i<rows(); i++)
513 for (size_type j=0; j<cols(); j++)
514 y[j] -= conjugateComplex((*this)[i][j])*x[i];
515 }
516
518 template<class X, class Y>
519 void usmv (const field_type& alpha, const X& x, Y& y) const
520 {
521#ifdef DUNE_FMatrix_WITH_CHECKING
522 if (x.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
523 if (y.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
524#endif
525 for (size_type i=0; i<rows(); i++)
526 for (size_type j=0; j<cols(); j++)
527 y[i] += alpha * (*this)[i][j] * x[j];
528 }
529
531 template<class X, class Y>
532 void usmtv (const field_type& alpha, const X& x, Y& y) const
533 {
534#ifdef DUNE_FMatrix_WITH_CHECKING
535 if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
536 if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
537#endif
538
539 for (size_type i=0; i<rows(); i++)
540 for (size_type j=0; j<cols(); j++)
541 y[j] += alpha*(*this)[i][j]*x[i];
542 }
543
545 template<class X, class Y>
546 void usmhv (const field_type& alpha, const X& x, Y& y) const
547 {
548#ifdef DUNE_FMatrix_WITH_CHECKING
549 if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
550 if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
551#endif
552
553 for (size_type i=0; i<rows(); i++)
554 for (size_type j=0; j<cols(); j++)
555 y[j] += alpha*conjugateComplex((*this)[i][j])*x[i];
556 }
557
558 //===== norms
559
561 typename FieldTraits<value_type>::real_type frobenius_norm () const
562 {
563 typename FieldTraits<value_type>::real_type sum=(0.0);
564 for (size_type i=0; i<rows(); ++i) sum += (*this)[i].two_norm2();
565 return fvmeta::sqrt(sum);
566 }
567
569 typename FieldTraits<value_type>::real_type frobenius_norm2 () const
570 {
571 typename FieldTraits<value_type>::real_type sum=(0.0);
572 for (size_type i=0; i<rows(); ++i) sum += (*this)[i].two_norm2();
573 return sum;
574 }
575
577 typename FieldTraits<value_type>::real_type infinity_norm () const
578 {
579 if (size() == 0)
580 return 0.0;
581
582 ConstIterator it = begin();
583 typename remove_const< typename FieldTraits<value_type>::real_type >::type max = it->one_norm();
584 for (it = it + 1; it != end(); ++it)
585 max = std::max(max, it->one_norm());
586
587 return max;
588 }
589
591 typename FieldTraits<value_type>::real_type infinity_norm_real () const
592 {
593 if (size() == 0)
594 return 0.0;
595
596 ConstIterator it = begin();
597 typename remove_const< typename FieldTraits<value_type>::real_type >::type max = it->one_norm_real();
598 for (it = it + 1; it != end(); ++it)
599 max = std::max(max, it->one_norm_real());
600
601 return max;
602 }
603
604 //===== solve
605
610 template <class V>
611 void solve (V& x, const V& b) const;
612
617 void invert();
618
621
623 template<typename M2>
625 {
626 assert(M.rows() == M.cols() && M.rows() == rows());
627 MAT C(asImp());
628
629 for (size_type i=0; i<rows(); i++)
630 for (size_type j=0; j<cols(); j++) {
631 (*this)[i][j] = 0;
632 for (size_type k=0; k<rows(); k++)
633 (*this)[i][j] += M[i][k]*C[k][j];
634 }
635
636 return asImp();
637 }
638
640 template<typename M2>
642 {
643 assert(M.rows() == M.cols() && M.cols() == cols());
644 MAT C(asImp());
645
646 for (size_type i=0; i<rows(); i++)
647 for (size_type j=0; j<cols(); j++) {
648 (*this)[i][j] = 0;
649 for (size_type k=0; k<cols(); k++)
650 (*this)[i][j] += C[i][k]*M[k][j];
651 }
652 return asImp();
653 }
654
655#if 0
657 template<int l>
658 DenseMatrix<K,l,cols> leftmultiplyany (const FieldMatrix<K,l,rows>& M) const
659 {
661
662 for (size_type i=0; i<l; i++) {
663 for (size_type j=0; j<cols(); j++) {
664 C[i][j] = 0;
665 for (size_type k=0; k<rows(); k++)
666 C[i][j] += M[i][k]*(*this)[k][j];
667 }
668 }
669 return C;
670 }
671
673 template<int l>
674 FieldMatrix<K,rows,l> rightmultiplyany (const FieldMatrix<K,cols,l>& M) const
675 {
676 FieldMatrix<K,rows,l> C;
677
678 for (size_type i=0; i<rows(); i++) {
679 for (size_type j=0; j<l; j++) {
680 C[i][j] = 0;
681 for (size_type k=0; k<cols(); k++)
682 C[i][j] += (*this)[i][k]*M[k][j];
683 }
684 }
685 return C;
686 }
687#endif
688
689 //===== sizes
690
692 size_type N () const
693 {
694 return rows();
695 }
696
698 size_type M () const
699 {
700 return cols();
701 }
702
705 {
706 return asImp().mat_rows();
707 }
708
711 {
712 return asImp().mat_cols();
713 }
714
715 //===== query
716
718 bool exists (size_type i, size_type j) const
719 {
720#ifdef DUNE_FMatrix_WITH_CHECKING
721 if (i<0 || i>=rows()) DUNE_THROW(FMatrixError,"row index out of range");
722 if (j<0 || j>=cols()) DUNE_THROW(FMatrixError,"column index out of range");
723#else
726#endif
727 return true;
728 }
729
730 private:
731
732#ifndef DOXYGEN
733 struct ElimPivot
734 {
735 ElimPivot(std::vector<size_type> & pivot);
736
737 void swap(int i, int j);
738
739 template<typename T>
740 void operator()(const T&, int, int)
741 {}
742
743 std::vector<size_type> & pivot_;
744 };
745
746 template<typename V>
747 struct Elim
748 {
749 Elim(V& rhs);
750
751 void swap(int i, int j);
752
753 void operator()(const typename V::field_type& factor, int k, int i);
754
755 V* rhs_;
756 };
757
758 struct ElimDet
759 {
760 ElimDet(field_type& sign) : sign_(sign)
761 { sign_ = 1; }
762
763 void swap(int, int)
764 { sign_ *= -1; }
765
766 void operator()(const field_type&, int, int)
767 {}
768
769 field_type& sign_;
770 };
771#endif // DOXYGEN
772
773 template<class Func>
774 void luDecomposition(DenseMatrix<MAT>& A, Func func) const;
775 };
776
777#ifndef DOXYGEN
778 template<typename MAT>
779 DenseMatrix<MAT>::ElimPivot::ElimPivot(std::vector<size_type> & pivot)
780 : pivot_(pivot)
781 {
782 typedef typename std::vector<size_type>::size_type size_type;
783 for(size_type i=0; i < pivot_.size(); ++i) pivot_[i]=i;
784 }
785
786 template<typename MAT>
787 void DenseMatrix<MAT>::ElimPivot::swap(int i, int j)
788 {
789 pivot_[i]=j;
790 }
791
792 template<typename MAT>
793 template<typename V>
794 DenseMatrix<MAT>::Elim<V>::Elim(V& rhs)
795 : rhs_(&rhs)
796 {}
797
798 template<typename MAT>
799 template<typename V>
800 void DenseMatrix<MAT>::Elim<V>::swap(int i, int j)
801 {
802 std::swap((*rhs_)[i], (*rhs_)[j]);
803 }
804
805 template<typename MAT>
806 template<typename V>
807 void DenseMatrix<MAT>::
808 Elim<V>::operator()(const typename V::field_type& factor, int k, int i)
809 {
810 (*rhs_)[k] -= factor*(*rhs_)[i];
811 }
812 template<typename MAT>
813 template<typename Func>
814 inline void DenseMatrix<MAT>::luDecomposition(DenseMatrix<MAT>& A, Func func) const
815 {
816 typedef typename FieldTraits<value_type>::real_type
817 real_type;
818 real_type norm = A.infinity_norm_real(); // for relative thresholds
821
822 // LU decomposition of A in A
823 for (size_type i=0; i<rows(); i++) // loop over all rows
824 {
825 typename FieldTraits<value_type>::real_type pivmax=fvmeta::absreal(A[i][i]);
826
827 // pivoting ?
828 if (pivmax<pivthres)
829 {
830 // compute maximum of column
831 size_type imax=i;
832 typename FieldTraits<value_type>::real_type abs(0.0);
833 for (size_type k=i+1; k<rows(); k++)
834 if ((abs=fvmeta::absreal(A[k][i]))>pivmax)
835 {
836 pivmax = abs; imax = k;
837 }
838 // swap rows
839 if (imax!=i) {
840 for (size_type j=0; j<rows(); j++)
841 std::swap(A[i][j],A[imax][j]);
842 func.swap(i, imax); // swap the pivot or rhs
843 }
844 }
845
846 // singular ?
847 if (pivmax<singthres)
848 DUNE_THROW(FMatrixError,"matrix is singular");
849
850 // eliminate
851 for (size_type k=i+1; k<rows(); k++)
852 {
853 field_type factor = A[k][i]/A[i][i];
854 A[k][i] = factor;
855 for (size_type j=i+1; j<rows(); j++)
856 A[k][j] -= factor*A[i][j];
857 func(factor, k, i);
858 }
859 }
860 }
861
862 template<typename MAT>
863 template <class V>
864 inline void DenseMatrix<MAT>::solve(V& x, const V& b) const
865 {
866 // never mind those ifs, because they get optimized away
867 if (rows()!=cols())
868 DUNE_THROW(FMatrixError, "Can't solve for a " << rows() << "x" << cols() << " matrix!");
869
870 if (rows()==1) {
871
872#ifdef DUNE_FMatrix_WITH_CHECKING
873 if (fvmeta::absreal((*this)[0][0])<FMatrixPrecision<>::absolute_limit())
874 DUNE_THROW(FMatrixError,"matrix is singular");
875#endif
876 x[0] = b[0]/(*this)[0][0];
877
878 }
879 else if (rows()==2) {
880
881 field_type detinv = (*this)[0][0]*(*this)[1][1]-(*this)[0][1]*(*this)[1][0];
882#ifdef DUNE_FMatrix_WITH_CHECKING
883 if (fvmeta::absreal(detinv)<FMatrixPrecision<>::absolute_limit())
884 DUNE_THROW(FMatrixError,"matrix is singular");
885#endif
886 detinv = 1.0/detinv;
887
888 x[0] = detinv*((*this)[1][1]*b[0]-(*this)[0][1]*b[1]);
889 x[1] = detinv*((*this)[0][0]*b[1]-(*this)[1][0]*b[0]);
890
891 }
892 else if (rows()==3) {
893
894 field_type d = determinant();
895#ifdef DUNE_FMatrix_WITH_CHECKING
896 if (fvmeta::absreal(d)<FMatrixPrecision<>::absolute_limit())
897 DUNE_THROW(FMatrixError,"matrix is singular");
898#endif
899
900 x[0] = (b[0]*(*this)[1][1]*(*this)[2][2] - b[0]*(*this)[2][1]*(*this)[1][2]
901 - b[1] *(*this)[0][1]*(*this)[2][2] + b[1]*(*this)[2][1]*(*this)[0][2]
902 + b[2] *(*this)[0][1]*(*this)[1][2] - b[2]*(*this)[1][1]*(*this)[0][2]) / d;
903
904 x[1] = ((*this)[0][0]*b[1]*(*this)[2][2] - (*this)[0][0]*b[2]*(*this)[1][2]
905 - (*this)[1][0] *b[0]*(*this)[2][2] + (*this)[1][0]*b[2]*(*this)[0][2]
906 + (*this)[2][0] *b[0]*(*this)[1][2] - (*this)[2][0]*b[1]*(*this)[0][2]) / d;
907
908 x[2] = ((*this)[0][0]*(*this)[1][1]*b[2] - (*this)[0][0]*(*this)[2][1]*b[1]
909 - (*this)[1][0] *(*this)[0][1]*b[2] + (*this)[1][0]*(*this)[2][1]*b[0]
910 + (*this)[2][0] *(*this)[0][1]*b[1] - (*this)[2][0]*(*this)[1][1]*b[0]) / d;
911
912 }
913 else {
914
915 V& rhs = x; // use x to store rhs
916 rhs = b; // copy data
917 Elim<V> elim(rhs);
918 MAT A(asImp());
919
920 luDecomposition(A, elim);
921
922 // backsolve
923 for(int i=rows()-1; i>=0; i--) {
924 for (size_type j=i+1; j<rows(); j++)
925 rhs[i] -= A[i][j]*x[j];
926 x[i] = rhs[i]/A[i][i];
927 }
928 }
929 }
930
931 template<typename MAT>
932 inline void DenseMatrix<MAT>::invert()
933 {
934 // never mind those ifs, because they get optimized away
935 if (rows()!=cols())
936 DUNE_THROW(FMatrixError, "Can't invert a " << rows() << "x" << cols() << " matrix!");
937
938 if (rows()==1) {
939
940#ifdef DUNE_FMatrix_WITH_CHECKING
941 if (fvmeta::absreal((*this)[0][0])<FMatrixPrecision<>::absolute_limit())
942 DUNE_THROW(FMatrixError,"matrix is singular");
943#endif
944 (*this)[0][0] = field_type( 1 ) / (*this)[0][0];
945
946 }
947 else if (rows()==2) {
948
949 field_type detinv = (*this)[0][0]*(*this)[1][1]-(*this)[0][1]*(*this)[1][0];
950#ifdef DUNE_FMatrix_WITH_CHECKING
951 if (fvmeta::absreal(detinv)<FMatrixPrecision<>::absolute_limit())
952 DUNE_THROW(FMatrixError,"matrix is singular");
953#endif
954 detinv = field_type( 1 ) / detinv;
955
956 field_type temp=(*this)[0][0];
957 (*this)[0][0] = (*this)[1][1]*detinv;
958 (*this)[0][1] = -(*this)[0][1]*detinv;
959 (*this)[1][0] = -(*this)[1][0]*detinv;
960 (*this)[1][1] = temp*detinv;
961
962 }
963 else {
964
965 MAT A(asImp());
966 std::vector<size_type> pivot(rows());
967 luDecomposition(A, ElimPivot(pivot));
968 DenseMatrix<MAT>& L=A;
969 DenseMatrix<MAT>& U=A;
970
971 // initialize inverse
972 *this=field_type();
973
974 for(size_type i=0; i<rows(); ++i)
975 (*this)[i][i]=1;
976
977 // L Y = I; multiple right hand sides
978 for (size_type i=0; i<rows(); i++)
979 for (size_type j=0; j<i; j++)
980 for (size_type k=0; k<rows(); k++)
981 (*this)[i][k] -= L[i][j]*(*this)[j][k];
982
983 // U A^{-1} = Y
984 for (size_type i=rows(); i>0;) {
985 --i;
986 for (size_type k=0; k<rows(); k++) {
987 for (size_type j=i+1; j<rows(); j++)
988 (*this)[i][k] -= U[i][j]*(*this)[j][k];
989 (*this)[i][k] /= U[i][i];
990 }
991 }
992
993 for(size_type i=rows(); i>0; ) {
994 --i;
995 if(i!=pivot[i])
996 for(size_type j=0; j<rows(); ++j)
997 std::swap((*this)[j][pivot[i]], (*this)[j][i]);
998 }
999 }
1000 }
1001
1002 // implementation of the determinant
1003 template<typename MAT>
1004 inline typename DenseMatrix<MAT>::field_type
1006 {
1007 // never mind those ifs, because they get optimized away
1008 if (rows()!=cols())
1009 DUNE_THROW(FMatrixError, "There is no determinant for a " << rows() << "x" << cols() << " matrix!");
1010
1011 if (rows()==1)
1012 return (*this)[0][0];
1013
1014 if (rows()==2)
1015 return (*this)[0][0]*(*this)[1][1] - (*this)[0][1]*(*this)[1][0];
1016
1017 if (rows()==3) {
1018 // code generated by maple
1019 field_type t4 = (*this)[0][0] * (*this)[1][1];
1020 field_type t6 = (*this)[0][0] * (*this)[1][2];
1021 field_type t8 = (*this)[0][1] * (*this)[1][0];
1022 field_type t10 = (*this)[0][2] * (*this)[1][0];
1023 field_type t12 = (*this)[0][1] * (*this)[2][0];
1024 field_type t14 = (*this)[0][2] * (*this)[2][0];
1025
1026 return (t4*(*this)[2][2]-t6*(*this)[2][1]-t8*(*this)[2][2]+
1027 t10*(*this)[2][1]+t12*(*this)[1][2]-t14*(*this)[1][1]);
1028
1029 }
1030
1031 MAT A(asImp());
1032 field_type det;
1033 try
1034 {
1035 luDecomposition(A, ElimDet(det));
1036 }
1037 catch (FMatrixError&)
1038 {
1039 return 0;
1040 }
1041 for (size_type i = 0; i < rows(); ++i)
1042 det *= A[i][i];
1043 return det;
1044 }
1045
1046#endif // DOXYGEN
1047
1048 namespace DenseMatrixHelp {
1049#if 0
1051 template <typename K>
1052 static inline K invertMatrix (const FieldMatrix<K,1,1> &matrix, FieldMatrix<K,1,1> &inverse)
1053 {
1054 inverse[0][0] = 1.0/matrix[0][0];
1055 return matrix[0][0];
1056 }
1057
1059 template <typename K>
1060 static inline K invertMatrix_retTransposed (const FieldMatrix<K,1,1> &matrix, FieldMatrix<K,1,1> &inverse)
1061 {
1062 return invertMatrix(matrix,inverse);
1063 }
1064
1065
1067 template <typename K>
1068 static inline K invertMatrix (const FieldMatrix<K,2,2> &matrix, FieldMatrix<K,2,2> &inverse)
1069 {
1070 // code generated by maple
1071 field_type det = (matrix[0][0]*matrix[1][1] - matrix[0][1]*matrix[1][0]);
1072 field_type det_1 = 1.0/det;
1073 inverse[0][0] = matrix[1][1] * det_1;
1074 inverse[0][1] = - matrix[0][1] * det_1;
1075 inverse[1][0] = - matrix[1][0] * det_1;
1076 inverse[1][1] = matrix[0][0] * det_1;
1077 return det;
1078 }
1079
1082 template <typename K>
1083 static inline K invertMatrix_retTransposed (const FieldMatrix<K,2,2> &matrix, FieldMatrix<K,2,2> &inverse)
1084 {
1085 // code generated by maple
1086 field_type det = (matrix[0][0]*matrix[1][1] - matrix[0][1]*matrix[1][0]);
1087 field_type det_1 = 1.0/det;
1088 inverse[0][0] = matrix[1][1] * det_1;
1089 inverse[1][0] = - matrix[0][1] * det_1;
1090 inverse[0][1] = - matrix[1][0] * det_1;
1091 inverse[1][1] = matrix[0][0] * det_1;
1092 return det;
1093 }
1094
1096 template <typename K>
1097 static inline K invertMatrix (const FieldMatrix<K,3,3> &matrix, FieldMatrix<K,3,3> &inverse)
1098 {
1099 // code generated by maple
1100 field_type t4 = matrix[0][0] * matrix[1][1];
1101 field_type t6 = matrix[0][0] * matrix[1][2];
1102 field_type t8 = matrix[0][1] * matrix[1][0];
1103 field_type t10 = matrix[0][2] * matrix[1][0];
1104 field_type t12 = matrix[0][1] * matrix[2][0];
1105 field_type t14 = matrix[0][2] * matrix[2][0];
1106
1107 field_type det = (t4*matrix[2][2]-t6*matrix[2][1]-t8*matrix[2][2]+
1108 t10*matrix[2][1]+t12*matrix[1][2]-t14*matrix[1][1]);
1109 field_type t17 = 1.0/det;
1110
1111 inverse[0][0] = (matrix[1][1] * matrix[2][2] - matrix[1][2] * matrix[2][1])*t17;
1112 inverse[0][1] = -(matrix[0][1] * matrix[2][2] - matrix[0][2] * matrix[2][1])*t17;
1113 inverse[0][2] = (matrix[0][1] * matrix[1][2] - matrix[0][2] * matrix[1][1])*t17;
1114 inverse[1][0] = -(matrix[1][0] * matrix[2][2] - matrix[1][2] * matrix[2][0])*t17;
1115 inverse[1][1] = (matrix[0][0] * matrix[2][2] - t14) * t17;
1116 inverse[1][2] = -(t6-t10) * t17;
1117 inverse[2][0] = (matrix[1][0] * matrix[2][1] - matrix[1][1] * matrix[2][0]) * t17;
1118 inverse[2][1] = -(matrix[0][0] * matrix[2][1] - t12) * t17;
1119 inverse[2][2] = (t4-t8) * t17;
1120
1121 return det;
1122 }
1123
1125 template <typename K>
1126 static inline K invertMatrix_retTransposed (const FieldMatrix<K,3,3> &matrix, FieldMatrix<K,3,3> &inverse)
1127 {
1128 // code generated by maple
1129 field_type t4 = matrix[0][0] * matrix[1][1];
1130 field_type t6 = matrix[0][0] * matrix[1][2];
1131 field_type t8 = matrix[0][1] * matrix[1][0];
1132 field_type t10 = matrix[0][2] * matrix[1][0];
1133 field_type t12 = matrix[0][1] * matrix[2][0];
1134 field_type t14 = matrix[0][2] * matrix[2][0];
1135
1136 field_type det = (t4*matrix[2][2]-t6*matrix[2][1]-t8*matrix[2][2]+
1137 t10*matrix[2][1]+t12*matrix[1][2]-t14*matrix[1][1]);
1138 field_type t17 = 1.0/det;
1139
1140 inverse[0][0] = (matrix[1][1] * matrix[2][2] - matrix[1][2] * matrix[2][1])*t17;
1141 inverse[1][0] = -(matrix[0][1] * matrix[2][2] - matrix[0][2] * matrix[2][1])*t17;
1142 inverse[2][0] = (matrix[0][1] * matrix[1][2] - matrix[0][2] * matrix[1][1])*t17;
1143 inverse[0][1] = -(matrix[1][0] * matrix[2][2] - matrix[1][2] * matrix[2][0])*t17;
1144 inverse[1][1] = (matrix[0][0] * matrix[2][2] - t14) * t17;
1145 inverse[2][1] = -(t6-t10) * t17;
1146 inverse[0][2] = (matrix[1][0] * matrix[2][1] - matrix[1][1] * matrix[2][0]) * t17;
1147 inverse[1][2] = -(matrix[0][0] * matrix[2][1] - t12) * t17;
1148 inverse[2][2] = (t4-t8) * t17;
1149
1150 return det;
1151 }
1152
1154 template< class K, int m, int n, int p >
1155 static inline void multMatrix ( const FieldMatrix< K, m, n > &A,
1156 const FieldMatrix< K, n, p > &B,
1157 FieldMatrix< K, m, p > &ret )
1158 {
1159 typedef typename FieldMatrix< K, m, p > :: size_type size_type;
1160
1161 for( size_type i = 0; i < m; ++i )
1162 {
1163 for( size_type j = 0; j < p; ++j )
1164 {
1165 ret[ i ][ j ] = K( 0 );
1166 for( size_type k = 0; k < n; ++k )
1167 ret[ i ][ j ] += A[ i ][ k ] * B[ k ][ j ];
1168 }
1169 }
1170 }
1171
1173 template <typename K, int rows, int cols>
1174 static inline void multTransposedMatrix(const FieldMatrix<K,rows,cols> &matrix, FieldMatrix<K,cols,cols>& ret)
1175 {
1176 typedef typename FieldMatrix<K,rows,cols>::size_type size_type;
1177
1178 for(size_type i=0; i<cols(); i++)
1179 for(size_type j=0; j<cols(); j++)
1180 {
1181 ret[i][j]=0.0;
1182 for(size_type k=0; k<rows(); k++)
1183 ret[i][j]+=matrix[k][i]*matrix[k][j];
1184 }
1185 }
1186#endif
1187
1189 template <typename MAT, typename V1, typename V2>
1190 static inline void multAssign(const DenseMatrix<MAT> &matrix, const DenseVector<V1> & x, DenseVector<V2> & ret)
1191 {
1192 assert(x.size() == matrix.cols());
1193 assert(ret.size() == matrix.rows());
1194 typedef typename DenseMatrix<MAT>::size_type size_type;
1195
1196 for(size_type i=0; i<matrix.rows(); ++i)
1197 {
1198 ret[i] = 0.0;
1199 for(size_type j=0; j<matrix.cols(); ++j)
1200 {
1201 ret[i] += matrix[i][j]*x[j];
1202 }
1203 }
1204 }
1205
1206#if 0
1208 template <typename K, int rows, int cols>
1209 static inline void multAssignTransposed( const FieldMatrix<K,rows,cols> &matrix, const FieldVector<K,rows> & x, FieldVector<K,cols> & ret)
1210 {
1211 typedef typename FieldMatrix<K,rows,cols>::size_type size_type;
1212
1213 for(size_type i=0; i<cols(); ++i)
1214 {
1215 ret[i] = 0.0;
1216 for(size_type j=0; j<rows(); ++j)
1217 ret[i] += matrix[j][i]*x[j];
1218 }
1219 }
1220
1222 template <typename K, int rows, int cols>
1223 static inline FieldVector<K,rows> mult(const FieldMatrix<K,rows,cols> &matrix, const FieldVector<K,cols> & x)
1224 {
1225 FieldVector<K,rows> ret;
1226 multAssign(matrix,x,ret);
1227 return ret;
1228 }
1229
1231 template <typename K, int rows, int cols>
1232 static inline FieldVector<K,cols> multTransposed(const FieldMatrix<K,rows,cols> &matrix, const FieldVector<K,rows> & x)
1233 {
1234 FieldVector<K,cols> ret;
1235 multAssignTransposed( matrix, x, ret );
1236 return ret;
1237 }
1238#endif
1239
1240 } // end namespace DenseMatrixHelp
1241
1243 template<typename MAT>
1244 std::ostream& operator<< (std::ostream& s, const DenseMatrix<MAT>& a)
1245 {
1246 for (typename DenseMatrix<MAT>::size_type i=0; i<a.rows(); i++)
1247 s << a[i] << std::endl;
1248 return s;
1249 }
1250
1253} // end namespace Dune
1254
1255#endif
Checks wether a type is convertible to another.
Definition: typetraits.hh:213
@ exists
True if the conversion exists.
Definition: typetraits.hh:223
Generic iterator class for dense vector and matrix implementations.
Definition: densevector.hh:122
A dense n x m matrix.
Definition: densematrix.hh:192
DenseMatrix & operator/=(const field_type &k)
vector space division by scalar
Definition: densematrix.hh:360
ConstIterator const_iterator
typedef for stl compliant access
Definition: densematrix.hh:290
void usmtv(const field_type &alpha, const X &x, Y &y) const
y += alpha A^T x
Definition: densematrix.hh:532
bool operator!=(const DenseMatrix< Other > &y) const
Binary matrix incomparison.
Definition: densematrix.hh:387
FieldTraits< value_type >::real_type infinity_norm() const
infinity norm (row sum norm, how to generalize for blocks?)
Definition: densematrix.hh:577
void mv(const X &x, Y &y) const
y = A x
Definition: densematrix.hh:397
Traits::value_type field_type
export the type representing the field
Definition: densematrix.hh:209
void solve(V &x, const V &b) const
Solve system A x = b.
size_type cols() const
number of columns
Definition: densematrix.hh:710
ConstIterator beforeEnd() const
Definition: densematrix.hh:310
bool operator==(const DenseMatrix< Other > &y) const
Binary matrix comparison.
Definition: densematrix.hh:378
void mmtv(const X &x, Y &y) const
y -= A^T x
Definition: densematrix.hh:491
DenseMatrix & operator+=(const DenseMatrix< Other > &y)
vector space addition
Definition: densematrix.hh:335
ConstIterator beforeBegin() const
Definition: densematrix.hh:317
DenseIterator< DenseMatrix, row_type > Iterator
Iterator class for sequential access.
Definition: densematrix.hh:253
Traits::value_type block_type
export the type representing the components
Definition: densematrix.hh:212
void mtv(const X &x, Y &y) const
y = A^T x
Definition: densematrix.hh:413
DenseMatrix & operator*=(const field_type &k)
vector space multiplication with scalar
Definition: densematrix.hh:352
void usmhv(const field_type &alpha, const X &x, Y &y) const
y += alpha A^H x
Definition: densematrix.hh:546
DenseMatrix & operator-=(const DenseMatrix< Other > &y)
vector space subtraction
Definition: densematrix.hh:344
FieldTraits< value_type >::real_type infinity_norm_real() const
simplified infinity norm (uses Manhattan norm for complex values)
Definition: densematrix.hh:591
row_type::Iterator ColIterator
rename the iterators for easier access
Definition: densematrix.hh:259
size_type size() const
size method (number of rows)
Definition: densematrix.hh:246
size_type M() const
number of columns
Definition: densematrix.hh:698
MAT & rightmultiply(const DenseMatrix< M2 > &M)
Multiplies M from the right to this matrix.
Definition: densematrix.hh:641
void usmv(const field_type &alpha, const X &x, Y &y) const
y += alpha A x
Definition: densematrix.hh:519
Iterator end()
end iterator
Definition: densematrix.hh:268
Iterator beforeBegin()
Definition: densematrix.hh:282
Iterator beforeEnd()
Definition: densematrix.hh:275
Traits::value_type value_type
export the type representing the field
Definition: densematrix.hh:206
@ blocklevel
The number of block levels we contain. This is 1.
Definition: densematrix.hh:229
size_type rows() const
number of rows
Definition: densematrix.hh:704
void mmv(const X &x, Y &y) const
y -= A x
Definition: densematrix.hh:478
void invert()
Compute inverse.
MAT & leftmultiply(const DenseMatrix< M2 > &M)
Multiplies M from the left to this matrix.
Definition: densematrix.hh:624
void mmhv(const X &x, Y &y) const
y -= A^H x
Definition: densematrix.hh:505
Traits::derived_type derived_type
type of derived matrix class
Definition: densematrix.hh:203
row_reference operator[](size_type i)
random access
Definition: densematrix.hh:235
row_type::ConstIterator ConstColIterator
rename the iterators for easier access
Definition: densematrix.hh:294
bool exists(size_type i, size_type j) const
return true when (i,j) is in pattern
Definition: densematrix.hh:718
DenseMatrix & axpy(const field_type &k, const DenseMatrix< Other > &y)
vector space axpy operation (*this += k y)
Definition: densematrix.hh:369
Iterator RowIterator
rename the iterators for easier access
Definition: densematrix.hh:257
FieldTraits< value_type >::real_type frobenius_norm() const
frobenius norm: sqrt(sum over squared values of entries)
Definition: densematrix.hh:561
void umv(const X &x, Y &y) const
y += A x
Definition: densematrix.hh:435
size_type N() const
number of rows
Definition: densematrix.hh:692
Traits::row_type row_type
The type used to represent a row (must fulfill the Dune::DenseVector interface)
Definition: densematrix.hh:218
Traits::size_type size_type
The type used for the index access and size operation.
Definition: densematrix.hh:215
Traits::const_row_reference const_row_reference
The type used to represent a reference to a constant row (usually const row_type &)
Definition: densematrix.hh:224
FieldTraits< value_type >::real_type frobenius_norm2() const
square of frobenius norm, need for block recursion
Definition: densematrix.hh:569
Traits::row_reference row_reference
The type used to represent a reference to a row (usually row_type &)
Definition: densematrix.hh:221
Iterator iterator
typedef for stl compliant access
Definition: densematrix.hh:255
ConstIterator ConstRowIterator
rename the iterators for easier access
Definition: densematrix.hh:292
field_type determinant() const
calculates the determinant of this matrix
void umtv(const X &x, Y &y) const
y += A^T x
Definition: densematrix.hh:450
ConstIterator begin() const
begin iterator
Definition: densematrix.hh:297
DenseIterator< const DenseMatrix, const row_type > ConstIterator
Iterator class for sequential access.
Definition: densematrix.hh:288
Iterator begin()
begin iterator
Definition: densematrix.hh:262
void umhv(const X &x, Y &y) const
y += A^H x
Definition: densematrix.hh:464
ConstIterator end() const
end iterator
Definition: densematrix.hh:303
Interface for a class of dense vectors over a given field.
Definition: densevector.hh:223
size_type size() const
size method
Definition: densevector.hh:285
size_type N() const
number of blocks in the vector (are of size 1 here)
Definition: densevector.hh:569
Base class for Dune-Exceptions.
Definition: exceptions.hh:92
Error thrown if operations of a FieldMatrix fail.
Definition: densematrix.hh:178
static ctype pivoting_limit()
return threshold to do pivoting
Definition: precision.hh:26
static ctype absolute_limit()
return threshold to declare matrix singular
Definition: precision.hh:50
static ctype singular_limit()
return threshold to declare matrix singular
Definition: precision.hh:38
A dense n x m matrix.
Definition: fmatrix.hh:67
vector space out of a tensor product of fields.
Definition: fvector.hh:92
FieldVector()
Constructor making default-initialized vector.
Definition: fvector.hh:106
FieldTraits< field_type >::real_type infinity_norm_real() const
simplified infinity norm (uses Manhattan norm for complex values)
Definition: matrix.hh:518
A free function to provide the demangled class name of a given object or type as a string.
static void multAssign(const DenseMatrix< MAT > &matrix, const DenseVector< V1 > &x, DenseVector< V2 > &ret)
calculates ret = matrix * x
Definition: densematrix.hh:1190
A few common exception classes.
Implements a vector constructed from a given type representing a field and a compile-time given size.
#define dune_static_assert(COND, MSG)
Helper template so that compilation fails if condition is not true.
Definition: static_assert.hh:79
std::ostream & operator<<(std::ostream &s, const array< T, N > &e)
Output operator for array.
Definition: array.hh:159
#define DUNE_THROW(E, m)
Definition: exceptions.hh:244
Some useful basic math stuff.
Dune namespace.
Definition: alignment.hh:14
int sign(const T &val)
Return the sign of the value.
Definition: math.hh:119
K conjugateComplex(const K &x)
compute conjugate complex of x
Definition: math.hh:103
Various precision settings for calculations with FieldMatrix and FieldVector.
Fallback implementation of the C++0x static_assert feature.
you have to specialize this structure for any type that should be assignable to a DenseMatrix
Definition: densematrix.hh:73
Removes a const qualifier while preserving others.
Definition: typetraits.hh:176
Definition of the DUNE_UNUSED macro for the case that config.h is not available.
#define DUNE_UNUSED_PARAMETER(parm)
A macro to mark intentional unused function parameters with.
Definition: unused.hh:18
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 12, 23:30, 2024)