dune-common  2.2.1
densematrix.hh
Go to the documentation of this file.
1 // -*- tab-width: 8; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2 // vi: set et ts=8 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 
12 #include <dune/common/misc.hh>
14 #include <dune/common/fvector.hh>
15 #include <dune/common/precision.hh>
17 #include <dune/common/classname.hh>
18 
19 
20 namespace Dune
21 {
22 
23  template<typename M> class DenseMatrix;
24 
25  template<typename M>
26  struct FieldTraits< DenseMatrix<M> >
27  {
30  };
31 
32  /*
33  work around a problem of FieldMatrix/FieldVector,
34  there is no unique way to obtain the size of a class
35  */
36  template<class K, int N, int M> class FieldMatrix;
37  template<class K, int N> class FieldVector;
38  namespace {
39  template<class V>
40  struct VectorSize
41  {
42  static typename V::size_type size(const V & v) { return v.size(); }
43  };
44 
45  template<class K, int N>
46  struct VectorSize< const FieldVector<K,N> >
47  {
48  typedef FieldVector<K,N> V;
49  static typename V::size_type size(const V & v) { return N; }
50  };
51  }
52 
68  template<typename M, typename T>
70  {
71  DUNE_THROW(NotImplemented, "You need to specialise the method istl_assign_to_fmatrix(DenseMatrix<M>& f, const T& t) "
72  << "(with M being " << className<M>() << ") "
73  << "for T == " << className<T>() << "!");
74  }
75 
76  namespace
77  {
78  template<bool b>
79  struct DenseMatrixAssigner
80  {
81  template<typename M, typename T>
82  static void assign(DenseMatrix<M>& fm, const T& t)
83  {
85  }
86 
87  };
88 
89 
90  template<>
91  struct DenseMatrixAssigner<true>
92  {
93  template<typename M, typename T>
94  static void assign(DenseMatrix<M>& fm, const T& t)
95  {
96  fm = static_cast<const typename DenseMatVecTraits<M>::value_type>(t);
97  }
98  };
99  }
100 
102  class FMatrixError : public Exception {};
103 
114  template<typename MAT>
116  {
118 
119  // Curiously recurring template pattern
120  MAT & asImp() { return static_cast<MAT&>(*this); }
121  const MAT & asImp() const { return static_cast<const MAT&>(*this); }
122 
123  public:
124  //===== type definitions and constants
125 
128 
130  typedef typename Traits::value_type value_type;
131 
133  typedef typename Traits::value_type field_type;
134 
136  typedef typename Traits::value_type block_type;
137 
139  typedef typename Traits::size_type size_type;
140 
142  typedef typename Traits::row_type row_type;
143 
146 
149 
151  enum {
154  };
155 
156  //===== access to components
157 
160  {
161  return asImp().mat_access(i);
162  }
163 
165  {
166  return asImp().mat_access(i);
167  }
168 
170  size_type size() const
171  {
172  return rows();
173  }
174 
175  //===== iterator interface to rows of the matrix
183  typedef typename row_type::Iterator ColIterator;
184 
187  {
188  return Iterator(*this,0);
189  }
190 
193  {
194  return Iterator(*this,rows());
195  }
196 
200  {
201  return Iterator(*this,rows()-1);
202  }
203 
207  {
208  return Iterator(*this,-1);
209  }
210 
219 
222  {
223  return ConstIterator(*this,0);
224  }
225 
228  {
229  return ConstIterator(*this,rows());
230  }
231 
235  {
236  return ConstIterator(*this,rows()-1);
237  }
238 
242  {
243  return ConstIterator(*this,-1);
244  }
245 
246  //===== assignment from scalar
248  {
249  for (size_type i=0; i<rows(); i++)
250  (*this)[i] = f;
251  return *this;
252  }
253 
254  template<typename T>
256  {
257  DenseMatrixAssigner<Conversion<T,field_type>::exists>::assign(*this, t);
258  return *this;
259  }
260  //===== vector space arithmetic
261 
263  template <class Other>
265  {
266  for (size_type i=0; i<rows(); i++)
267  (*this)[i] += y[i];
268  return *this;
269  }
270 
272  template <class Other>
274  {
275  for (size_type i=0; i<rows(); i++)
276  (*this)[i] -= y[i];
277  return *this;
278  }
279 
282  {
283  for (size_type i=0; i<rows(); i++)
284  (*this)[i] *= k;
285  return *this;
286  }
287 
290  {
291  for (size_type i=0; i<rows(); i++)
292  (*this)[i] /= k;
293  return *this;
294  }
295 
297  template <class Other>
299  {
300  for( size_type i = 0; i < rows(); ++i )
301  (*this)[ i ].axpy( k, y[ i ] );
302  return *this;
303  }
304 
306  template <class Other>
307  bool operator== (const DenseMatrix<Other>& y) const
308  {
309  for (size_type i=0; i<rows(); i++)
310  if ((*this)[i]!=y[i])
311  return false;
312  return true;
313  }
315  template <class Other>
316  bool operator!= (const DenseMatrix<Other>& y) const
317  {
318  return !operator==(y);
319  }
320 
321 
322  //===== linear maps
323 
325  template<class X, class Y>
326  void mv (const X& x, Y& y) const
327  {
328 #ifdef DUNE_FMatrix_WITH_CHECKING
329  if (x.N()!=M()) DUNE_THROW(FMatrixError,"Index out of range");
330  if (y.N()!=N()) DUNE_THROW(FMatrixError,"Index out of range");
331 #endif
332  for (size_type i=0; i<rows(); ++i)
333  {
334  y[i] = 0;
335  for (size_type j=0; j<cols(); j++)
336  y[i] += (*this)[i][j] * x[j];
337  }
338  }
339 
341  template< class X, class Y >
342  void mtv ( const X &x, Y &y ) const
343  {
344 #ifdef DUNE_FMatrix_WITH_CHECKING
345  //assert( &x != &y );
346  //This assert did not work for me. Compile error:
347  // comparison between distinct pointer types ‘const
348  // Dune::FieldVector<double, 3>*’ and ‘Dune::FieldVector<double, 2>*’ lacks a cast
349  if( x.N() != N() )
350  DUNE_THROW( FMatrixError, "Index out of range." );
351  if( y.N() != M() )
352  DUNE_THROW( FMatrixError, "Index out of range." );
353 #endif
354  for( size_type i = 0; i < cols(); ++i )
355  {
356  y[ i ] = 0;
357  for( size_type j = 0; j < rows(); ++j )
358  y[ i ] += (*this)[ j ][ i ] * x[ j ];
359  }
360  }
361 
363  template<class X, class Y>
364  void umv (const X& x, Y& y) const
365  {
366 #ifdef DUNE_FMatrix_WITH_CHECKING
367  if (x.N()!=M())
368  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);
369  if (y.N()!=N())
370  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);
371 #endif
372  for (size_type i=0; i<rows(); i++)
373  for (size_type j=0; j<cols(); j++)
374  y[i] += (*this)[i][j] * x[j];
375  }
376 
378  template<class X, class Y>
379  void umtv (const X& x, Y& y) const
380  {
381 #ifdef DUNE_FMatrix_WITH_CHECKING
382  if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
383  if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
384 #endif
385 
386  for (size_type i=0; i<rows(); i++)
387  for (size_type j=0; j<cols(); j++)
388  y[j] += (*this)[i][j]*x[i];
389  }
390 
392  template<class X, class Y>
393  void umhv (const X& x, Y& y) const
394  {
395 #ifdef DUNE_FMatrix_WITH_CHECKING
396  if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
397  if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
398 #endif
399 
400  for (size_type i=0; i<rows(); i++)
401  for (size_type j=0; j<cols(); j++)
402  y[j] += conjugateComplex((*this)[i][j])*x[i];
403  }
404 
406  template<class X, class Y>
407  void mmv (const X& x, Y& y) const
408  {
409 #ifdef DUNE_FMatrix_WITH_CHECKING
410  if (x.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
411  if (y.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
412 #endif
413  for (size_type i=0; i<rows(); i++)
414  for (size_type j=0; j<cols(); j++)
415  y[i] -= (*this)[i][j] * x[j];
416  }
417 
419  template<class X, class Y>
420  void mmtv (const X& x, Y& y) const
421  {
422 #ifdef DUNE_FMatrix_WITH_CHECKING
423  if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
424  if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
425 #endif
426 
427  for (size_type i=0; i<rows(); i++)
428  for (size_type j=0; j<cols(); j++)
429  y[j] -= (*this)[i][j]*x[i];
430  }
431 
433  template<class X, class Y>
434  void mmhv (const X& x, Y& y) const
435  {
436 #ifdef DUNE_FMatrix_WITH_CHECKING
437  if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
438  if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
439 #endif
440 
441  for (size_type i=0; i<rows(); i++)
442  for (size_type j=0; j<cols(); j++)
443  y[j] -= conjugateComplex((*this)[i][j])*x[i];
444  }
445 
447  template<class X, class Y>
448  void usmv (const field_type& alpha, const X& x, Y& y) const
449  {
450 #ifdef DUNE_FMatrix_WITH_CHECKING
451  if (x.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
452  if (y.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
453 #endif
454  for (size_type i=0; i<rows(); i++)
455  for (size_type j=0; j<cols(); j++)
456  y[i] += alpha * (*this)[i][j] * x[j];
457  }
458 
460  template<class X, class Y>
461  void usmtv (const field_type& alpha, const X& x, Y& y) const
462  {
463 #ifdef DUNE_FMatrix_WITH_CHECKING
464  if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
465  if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
466 #endif
467 
468  for (size_type i=0; i<rows(); i++)
469  for (size_type j=0; j<cols(); j++)
470  y[j] += alpha*(*this)[i][j]*x[i];
471  }
472 
474  template<class X, class Y>
475  void usmhv (const field_type& alpha, const X& x, Y& y) const
476  {
477 #ifdef DUNE_FMatrix_WITH_CHECKING
478  if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
479  if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
480 #endif
481 
482  for (size_type i=0; i<rows(); i++)
483  for (size_type j=0; j<cols(); j++)
484  y[j] += alpha*conjugateComplex((*this)[i][j])*x[i];
485  }
486 
487  //===== norms
488 
491  {
492  typename FieldTraits<value_type>::real_type sum=(0.0);
493  for (size_type i=0; i<rows(); ++i) sum += (*this)[i].two_norm2();
494  return fvmeta::sqrt(sum);
495  }
496 
499  {
500  typename FieldTraits<value_type>::real_type sum=(0.0);
501  for (size_type i=0; i<rows(); ++i) sum += (*this)[i].two_norm2();
502  return sum;
503  }
504 
507  {
508  if (size() == 0)
509  return 0.0;
510 
511  ConstIterator it = begin();
512  typename remove_const< typename FieldTraits<value_type>::real_type >::type max = it->one_norm();
513  for (it = it + 1; it != end(); ++it)
514  max = std::max(max, it->one_norm());
515 
516  return max;
517  }
518 
521  {
522  if (size() == 0)
523  return 0.0;
524 
525  ConstIterator it = begin();
526  typename remove_const< typename FieldTraits<value_type>::real_type >::type max = it->one_norm_real();
527  for (it = it + 1; it != end(); ++it)
528  max = std::max(max, it->one_norm_real());
529 
530  return max;
531  }
532 
533  //===== solve
534 
539  template <class V>
540  void solve (V& x, const V& b) const;
541 
546  void invert();
547 
549  field_type determinant () const;
550 
552  template<typename M2>
554  {
555  assert(M.rows() == M.cols() && M.rows() == rows());
556  MAT C(asImp());
557 
558  for (size_type i=0; i<rows(); i++)
559  for (size_type j=0; j<cols(); j++) {
560  (*this)[i][j] = 0;
561  for (size_type k=0; k<rows(); k++)
562  (*this)[i][j] += M[i][k]*C[k][j];
563  }
564 
565  return asImp();
566  }
567 
569  template<typename M2>
571  {
572  assert(M.rows() == M.cols() && M.cols() == cols());
573  MAT C(asImp());
574 
575  for (size_type i=0; i<rows(); i++)
576  for (size_type j=0; j<cols(); j++) {
577  (*this)[i][j] = 0;
578  for (size_type k=0; k<cols(); k++)
579  (*this)[i][j] += C[i][k]*M[k][j];
580  }
581  return asImp();
582  }
583 
584 #if 0
585 
586  template<int l>
587  DenseMatrix<K,l,cols> leftmultiplyany (const FieldMatrix<K,l,rows>& M) const
588  {
590 
591  for (size_type i=0; i<l; i++) {
592  for (size_type j=0; j<cols(); j++) {
593  C[i][j] = 0;
594  for (size_type k=0; k<rows(); k++)
595  C[i][j] += M[i][k]*(*this)[k][j];
596  }
597  }
598  return C;
599  }
600 
602  template<int l>
603  FieldMatrix<K,rows,l> rightmultiplyany (const FieldMatrix<K,cols,l>& M) const
604  {
605  FieldMatrix<K,rows,l> C;
606 
607  for (size_type i=0; i<rows(); i++) {
608  for (size_type j=0; j<l; j++) {
609  C[i][j] = 0;
610  for (size_type k=0; k<cols(); k++)
611  C[i][j] += (*this)[i][k]*M[k][j];
612  }
613  }
614  return C;
615  }
616 #endif
617 
618  //===== sizes
619 
621  size_type N () const
622  {
623  return rows();
624  }
625 
627  size_type M () const
628  {
629  return cols();
630  }
631 
633  size_type rows() const
634  {
635  return asImp().mat_rows();
636  }
637 
639  size_type cols() const
640  {
641  return asImp().mat_cols();
642  }
643 
644  //===== query
645 
647  bool exists (size_type i, size_type j) const
648  {
649 #ifdef DUNE_FMatrix_WITH_CHECKING
650  if (i<0 || i>=rows()) DUNE_THROW(FMatrixError,"row index out of range");
651  if (j<0 || j>=cols()) DUNE_THROW(FMatrixError,"column index out of range");
652 #endif
653  return true;
654  }
655 
656  private:
657 
658 #ifndef DOXYGEN
659  struct ElimPivot
660  {
661  ElimPivot(std::vector<size_type> & pivot);
662 
663  void swap(int i, int j);
664 
665  template<typename T>
666  void operator()(const T&, int k, int i)
667  {}
668 
669  std::vector<size_type> & pivot_;
670  };
671 
672  template<typename V>
673  struct Elim
674  {
675  Elim(V& rhs);
676 
677  void swap(int i, int j);
678 
679  void operator()(const typename V::field_type& factor, int k, int i);
680 
681  V* rhs_;
682  };
683 
684  struct ElimDet
685  {
686  ElimDet(field_type& sign) : sign_(sign)
687  { sign_ = 1; }
688 
689  void swap(int i, int j)
690  { sign_ *= -1; }
691 
692  void operator()(const field_type&, int k, int i)
693  {}
694 
695  field_type& sign_;
696  };
697 #endif // DOXYGEN
698 
699  template<class Func>
700  void luDecomposition(DenseMatrix<MAT>& A, Func func) const;
701  };
702 
703 #ifndef DOXYGEN
704  template<typename MAT>
705  DenseMatrix<MAT>::ElimPivot::ElimPivot(std::vector<size_type> & pivot)
706  : pivot_(pivot)
707  {
708  typedef typename std::vector<size_type>::size_type size_type;
709  for(size_type i=0; i < pivot_.size(); ++i) pivot_[i]=i;
710  }
711 
712  template<typename MAT>
713  void DenseMatrix<MAT>::ElimPivot::swap(int i, int j)
714  {
715  pivot_[i]=j;
716  }
717 
718  template<typename MAT>
719  template<typename V>
720  DenseMatrix<MAT>::Elim<V>::Elim(V& rhs)
721  : rhs_(&rhs)
722  {}
723 
724  template<typename MAT>
725  template<typename V>
726  void DenseMatrix<MAT>::Elim<V>::swap(int i, int j)
727  {
728  std::swap((*rhs_)[i], (*rhs_)[j]);
729  }
730 
731  template<typename MAT>
732  template<typename V>
733  void DenseMatrix<MAT>::
734  Elim<V>::operator()(const typename V::field_type& factor, int k, int i)
735  {
736  (*rhs_)[k] -= factor*(*rhs_)[i];
737  }
738  template<typename MAT>
739  template<typename Func>
740  inline void DenseMatrix<MAT>::luDecomposition(DenseMatrix<MAT>& A, Func func) const
741  {
742  typedef typename FieldTraits<value_type>::real_type
743  real_type;
744  real_type norm = A.infinity_norm_real(); // for relative thresholds
745  real_type pivthres = std::max( FMatrixPrecision< real_type >::absolute_limit(), norm * FMatrixPrecision< real_type >::pivoting_limit() );
746  real_type singthres = std::max( FMatrixPrecision< real_type >::absolute_limit(), norm * FMatrixPrecision< real_type >::singular_limit() );
747 
748  // LU decomposition of A in A
749  for (size_type i=0; i<rows(); i++) // loop over all rows
750  {
751  typename FieldTraits<value_type>::real_type pivmax=fvmeta::absreal(A[i][i]);
752 
753  // pivoting ?
754  if (pivmax<pivthres)
755  {
756  // compute maximum of column
757  size_type imax=i;
758  typename FieldTraits<value_type>::real_type abs(0.0);
759  for (size_type k=i+1; k<rows(); k++)
760  if ((abs=fvmeta::absreal(A[k][i]))>pivmax)
761  {
762  pivmax = abs; imax = k;
763  }
764  // swap rows
765  if (imax!=i){
766  for (size_type j=0; j<rows(); j++)
767  std::swap(A[i][j],A[imax][j]);
768  func.swap(i, imax); // swap the pivot or rhs
769  }
770  }
771 
772  // singular ?
773  if (pivmax<singthres)
774  DUNE_THROW(FMatrixError,"matrix is singular");
775 
776  // eliminate
777  for (size_type k=i+1; k<rows(); k++)
778  {
779  field_type factor = A[k][i]/A[i][i];
780  A[k][i] = factor;
781  for (size_type j=i+1; j<rows(); j++)
782  A[k][j] -= factor*A[i][j];
783  func(factor, k, i);
784  }
785  }
786  }
787 
788  template<typename MAT>
789  template <class V>
790  inline void DenseMatrix<MAT>::solve(V& x, const V& b) const
791  {
792  // never mind those ifs, because they get optimized away
793  if (rows()!=cols())
794  DUNE_THROW(FMatrixError, "Can't solve for a " << rows() << "x" << cols() << " matrix!");
795 
796  if (rows()==1) {
797 
798 #ifdef DUNE_FMatrix_WITH_CHECKING
799  if (fvmeta::absreal((*this)[0][0])<FMatrixPrecision<>::absolute_limit())
800  DUNE_THROW(FMatrixError,"matrix is singular");
801 #endif
802  x[0] = b[0]/(*this)[0][0];
803 
804  }
805  else if (rows()==2) {
806 
807  field_type detinv = (*this)[0][0]*(*this)[1][1]-(*this)[0][1]*(*this)[1][0];
808 #ifdef DUNE_FMatrix_WITH_CHECKING
809  if (fvmeta::absreal(detinv)<FMatrixPrecision<>::absolute_limit())
810  DUNE_THROW(FMatrixError,"matrix is singular");
811 #endif
812  detinv = 1.0/detinv;
813 
814  x[0] = detinv*((*this)[1][1]*b[0]-(*this)[0][1]*b[1]);
815  x[1] = detinv*((*this)[0][0]*b[1]-(*this)[1][0]*b[0]);
816 
817  }
818  else if (rows()==3) {
819 
820  field_type d = determinant();
821 #ifdef DUNE_FMatrix_WITH_CHECKING
822  if (fvmeta::absreal(d)<FMatrixPrecision<>::absolute_limit())
823  DUNE_THROW(FMatrixError,"matrix is singular");
824 #endif
825 
826  x[0] = (b[0]*(*this)[1][1]*(*this)[2][2] - b[0]*(*this)[2][1]*(*this)[1][2]
827  - b[1] *(*this)[0][1]*(*this)[2][2] + b[1]*(*this)[2][1]*(*this)[0][2]
828  + b[2] *(*this)[0][1]*(*this)[1][2] - b[2]*(*this)[1][1]*(*this)[0][2]) / d;
829 
830  x[1] = ((*this)[0][0]*b[1]*(*this)[2][2] - (*this)[0][0]*b[2]*(*this)[1][2]
831  - (*this)[1][0] *b[0]*(*this)[2][2] + (*this)[1][0]*b[2]*(*this)[0][2]
832  + (*this)[2][0] *b[0]*(*this)[1][2] - (*this)[2][0]*b[1]*(*this)[0][2]) / d;
833 
834  x[2] = ((*this)[0][0]*(*this)[1][1]*b[2] - (*this)[0][0]*(*this)[2][1]*b[1]
835  - (*this)[1][0] *(*this)[0][1]*b[2] + (*this)[1][0]*(*this)[2][1]*b[0]
836  + (*this)[2][0] *(*this)[0][1]*b[1] - (*this)[2][0]*(*this)[1][1]*b[0]) / d;
837 
838  }
839  else {
840 
841  V& rhs = x; // use x to store rhs
842  rhs = b; // copy data
843  Elim<V> elim(rhs);
844  MAT A(asImp());
845 
846  luDecomposition(A, elim);
847 
848  // backsolve
849  for(int i=rows()-1; i>=0; i--){
850  for (size_type j=i+1; j<rows(); j++)
851  rhs[i] -= A[i][j]*x[j];
852  x[i] = rhs[i]/A[i][i];
853  }
854  }
855  }
856 
857  template<typename MAT>
858  inline void DenseMatrix<MAT>::invert()
859  {
860  // never mind those ifs, because they get optimized away
861  if (rows()!=cols())
862  DUNE_THROW(FMatrixError, "Can't invert a " << rows() << "x" << cols() << " matrix!");
863 
864  if (rows()==1) {
865 
866 #ifdef DUNE_FMatrix_WITH_CHECKING
867  if (fvmeta::absreal((*this)[0][0])<FMatrixPrecision<>::absolute_limit())
868  DUNE_THROW(FMatrixError,"matrix is singular");
869 #endif
870  (*this)[0][0] = 1.0/(*this)[0][0];
871 
872  }
873  else if (rows()==2) {
874 
875  field_type detinv = (*this)[0][0]*(*this)[1][1]-(*this)[0][1]*(*this)[1][0];
876 #ifdef DUNE_FMatrix_WITH_CHECKING
877  if (fvmeta::absreal(detinv)<FMatrixPrecision<>::absolute_limit())
878  DUNE_THROW(FMatrixError,"matrix is singular");
879 #endif
880  detinv = 1.0/detinv;
881 
882  field_type temp=(*this)[0][0];
883  (*this)[0][0] = (*this)[1][1]*detinv;
884  (*this)[0][1] = -(*this)[0][1]*detinv;
885  (*this)[1][0] = -(*this)[1][0]*detinv;
886  (*this)[1][1] = temp*detinv;
887 
888  }
889  else {
890 
891  MAT A(asImp());
892  std::vector<size_type> pivot(rows());
893  luDecomposition(A, ElimPivot(pivot));
894  DenseMatrix<MAT>& L=A;
895  DenseMatrix<MAT>& U=A;
896 
897  // initialize inverse
898  *this=field_type();
899 
900  for(size_type i=0; i<rows(); ++i)
901  (*this)[i][i]=1;
902 
903  // L Y = I; multiple right hand sides
904  for (size_type i=0; i<rows(); i++)
905  for (size_type j=0; j<i; j++)
906  for (size_type k=0; k<rows(); k++)
907  (*this)[i][k] -= L[i][j]*(*this)[j][k];
908 
909  // U A^{-1} = Y
910  for (size_type i=rows(); i>0;){
911  --i;
912  for (size_type k=0; k<rows(); k++){
913  for (size_type j=i+1; j<rows(); j++)
914  (*this)[i][k] -= U[i][j]*(*this)[j][k];
915  (*this)[i][k] /= U[i][i];
916  }
917  }
918 
919  for(size_type i=rows(); i>0; ){
920  --i;
921  if(i!=pivot[i])
922  for(size_type j=0; j<rows(); ++j)
923  std::swap((*this)[j][pivot[i]], (*this)[j][i]);
924  }
925  }
926  }
927 
928  // implementation of the determinant
929  template<typename MAT>
930  inline typename DenseMatrix<MAT>::field_type
931  DenseMatrix<MAT>::determinant() const
932  {
933  // never mind those ifs, because they get optimized away
934  if (rows()!=cols())
935  DUNE_THROW(FMatrixError, "There is no determinant for a " << rows() << "x" << cols() << " matrix!");
936 
937  if (rows()==1)
938  return (*this)[0][0];
939 
940  if (rows()==2)
941  return (*this)[0][0]*(*this)[1][1] - (*this)[0][1]*(*this)[1][0];
942 
943  if (rows()==3) {
944  // code generated by maple
945  field_type t4 = (*this)[0][0] * (*this)[1][1];
946  field_type t6 = (*this)[0][0] * (*this)[1][2];
947  field_type t8 = (*this)[0][1] * (*this)[1][0];
948  field_type t10 = (*this)[0][2] * (*this)[1][0];
949  field_type t12 = (*this)[0][1] * (*this)[2][0];
950  field_type t14 = (*this)[0][2] * (*this)[2][0];
951 
952  return (t4*(*this)[2][2]-t6*(*this)[2][1]-t8*(*this)[2][2]+
953  t10*(*this)[2][1]+t12*(*this)[1][2]-t14*(*this)[1][1]);
954 
955  }
956 
957  MAT A(asImp());
958  field_type det;
959  try
960  {
961  luDecomposition(A, ElimDet(det));
962  }
963  catch (FMatrixError&)
964  {
965  return 0;
966  }
967  for (size_type i = 0; i < rows(); ++i)
968  det *= A[i][i];
969  return det;
970  }
971 
972 #endif // DOXYGEN
973 
974  namespace DenseMatrixHelp {
975 #if 0
976 
977  template <typename K>
978  static inline K invertMatrix (const FieldMatrix<K,1,1> &matrix, FieldMatrix<K,1,1> &inverse)
979  {
980  inverse[0][0] = 1.0/matrix[0][0];
981  return matrix[0][0];
982  }
983 
985  template <typename K>
986  static inline K invertMatrix_retTransposed (const FieldMatrix<K,1,1> &matrix, FieldMatrix<K,1,1> &inverse)
987  {
988  return invertMatrix(matrix,inverse);
989  }
990 
991 
993  template <typename K>
994  static inline K invertMatrix (const FieldMatrix<K,2,2> &matrix, FieldMatrix<K,2,2> &inverse)
995  {
996  // code generated by maple
997  field_type det = (matrix[0][0]*matrix[1][1] - matrix[0][1]*matrix[1][0]);
998  field_type det_1 = 1.0/det;
999  inverse[0][0] = matrix[1][1] * det_1;
1000  inverse[0][1] = - matrix[0][1] * det_1;
1001  inverse[1][0] = - matrix[1][0] * det_1;
1002  inverse[1][1] = matrix[0][0] * det_1;
1003  return det;
1004  }
1005 
1008  template <typename K>
1009  static inline K invertMatrix_retTransposed (const FieldMatrix<K,2,2> &matrix, FieldMatrix<K,2,2> &inverse)
1010  {
1011  // code generated by maple
1012  field_type det = (matrix[0][0]*matrix[1][1] - matrix[0][1]*matrix[1][0]);
1013  field_type det_1 = 1.0/det;
1014  inverse[0][0] = matrix[1][1] * det_1;
1015  inverse[1][0] = - matrix[0][1] * det_1;
1016  inverse[0][1] = - matrix[1][0] * det_1;
1017  inverse[1][1] = matrix[0][0] * det_1;
1018  return det;
1019  }
1020 
1022  template <typename K>
1023  static inline K invertMatrix (const FieldMatrix<K,3,3> &matrix, FieldMatrix<K,3,3> &inverse)
1024  {
1025  // code generated by maple
1026  field_type t4 = matrix[0][0] * matrix[1][1];
1027  field_type t6 = matrix[0][0] * matrix[1][2];
1028  field_type t8 = matrix[0][1] * matrix[1][0];
1029  field_type t10 = matrix[0][2] * matrix[1][0];
1030  field_type t12 = matrix[0][1] * matrix[2][0];
1031  field_type t14 = matrix[0][2] * matrix[2][0];
1032 
1033  field_type det = (t4*matrix[2][2]-t6*matrix[2][1]-t8*matrix[2][2]+
1034  t10*matrix[2][1]+t12*matrix[1][2]-t14*matrix[1][1]);
1035  field_type t17 = 1.0/det;
1036 
1037  inverse[0][0] = (matrix[1][1] * matrix[2][2] - matrix[1][2] * matrix[2][1])*t17;
1038  inverse[0][1] = -(matrix[0][1] * matrix[2][2] - matrix[0][2] * matrix[2][1])*t17;
1039  inverse[0][2] = (matrix[0][1] * matrix[1][2] - matrix[0][2] * matrix[1][1])*t17;
1040  inverse[1][0] = -(matrix[1][0] * matrix[2][2] - matrix[1][2] * matrix[2][0])*t17;
1041  inverse[1][1] = (matrix[0][0] * matrix[2][2] - t14) * t17;
1042  inverse[1][2] = -(t6-t10) * t17;
1043  inverse[2][0] = (matrix[1][0] * matrix[2][1] - matrix[1][1] * matrix[2][0]) * t17;
1044  inverse[2][1] = -(matrix[0][0] * matrix[2][1] - t12) * t17;
1045  inverse[2][2] = (t4-t8) * t17;
1046 
1047  return det;
1048  }
1049 
1051  template <typename K>
1052  static inline K invertMatrix_retTransposed (const FieldMatrix<K,3,3> &matrix, FieldMatrix<K,3,3> &inverse)
1053  {
1054  // code generated by maple
1055  field_type t4 = matrix[0][0] * matrix[1][1];
1056  field_type t6 = matrix[0][0] * matrix[1][2];
1057  field_type t8 = matrix[0][1] * matrix[1][0];
1058  field_type t10 = matrix[0][2] * matrix[1][0];
1059  field_type t12 = matrix[0][1] * matrix[2][0];
1060  field_type t14 = matrix[0][2] * matrix[2][0];
1061 
1062  field_type det = (t4*matrix[2][2]-t6*matrix[2][1]-t8*matrix[2][2]+
1063  t10*matrix[2][1]+t12*matrix[1][2]-t14*matrix[1][1]);
1064  field_type t17 = 1.0/det;
1065 
1066  inverse[0][0] = (matrix[1][1] * matrix[2][2] - matrix[1][2] * matrix[2][1])*t17;
1067  inverse[1][0] = -(matrix[0][1] * matrix[2][2] - matrix[0][2] * matrix[2][1])*t17;
1068  inverse[2][0] = (matrix[0][1] * matrix[1][2] - matrix[0][2] * matrix[1][1])*t17;
1069  inverse[0][1] = -(matrix[1][0] * matrix[2][2] - matrix[1][2] * matrix[2][0])*t17;
1070  inverse[1][1] = (matrix[0][0] * matrix[2][2] - t14) * t17;
1071  inverse[2][1] = -(t6-t10) * t17;
1072  inverse[0][2] = (matrix[1][0] * matrix[2][1] - matrix[1][1] * matrix[2][0]) * t17;
1073  inverse[1][2] = -(matrix[0][0] * matrix[2][1] - t12) * t17;
1074  inverse[2][2] = (t4-t8) * t17;
1075 
1076  return det;
1077  }
1078 
1080  template< class K, int m, int n, int p >
1081  static inline void multMatrix ( const FieldMatrix< K, m, n > &A,
1082  const FieldMatrix< K, n, p > &B,
1083  FieldMatrix< K, m, p > &ret )
1084  {
1085  typedef typename FieldMatrix< K, m, p > :: size_type size_type;
1086 
1087  for( size_type i = 0; i < m; ++i )
1088  {
1089  for( size_type j = 0; j < p; ++j )
1090  {
1091  ret[ i ][ j ] = K( 0 );
1092  for( size_type k = 0; k < n; ++k )
1093  ret[ i ][ j ] += A[ i ][ k ] * B[ k ][ j ];
1094  }
1095  }
1096  }
1097 
1099  template <typename K, int rows, int cols>
1100  static inline void multTransposedMatrix(const FieldMatrix<K,rows,cols> &matrix, FieldMatrix<K,cols,cols>& ret)
1101  {
1102  typedef typename FieldMatrix<K,rows,cols>::size_type size_type;
1103 
1104  for(size_type i=0; i<cols(); i++)
1105  for(size_type j=0; j<cols(); j++)
1106  {
1107  ret[i][j]=0.0;
1108  for(size_type k=0; k<rows(); k++)
1109  ret[i][j]+=matrix[k][i]*matrix[k][j];
1110  }
1111  }
1112 #endif
1113 
1115  template <typename MAT, typename V1, typename V2>
1116  static inline void multAssign(const DenseMatrix<MAT> &matrix, const DenseVector<V1> & x, DenseVector<V2> & ret)
1117  {
1118  assert(x.size() == matrix.cols());
1119  assert(ret.size() == matrix.rows());
1120  typedef typename DenseMatrix<MAT>::size_type size_type;
1121 
1122  for(size_type i=0; i<matrix.rows(); ++i)
1123  {
1124  ret[i] = 0.0;
1125  for(size_type j=0; j<matrix.cols(); ++j)
1126  {
1127  ret[i] += matrix[i][j]*x[j];
1128  }
1129  }
1130  }
1131 
1132 #if 0
1133 
1134  template <typename K, int rows, int cols>
1135  static inline void multAssignTransposed( const FieldMatrix<K,rows,cols> &matrix, const FieldVector<K,rows> & x, FieldVector<K,cols> & ret)
1136  {
1137  typedef typename FieldMatrix<K,rows,cols>::size_type size_type;
1138 
1139  for(size_type i=0; i<cols(); ++i)
1140  {
1141  ret[i] = 0.0;
1142  for(size_type j=0; j<rows(); ++j)
1143  ret[i] += matrix[j][i]*x[j];
1144  }
1145  }
1146 
1148  template <typename K, int rows, int cols>
1149  static inline FieldVector<K,rows> mult(const FieldMatrix<K,rows,cols> &matrix, const FieldVector<K,cols> & x)
1150  {
1151  FieldVector<K,rows> ret;
1152  multAssign(matrix,x,ret);
1153  return ret;
1154  }
1155 
1157  template <typename K, int rows, int cols>
1158  static inline FieldVector<K,cols> multTransposed(const FieldMatrix<K,rows,cols> &matrix, const FieldVector<K,rows> & x)
1159  {
1160  FieldVector<K,cols> ret;
1161  multAssignTransposed( matrix, x, ret );
1162  return ret;
1163  }
1164 #endif
1165 
1166  } // end namespace DenseMatrixHelp
1167 
1169  template<typename MAT>
1170  std::ostream& operator<< (std::ostream& s, const DenseMatrix<MAT>& a)
1171  {
1172  for (typename DenseMatrix<MAT>::size_type i=0; i<a.rows(); i++)
1173  s << a[i] << std::endl;
1174  return s;
1175  }
1176 
1179 } // end namespace Dune
1180 
1181 #endif