Dune Core Modules (2.5.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 #ifndef DUNE_DENSEMATRIX_HH
4 #define DUNE_DENSEMATRIX_HH
5 
6 #include <cmath>
7 #include <cstddef>
8 #include <iostream>
9 #include <type_traits>
10 #include <vector>
11 
14 #include <dune/common/fvector.hh>
15 #include <dune/common/precision.hh>
16 #include <dune/common/classname.hh>
17 #include <dune/common/math.hh>
19 #include <dune/common/unused.hh>
20 
21 
22 namespace Dune
23 {
24 
25  template<typename M> class DenseMatrix;
26 
27  template<typename M>
28  struct FieldTraits< DenseMatrix<M> >
29  {
30  typedef const typename FieldTraits< typename DenseMatVecTraits<M>::value_type >::field_type field_type;
31  typedef const typename FieldTraits< typename DenseMatVecTraits<M>::value_type >::real_type real_type;
32  };
33 
34  /*
35  work around a problem of FieldMatrix/FieldVector,
36  there is no unique way to obtain the size of a class
37  */
38  template<class K, int N, int M> class FieldMatrix;
39  template<class K, int N> class FieldVector;
40  namespace {
41  template<class V>
42  struct VectorSize
43  {
44  static typename V::size_type size(const V & v) { return v.size(); }
45  };
46 
47  template<class K, int N>
48  struct VectorSize< const FieldVector<K,N> >
49  {
50  typedef FieldVector<K,N> V;
51  static typename V::size_type size(const V & v) { return N; }
52  };
53  }
54 
73  template< class DenseMatrix, class RHS >
75 
76 #ifndef DOXYGEN
77  namespace Impl
78  {
79 
80  template< class DenseMatrix, class RHS, class = void >
82  {};
83 
84  template< class DenseMatrix, class RHS >
85  class DenseMatrixAssigner< DenseMatrix, RHS, std::enable_if_t< Dune::IsNumber< RHS >::value > >
86  {
87  public:
88  static void apply ( DenseMatrix &denseMatrix, const RHS &rhs )
89  {
90  typedef typename DenseMatrix::field_type field_type;
91  std::fill( denseMatrix.begin(), denseMatrix.end(), static_cast< field_type >( rhs ) );
92  }
93  };
94 
95  template< class DenseMatrix, class RHS >
96  class DenseMatrixAssigner< DenseMatrix, RHS, std::enable_if_t< !std::is_same< typename RHS::const_iterator, void >::value > >
97  {
98  public:
99  static void apply ( DenseMatrix &denseMatrix, const RHS &rhs )
100  {
101  DUNE_ASSERT_BOUNDS(rhs.N() == denseMatrix.N());
102  DUNE_ASSERT_BOUNDS(rhs.M() == denseMatrix.M());
103  typename DenseMatrix::iterator tIt = std::begin(denseMatrix);
104  typename RHS::const_iterator sIt = std::begin(rhs);
105  for(; sIt != std::end(rhs); ++tIt, ++sIt)
106  std::copy(std::begin(*sIt), std::end(*sIt), std::begin(*tIt));
107  }
108  };
109 
110  } // namespace Impl
111 
112 
113 
114  template< class DenseMatrix, class RHS >
115  struct DenseMatrixAssigner
116  : public Impl::DenseMatrixAssigner< DenseMatrix, RHS >
117  {};
118 
119 
120  namespace Impl
121  {
122 
123  template< class DenseMatrix, class RHS >
124  std::true_type hasDenseMatrixAssigner ( DenseMatrix &, const RHS &, decltype( Dune::DenseMatrixAssigner< DenseMatrix, RHS >::apply( std::declval< DenseMatrix & >(), std::declval< const RHS & >() ) ) * = nullptr );
125 
126  std::false_type hasDenseMatrixAssigner ( ... );
127 
128  } // namespace Impl
129 
130  template< class DenseMatrix, class RHS >
131  struct HasDenseMatrixAssigner
132  : public decltype( Impl::hasDenseMatrixAssigner( std::declval< DenseMatrix & >(), std::declval< const RHS & >() ) )
133  {};
134 
135 #endif // #ifndef DOXYGEN
136 
137 
138 
140  class FMatrixError : public MathError {};
141 
152  template<typename MAT>
154  {
155  typedef DenseMatVecTraits<MAT> Traits;
156 
157  // Curiously recurring template pattern
158  MAT & asImp() { return static_cast<MAT&>(*this); }
159  const MAT & asImp() const { return static_cast<const MAT&>(*this); }
160 
161  public:
162  //===== type definitions and constants
163 
165  typedef typename Traits::derived_type derived_type;
166 
168  typedef typename Traits::value_type value_type;
169 
171  typedef typename Traits::value_type field_type;
172 
174  typedef typename Traits::value_type block_type;
175 
177  typedef typename Traits::size_type size_type;
178 
180  typedef typename Traits::row_type row_type;
181 
183  typedef typename Traits::row_reference row_reference;
184 
186  typedef typename Traits::const_row_reference const_row_reference;
187 
189  enum {
191  blocklevel = 1
192  };
193 
194  //===== access to components
195 
198  {
199  return asImp().mat_access(i);
200  }
201 
203  {
204  return asImp().mat_access(i);
205  }
206 
208  size_type size() const
209  {
210  return rows();
211  }
212 
213  //===== iterator interface to rows of the matrix
221  typedef typename std::remove_reference<row_reference>::type::Iterator ColIterator;
222 
225  {
226  return Iterator(*this,0);
227  }
228 
231  {
232  return Iterator(*this,rows());
233  }
234 
238  {
239  return Iterator(*this,rows()-1);
240  }
241 
245  {
246  return Iterator(*this,-1);
247  }
248 
256  typedef typename std::remove_reference<const_row_reference>::type::ConstIterator ConstColIterator;
257 
260  {
261  return ConstIterator(*this,0);
262  }
263 
266  {
267  return ConstIterator(*this,rows());
268  }
269 
273  {
274  return ConstIterator(*this,rows()-1);
275  }
276 
280  {
281  return ConstIterator(*this,-1);
282  }
283 
284  //===== assignment
285 
286  template< class RHS, class = std::enable_if_t< HasDenseMatrixAssigner< MAT, RHS >::value > >
287  DenseMatrix &operator= ( const RHS &rhs )
288  {
290  return *this;
291  }
292 
293  //===== vector space arithmetic
294 
296  template <class Other>
298  {
299  DUNE_ASSERT_BOUNDS(rows() == y.rows());
300  for (size_type i=0; i<rows(); i++)
301  (*this)[i] += y[i];
302  return *this;
303  }
304 
306  template <class Other>
308  {
309  DUNE_ASSERT_BOUNDS(rows() == y.rows());
310  for (size_type i=0; i<rows(); i++)
311  (*this)[i] -= y[i];
312  return *this;
313  }
314 
317  {
318  for (size_type i=0; i<rows(); i++)
319  (*this)[i] *= k;
320  return *this;
321  }
322 
325  {
326  for (size_type i=0; i<rows(); i++)
327  (*this)[i] /= k;
328  return *this;
329  }
330 
332  template <class Other>
334  {
335  DUNE_ASSERT_BOUNDS(rows() == y.rows());
336  for( size_type i = 0; i < rows(); ++i )
337  (*this)[ i ].axpy( k, y[ i ] );
338  return *this;
339  }
340 
342  template <class Other>
343  bool operator== (const DenseMatrix<Other>& y) const
344  {
345  DUNE_ASSERT_BOUNDS(rows() == y.rows());
346  for (size_type i=0; i<rows(); i++)
347  if ((*this)[i]!=y[i])
348  return false;
349  return true;
350  }
352  template <class Other>
353  bool operator!= (const DenseMatrix<Other>& y) const
354  {
355  return !operator==(y);
356  }
357 
358 
359  //===== linear maps
360 
362  template<class X, class Y>
363  void mv (const X& x, Y& y) const
364  {
365  DUNE_ASSERT_BOUNDS((void*)(&x) != (void*)(&y));
366  DUNE_ASSERT_BOUNDS(x.N() == M());
367  DUNE_ASSERT_BOUNDS(y.N() == N());
368 
369  using field_type = typename FieldTraits<Y>::field_type;
370  for (size_type i=0; i<rows(); ++i)
371  {
372  y[i] = field_type(0);
373  for (size_type j=0; j<cols(); j++)
374  y[i] += (*this)[i][j] * x[j];
375  }
376  }
377 
379  template< class X, class Y >
380  void mtv ( const X &x, Y &y ) const
381  {
382  DUNE_ASSERT_BOUNDS((void*)(&x) != (void*)(&y));
383  DUNE_ASSERT_BOUNDS(x.N() == N());
384  DUNE_ASSERT_BOUNDS(y.N() == M());
385 
386  using field_type = typename FieldTraits<Y>::field_type;
387  for(size_type i = 0; i < cols(); ++i)
388  {
389  y[i] = field_type(0);
390  for(size_type j = 0; j < rows(); ++j)
391  y[i] += (*this)[j][i] * x[j];
392  }
393  }
394 
396  template<class X, class Y>
397  void umv (const X& x, Y& y) const
398  {
399  DUNE_ASSERT_BOUNDS(x.N() == M());
400  DUNE_ASSERT_BOUNDS(y.N() == N());
401  for (size_type i=0; i<rows(); i++)
402  for (size_type j=0; j<cols(); j++)
403  y[i] += (*this)[i][j] * x[j];
404  }
405 
407  template<class X, class Y>
408  void umtv (const X& x, Y& y) const
409  {
410  DUNE_ASSERT_BOUNDS(x.N() == N());
411  DUNE_ASSERT_BOUNDS(y.N() == M());
412  for (size_type i=0; i<rows(); i++)
413  for (size_type j=0; j<cols(); j++)
414  y[j] += (*this)[i][j]*x[i];
415  }
416 
418  template<class X, class Y>
419  void umhv (const X& x, Y& y) const
420  {
421  DUNE_ASSERT_BOUNDS(x.N() == N());
422  DUNE_ASSERT_BOUNDS(y.N() == M());
423  for (size_type i=0; i<rows(); i++)
424  for (size_type j=0; j<cols(); j++)
425  y[j] += conjugateComplex((*this)[i][j])*x[i];
426  }
427 
429  template<class X, class Y>
430  void mmv (const X& x, Y& y) const
431  {
432  DUNE_ASSERT_BOUNDS(x.N() == M());
433  DUNE_ASSERT_BOUNDS(y.N() == N());
434  for (size_type i=0; i<rows(); i++)
435  for (size_type j=0; j<cols(); j++)
436  y[i] -= (*this)[i][j] * x[j];
437  }
438 
440  template<class X, class Y>
441  void mmtv (const X& x, Y& y) const
442  {
443  DUNE_ASSERT_BOUNDS(x.N() == N());
444  DUNE_ASSERT_BOUNDS(y.N() == M());
445  for (size_type i=0; i<rows(); i++)
446  for (size_type j=0; j<cols(); j++)
447  y[j] -= (*this)[i][j]*x[i];
448  }
449 
451  template<class X, class Y>
452  void mmhv (const X& x, Y& y) const
453  {
454  DUNE_ASSERT_BOUNDS(x.N() == N());
455  DUNE_ASSERT_BOUNDS(y.N() == M());
456  for (size_type i=0; i<rows(); i++)
457  for (size_type j=0; j<cols(); j++)
458  y[j] -= conjugateComplex((*this)[i][j])*x[i];
459  }
460 
462  template<class X, class Y>
463  void usmv (const typename FieldTraits<Y>::field_type & alpha,
464  const X& x, Y& y) const
465  {
466  DUNE_ASSERT_BOUNDS(x.N() == M());
467  DUNE_ASSERT_BOUNDS(y.N() == N());
468  for (size_type i=0; i<rows(); i++)
469  for (size_type j=0; j<cols(); j++)
470  y[i] += alpha * (*this)[i][j] * x[j];
471  }
472 
474  template<class X, class Y>
475  void usmtv (const typename FieldTraits<Y>::field_type & alpha,
476  const X& x, Y& y) const
477  {
478  DUNE_ASSERT_BOUNDS(x.N() == N());
479  DUNE_ASSERT_BOUNDS(y.N() == M());
480  for (size_type i=0; i<rows(); i++)
481  for (size_type j=0; j<cols(); j++)
482  y[j] += alpha*(*this)[i][j]*x[i];
483  }
484 
486  template<class X, class Y>
487  void usmhv (const typename FieldTraits<Y>::field_type & alpha,
488  const X& x, Y& y) const
489  {
490  DUNE_ASSERT_BOUNDS(x.N() == N());
491  DUNE_ASSERT_BOUNDS(y.N() == M());
492  for (size_type i=0; i<rows(); i++)
493  for (size_type j=0; j<cols(); j++)
494  y[j] += alpha*conjugateComplex((*this)[i][j])*x[i];
495  }
496 
497  //===== norms
498 
500  typename FieldTraits<value_type>::real_type frobenius_norm () const
501  {
502  typename FieldTraits<value_type>::real_type sum=(0.0);
503  for (size_type i=0; i<rows(); ++i) sum += (*this)[i].two_norm2();
504  return fvmeta::sqrt(sum);
505  }
506 
508  typename FieldTraits<value_type>::real_type frobenius_norm2 () const
509  {
510  typename FieldTraits<value_type>::real_type sum=(0.0);
511  for (size_type i=0; i<rows(); ++i) sum += (*this)[i].two_norm2();
512  return sum;
513  }
514 
516  template <typename vt = value_type,
517  typename std::enable_if<!has_nan<vt>::value, int>::type = 0>
518  typename FieldTraits<vt>::real_type infinity_norm() const {
519  using real_type = typename FieldTraits<vt>::real_type;
520  using std::max;
521 
522  real_type norm = 0;
523  for (auto const &x : *this) {
524  real_type const a = x.one_norm();
525  norm = max(a, norm);
526  }
527  return norm;
528  }
529 
531  template <typename vt = value_type,
532  typename std::enable_if<!has_nan<vt>::value, int>::type = 0>
533  typename FieldTraits<vt>::real_type infinity_norm_real() const {
534  using real_type = typename FieldTraits<vt>::real_type;
535  using std::max;
536 
537  real_type norm = 0;
538  for (auto const &x : *this) {
539  real_type const a = x.one_norm_real();
540  norm = max(a, norm);
541  }
542  return norm;
543  }
544 
546  template <typename vt = value_type,
547  typename std::enable_if<has_nan<vt>::value, int>::type = 0>
548  typename FieldTraits<vt>::real_type infinity_norm() const {
549  using real_type = typename FieldTraits<vt>::real_type;
550  using std::max;
551 
552  real_type norm = 0;
553  real_type isNaN = 1;
554  for (auto const &x : *this) {
555  real_type const a = x.one_norm();
556  norm = max(a, norm);
557  isNaN += a;
558  }
559  isNaN /= isNaN;
560  return norm * isNaN;
561  }
562 
564  template <typename vt = value_type,
565  typename std::enable_if<has_nan<vt>::value, int>::type = 0>
566  typename FieldTraits<vt>::real_type infinity_norm_real() const {
567  using real_type = typename FieldTraits<vt>::real_type;
568  using std::max;
569 
570  real_type norm = 0;
571  real_type isNaN = 1;
572  for (auto const &x : *this) {
573  real_type const a = x.one_norm_real();
574  norm = max(a, norm);
575  isNaN += a;
576  }
577  isNaN /= isNaN;
578  return norm * isNaN;
579  }
580 
581  //===== solve
582 
587  template <class V>
588  void solve (V& x, const V& b) const;
589 
594  void invert();
595 
598 
600  template<typename M2>
602  {
603  DUNE_ASSERT_BOUNDS(M.rows() == M.cols());
604  DUNE_ASSERT_BOUNDS(M.rows() == rows());
605  MAT C(asImp());
606 
607  for (size_type i=0; i<rows(); i++)
608  for (size_type j=0; j<cols(); j++) {
609  (*this)[i][j] = 0;
610  for (size_type k=0; k<rows(); k++)
611  (*this)[i][j] += M[i][k]*C[k][j];
612  }
613 
614  return asImp();
615  }
616 
618  template<typename M2>
620  {
621  DUNE_ASSERT_BOUNDS(M.rows() == M.cols());
622  DUNE_ASSERT_BOUNDS(M.cols() == cols());
623  MAT C(asImp());
624 
625  for (size_type i=0; i<rows(); i++)
626  for (size_type j=0; j<cols(); j++) {
627  (*this)[i][j] = 0;
628  for (size_type k=0; k<cols(); k++)
629  (*this)[i][j] += C[i][k]*M[k][j];
630  }
631  return asImp();
632  }
633 
634 #if 0
636  template<int l>
637  DenseMatrix<K,l,cols> leftmultiplyany (const FieldMatrix<K,l,rows>& M) const
638  {
640 
641  for (size_type i=0; i<l; i++) {
642  for (size_type j=0; j<cols(); j++) {
643  C[i][j] = 0;
644  for (size_type k=0; k<rows(); k++)
645  C[i][j] += M[i][k]*(*this)[k][j];
646  }
647  }
648  return C;
649  }
650 
652  template<int l>
653  FieldMatrix<K,rows,l> rightmultiplyany (const FieldMatrix<K,cols,l>& M) const
654  {
655  FieldMatrix<K,rows,l> C;
656 
657  for (size_type i=0; i<rows(); i++) {
658  for (size_type j=0; j<l; j++) {
659  C[i][j] = 0;
660  for (size_type k=0; k<cols(); k++)
661  C[i][j] += (*this)[i][k]*M[k][j];
662  }
663  }
664  return C;
665  }
666 #endif
667 
668  //===== sizes
669 
671  size_type N () const
672  {
673  return rows();
674  }
675 
677  size_type M () const
678  {
679  return cols();
680  }
681 
683  size_type rows() const
684  {
685  return asImp().mat_rows();
686  }
687 
689  size_type cols() const
690  {
691  return asImp().mat_cols();
692  }
693 
694  //===== query
695 
697  bool exists (size_type i, size_type j) const
698  {
701  DUNE_ASSERT_BOUNDS(i >= 0 && i < rows());
702  DUNE_ASSERT_BOUNDS(j >= 0 && j < cols());
703  return true;
704  }
705 
706  private:
707 
708 #ifndef DOXYGEN
709  struct ElimPivot
710  {
711  ElimPivot(std::vector<size_type> & pivot);
712 
713  void swap(int i, int j);
714 
715  template<typename T>
716  void operator()(const T&, int, int)
717  {}
718 
719  std::vector<size_type> & pivot_;
720  };
721 
722  template<typename V>
723  struct Elim
724  {
725  Elim(V& rhs);
726 
727  void swap(int i, int j);
728 
729  void operator()(const typename V::field_type& factor, int k, int i);
730 
731  V* rhs_;
732  };
733 
734  struct ElimDet
735  {
736  ElimDet(field_type& sign) : sign_(sign)
737  { sign_ = 1; }
738 
739  void swap(int, int)
740  { sign_ *= -1; }
741 
742  void operator()(const field_type&, int, int)
743  {}
744 
745  field_type& sign_;
746  };
747 #endif // DOXYGEN
748 
749  template<class Func>
750  void luDecomposition(DenseMatrix<MAT>& A, Func func) const;
751  };
752 
753 #ifndef DOXYGEN
754  template<typename MAT>
755  DenseMatrix<MAT>::ElimPivot::ElimPivot(std::vector<size_type> & pivot)
756  : pivot_(pivot)
757  {
758  typedef typename std::vector<size_type>::size_type size_type;
759  for(size_type i=0; i < pivot_.size(); ++i) pivot_[i]=i;
760  }
761 
762  template<typename MAT>
763  void DenseMatrix<MAT>::ElimPivot::swap(int i, int j)
764  {
765  pivot_[i]=j;
766  }
767 
768  template<typename MAT>
769  template<typename V>
770  DenseMatrix<MAT>::Elim<V>::Elim(V& rhs)
771  : rhs_(&rhs)
772  {}
773 
774  template<typename MAT>
775  template<typename V>
776  void DenseMatrix<MAT>::Elim<V>::swap(int i, int j)
777  {
778  std::swap((*rhs_)[i], (*rhs_)[j]);
779  }
780 
781  template<typename MAT>
782  template<typename V>
783  void DenseMatrix<MAT>::
784  Elim<V>::operator()(const typename V::field_type& factor, int k, int i)
785  {
786  (*rhs_)[k] -= factor*(*rhs_)[i];
787  }
788  template<typename MAT>
789  template<typename Func>
790  inline void DenseMatrix<MAT>::luDecomposition(DenseMatrix<MAT>& A, Func func) const
791  {
792  typedef typename FieldTraits<value_type>::real_type
793  real_type;
794  real_type norm = A.infinity_norm_real(); // for relative thresholds
797 
798  // LU decomposition of A in A
799  for (size_type i=0; i<rows(); i++) // loop over all rows
800  {
801  typename FieldTraits<value_type>::real_type pivmax=fvmeta::absreal(A[i][i]);
802 
803  // pivoting ?
804  if (pivmax<pivthres)
805  {
806  // compute maximum of column
807  size_type imax=i;
808  typename FieldTraits<value_type>::real_type abs(0.0);
809  for (size_type k=i+1; k<rows(); k++)
810  if ((abs=fvmeta::absreal(A[k][i]))>pivmax)
811  {
812  pivmax = abs; imax = k;
813  }
814  // swap rows
815  if (imax!=i) {
816  for (size_type j=0; j<rows(); j++)
817  std::swap(A[i][j],A[imax][j]);
818  func.swap(i, imax); // swap the pivot or rhs
819  }
820  }
821 
822  // singular ?
823  if (pivmax<singthres)
824  DUNE_THROW(FMatrixError,"matrix is singular");
825 
826  // eliminate
827  for (size_type k=i+1; k<rows(); k++)
828  {
829  field_type factor = A[k][i]/A[i][i];
830  A[k][i] = factor;
831  for (size_type j=i+1; j<rows(); j++)
832  A[k][j] -= factor*A[i][j];
833  func(factor, k, i);
834  }
835  }
836  }
837 
838  template<typename MAT>
839  template <class V>
840  inline void DenseMatrix<MAT>::solve(V& x, const V& b) const
841  {
842  // never mind those ifs, because they get optimized away
843  if (rows()!=cols())
844  DUNE_THROW(FMatrixError, "Can't solve for a " << rows() << "x" << cols() << " matrix!");
845 
846  if (rows()==1) {
847 
848 #ifdef DUNE_FMatrix_WITH_CHECKING
849  if (fvmeta::absreal((*this)[0][0])<FMatrixPrecision<>::absolute_limit())
850  DUNE_THROW(FMatrixError,"matrix is singular");
851 #endif
852  x[0] = b[0]/(*this)[0][0];
853 
854  }
855  else if (rows()==2) {
856 
857  field_type detinv = (*this)[0][0]*(*this)[1][1]-(*this)[0][1]*(*this)[1][0];
858 #ifdef DUNE_FMatrix_WITH_CHECKING
859  if (fvmeta::absreal(detinv)<FMatrixPrecision<>::absolute_limit())
860  DUNE_THROW(FMatrixError,"matrix is singular");
861 #endif
862  detinv = 1.0/detinv;
863 
864  x[0] = detinv*((*this)[1][1]*b[0]-(*this)[0][1]*b[1]);
865  x[1] = detinv*((*this)[0][0]*b[1]-(*this)[1][0]*b[0]);
866 
867  }
868  else if (rows()==3) {
869 
870  field_type d = determinant();
871 #ifdef DUNE_FMatrix_WITH_CHECKING
872  if (fvmeta::absreal(d)<FMatrixPrecision<>::absolute_limit())
873  DUNE_THROW(FMatrixError,"matrix is singular");
874 #endif
875 
876  x[0] = (b[0]*(*this)[1][1]*(*this)[2][2] - b[0]*(*this)[2][1]*(*this)[1][2]
877  - b[1] *(*this)[0][1]*(*this)[2][2] + b[1]*(*this)[2][1]*(*this)[0][2]
878  + b[2] *(*this)[0][1]*(*this)[1][2] - b[2]*(*this)[1][1]*(*this)[0][2]) / d;
879 
880  x[1] = ((*this)[0][0]*b[1]*(*this)[2][2] - (*this)[0][0]*b[2]*(*this)[1][2]
881  - (*this)[1][0] *b[0]*(*this)[2][2] + (*this)[1][0]*b[2]*(*this)[0][2]
882  + (*this)[2][0] *b[0]*(*this)[1][2] - (*this)[2][0]*b[1]*(*this)[0][2]) / d;
883 
884  x[2] = ((*this)[0][0]*(*this)[1][1]*b[2] - (*this)[0][0]*(*this)[2][1]*b[1]
885  - (*this)[1][0] *(*this)[0][1]*b[2] + (*this)[1][0]*(*this)[2][1]*b[0]
886  + (*this)[2][0] *(*this)[0][1]*b[1] - (*this)[2][0]*(*this)[1][1]*b[0]) / d;
887 
888  }
889  else {
890 
891  V& rhs = x; // use x to store rhs
892  rhs = b; // copy data
893  Elim<V> elim(rhs);
894  MAT A(asImp());
895 
896  luDecomposition(A, elim);
897 
898  // backsolve
899  for(int i=rows()-1; i>=0; i--) {
900  for (size_type j=i+1; j<rows(); j++)
901  rhs[i] -= A[i][j]*x[j];
902  x[i] = rhs[i]/A[i][i];
903  }
904  }
905  }
906 
907  template<typename MAT>
908  inline void DenseMatrix<MAT>::invert()
909  {
910  // never mind those ifs, because they get optimized away
911  if (rows()!=cols())
912  DUNE_THROW(FMatrixError, "Can't invert a " << rows() << "x" << cols() << " matrix!");
913 
914  if (rows()==1) {
915 
916 #ifdef DUNE_FMatrix_WITH_CHECKING
917  if (fvmeta::absreal((*this)[0][0])<FMatrixPrecision<>::absolute_limit())
918  DUNE_THROW(FMatrixError,"matrix is singular");
919 #endif
920  (*this)[0][0] = field_type( 1 ) / (*this)[0][0];
921 
922  }
923  else if (rows()==2) {
924 
925  field_type detinv = (*this)[0][0]*(*this)[1][1]-(*this)[0][1]*(*this)[1][0];
926 #ifdef DUNE_FMatrix_WITH_CHECKING
927  if (fvmeta::absreal(detinv)<FMatrixPrecision<>::absolute_limit())
928  DUNE_THROW(FMatrixError,"matrix is singular");
929 #endif
930  detinv = field_type( 1 ) / detinv;
931 
932  field_type temp=(*this)[0][0];
933  (*this)[0][0] = (*this)[1][1]*detinv;
934  (*this)[0][1] = -(*this)[0][1]*detinv;
935  (*this)[1][0] = -(*this)[1][0]*detinv;
936  (*this)[1][1] = temp*detinv;
937 
938  }
939  else if (rows()==3)
940  {
941  using K = field_type;
942  // code generated by maple
943  K t4 = (*this)[0][0] * (*this)[1][1];
944  K t6 = (*this)[0][0] * (*this)[1][2];
945  K t8 = (*this)[0][1] * (*this)[1][0];
946  K t10 = (*this)[0][2] * (*this)[1][0];
947  K t12 = (*this)[0][1] * (*this)[2][0];
948  K t14 = (*this)[0][2] * (*this)[2][0];
949 
950  K det = (t4*(*this)[2][2]-t6*(*this)[2][1]-t8*(*this)[2][2]+
951  t10*(*this)[2][1]+t12*(*this)[1][2]-t14*(*this)[1][1]);
952  K t17 = K(1.0)/det;
953 
954  K matrix01 = (*this)[0][1];
955  K matrix00 = (*this)[0][0];
956  K matrix10 = (*this)[1][0];
957  K matrix11 = (*this)[1][1];
958 
959  (*this)[0][0] = ((*this)[1][1] * (*this)[2][2] - (*this)[1][2] * (*this)[2][1])*t17;
960  (*this)[0][1] = -((*this)[0][1] * (*this)[2][2] - (*this)[0][2] * (*this)[2][1])*t17;
961  (*this)[0][2] = (matrix01 * (*this)[1][2] - (*this)[0][2] * (*this)[1][1])*t17;
962  (*this)[1][0] = -((*this)[1][0] * (*this)[2][2] - (*this)[1][2] * (*this)[2][0])*t17;
963  (*this)[1][1] = (matrix00 * (*this)[2][2] - t14) * t17;
964  (*this)[1][2] = -(t6-t10) * t17;
965  (*this)[2][0] = (matrix10 * (*this)[2][1] - matrix11 * (*this)[2][0]) * t17;
966  (*this)[2][1] = -(matrix00 * (*this)[2][1] - t12) * t17;
967  (*this)[2][2] = (t4-t8) * t17;
968  }
969  else {
970 
971  MAT A(asImp());
972  std::vector<size_type> pivot(rows());
973  luDecomposition(A, ElimPivot(pivot));
974  DenseMatrix<MAT>& L=A;
975  DenseMatrix<MAT>& U=A;
976 
977  // initialize inverse
978  *this=field_type();
979 
980  for(size_type i=0; i<rows(); ++i)
981  (*this)[i][i]=1;
982 
983  // L Y = I; multiple right hand sides
984  for (size_type i=0; i<rows(); i++)
985  for (size_type j=0; j<i; j++)
986  for (size_type k=0; k<rows(); k++)
987  (*this)[i][k] -= L[i][j]*(*this)[j][k];
988 
989  // U A^{-1} = Y
990  for (size_type i=rows(); i>0;) {
991  --i;
992  for (size_type k=0; k<rows(); k++) {
993  for (size_type j=i+1; j<rows(); j++)
994  (*this)[i][k] -= U[i][j]*(*this)[j][k];
995  (*this)[i][k] /= U[i][i];
996  }
997  }
998 
999  for(size_type i=rows(); i>0; ) {
1000  --i;
1001  if(i!=pivot[i])
1002  for(size_type j=0; j<rows(); ++j)
1003  std::swap((*this)[j][pivot[i]], (*this)[j][i]);
1004  }
1005  }
1006  }
1007 
1008  // implementation of the determinant
1009  template<typename MAT>
1010  inline typename DenseMatrix<MAT>::field_type
1012  {
1013  // never mind those ifs, because they get optimized away
1014  if (rows()!=cols())
1015  DUNE_THROW(FMatrixError, "There is no determinant for a " << rows() << "x" << cols() << " matrix!");
1016 
1017  if (rows()==1)
1018  return (*this)[0][0];
1019 
1020  if (rows()==2)
1021  return (*this)[0][0]*(*this)[1][1] - (*this)[0][1]*(*this)[1][0];
1022 
1023  if (rows()==3) {
1024  // code generated by maple
1025  field_type t4 = (*this)[0][0] * (*this)[1][1];
1026  field_type t6 = (*this)[0][0] * (*this)[1][2];
1027  field_type t8 = (*this)[0][1] * (*this)[1][0];
1028  field_type t10 = (*this)[0][2] * (*this)[1][0];
1029  field_type t12 = (*this)[0][1] * (*this)[2][0];
1030  field_type t14 = (*this)[0][2] * (*this)[2][0];
1031 
1032  return (t4*(*this)[2][2]-t6*(*this)[2][1]-t8*(*this)[2][2]+
1033  t10*(*this)[2][1]+t12*(*this)[1][2]-t14*(*this)[1][1]);
1034 
1035  }
1036 
1037  MAT A(asImp());
1038  field_type det;
1039  try
1040  {
1041  luDecomposition(A, ElimDet(det));
1042  }
1043  catch (FMatrixError&)
1044  {
1045  return 0;
1046  }
1047  for (size_type i = 0; i < rows(); ++i)
1048  det *= A[i][i];
1049  return det;
1050  }
1051 
1052 #endif // DOXYGEN
1053 
1054  namespace DenseMatrixHelp {
1055 
1057  template <typename MAT, typename V1, typename V2>
1058  static inline void multAssign(const DenseMatrix<MAT> &matrix, const DenseVector<V1> & x, DenseVector<V2> & ret)
1059  {
1060  DUNE_ASSERT_BOUNDS(x.size() == matrix.cols());
1061  DUNE_ASSERT_BOUNDS(ret.size() == matrix.rows());
1062  typedef typename DenseMatrix<MAT>::size_type size_type;
1063 
1064  for(size_type i=0; i<matrix.rows(); ++i)
1065  {
1066  ret[i] = 0.0;
1067  for(size_type j=0; j<matrix.cols(); ++j)
1068  {
1069  ret[i] += matrix[i][j]*x[j];
1070  }
1071  }
1072  }
1073 
1074 #if 0
1076  template <typename K, int rows, int cols>
1077  static inline void multAssignTransposed( const FieldMatrix<K,rows,cols> &matrix, const FieldVector<K,rows> & x, FieldVector<K,cols> & ret)
1078  {
1079  typedef typename FieldMatrix<K,rows,cols>::size_type size_type;
1080 
1081  for(size_type i=0; i<cols(); ++i)
1082  {
1083  ret[i] = 0.0;
1084  for(size_type j=0; j<rows(); ++j)
1085  ret[i] += matrix[j][i]*x[j];
1086  }
1087  }
1088 
1090  template <typename K, int rows, int cols>
1091  static inline FieldVector<K,rows> mult(const FieldMatrix<K,rows,cols> &matrix, const FieldVector<K,cols> & x)
1092  {
1093  FieldVector<K,rows> ret;
1094  multAssign(matrix,x,ret);
1095  return ret;
1096  }
1097 
1099  template <typename K, int rows, int cols>
1100  static inline FieldVector<K,cols> multTransposed(const FieldMatrix<K,rows,cols> &matrix, const FieldVector<K,rows> & x)
1101  {
1102  FieldVector<K,cols> ret;
1103  multAssignTransposed( matrix, x, ret );
1104  return ret;
1105  }
1106 #endif
1107 
1108  } // end namespace DenseMatrixHelp
1109 
1111  template<typename MAT>
1112  std::ostream& operator<< (std::ostream& s, const DenseMatrix<MAT>& a)
1113  {
1114  for (typename DenseMatrix<MAT>::size_type i=0; i<a.rows(); i++)
1115  s << a[i] << std::endl;
1116  return s;
1117  }
1118 
1121 } // end namespace Dune
1122 
1123 #endif
Macro for wrapping boundary checks.
Generic iterator class for dense vector and matrix implementations.
Definition: densevector.hh:128
A dense n x m matrix.
Definition: densematrix.hh:154
ConstIterator const_iterator
typedef for stl compliant access
Definition: densematrix.hh:252
bool operator!=(const DenseMatrix< Other > &y) const
Binary matrix incomparison.
Definition: densematrix.hh:353
void mv(const X &x, Y &y) const
y = A x
Definition: densematrix.hh:363
Traits::value_type field_type
export the type representing the field
Definition: densematrix.hh:171
void solve(V &x, const V &b) const
Solve system A x = b.
size_type cols() const
number of columns
Definition: densematrix.hh:689
ConstIterator beforeEnd() const
Definition: densematrix.hh:272
bool operator==(const DenseMatrix< Other > &y) const
Binary matrix comparison.
Definition: densematrix.hh:343
void mmtv(const X &x, Y &y) const
y -= A^T x
Definition: densematrix.hh:441
std::remove_reference< const_row_reference >::type::ConstIterator ConstColIterator
rename the iterators for easier access
Definition: densematrix.hh:256
ConstIterator beforeBegin() const
Definition: densematrix.hh:279
Traits::value_type block_type
export the type representing the components
Definition: densematrix.hh:174
DenseMatrix & operator+=(const DenseMatrix< Other > &y)
vector space addition
Definition: densematrix.hh:297
void mtv(const X &x, Y &y) const
y = A^T x
Definition: densematrix.hh:380
FieldTraits< value_type >::real_type frobenius_norm2() const
square of frobenius norm, need for block recursion
Definition: densematrix.hh:508
MAT & rightmultiply(const DenseMatrix< M2 > &M)
Multiplies M from the right to this matrix.
Definition: densematrix.hh:619
size_type size() const
size method (number of rows)
Definition: densematrix.hh:208
size_type M() const
number of columns
Definition: densematrix.hh:677
Iterator end()
end iterator
Definition: densematrix.hh:230
Iterator beforeBegin()
Definition: densematrix.hh:244
FieldTraits< vt >::real_type infinity_norm_real() const
simplified infinity norm (uses Manhattan norm for complex values)
Definition: densematrix.hh:533
DenseMatrix & operator-=(const DenseMatrix< Other > &y)
vector space subtraction
Definition: densematrix.hh:307
Iterator beforeEnd()
Definition: densematrix.hh:237
Traits::value_type value_type
export the type representing the field
Definition: densematrix.hh:168
void usmv(const typename FieldTraits< Y >::field_type &alpha, const X &x, Y &y) const
y += alpha A x
Definition: densematrix.hh:463
void usmhv(const typename FieldTraits< Y >::field_type &alpha, const X &x, Y &y) const
y += alpha A^H x
Definition: densematrix.hh:487
size_type rows() const
number of rows
Definition: densematrix.hh:683
void mmv(const X &x, Y &y) const
y -= A x
Definition: densematrix.hh:430
DenseMatrix & axpy(const field_type &k, const DenseMatrix< Other > &y)
vector space axpy operation (*this += k y)
Definition: densematrix.hh:333
void invert()
Compute inverse.
FieldTraits< value_type >::real_type frobenius_norm() const
frobenius norm: sqrt(sum over squared values of entries)
Definition: densematrix.hh:500
void usmtv(const typename FieldTraits< Y >::field_type &alpha, const X &x, Y &y) const
y += alpha A^T x
Definition: densematrix.hh:475
DenseMatrix & operator*=(const field_type &k)
vector space multiplication with scalar
Definition: densematrix.hh:316
void mmhv(const X &x, Y &y) const
y -= A^H x
Definition: densematrix.hh:452
Traits::derived_type derived_type
type of derived matrix class
Definition: densematrix.hh:165
row_reference operator[](size_type i)
random access
Definition: densematrix.hh:197
bool exists(size_type i, size_type j) const
return true when (i,j) is in pattern
Definition: densematrix.hh:697
Iterator RowIterator
rename the iterators for easier access
Definition: densematrix.hh:219
void umv(const X &x, Y &y) const
y += A x
Definition: densematrix.hh:397
DenseIterator< const DenseMatrix, const row_type, const_row_reference > ConstIterator
Iterator class for sequential access.
Definition: densematrix.hh:250
size_type N() const
number of rows
Definition: densematrix.hh:671
Traits::row_type row_type
The type used to represent a row (must fulfill the Dune::DenseVector interface)
Definition: densematrix.hh:180
DenseMatrix & operator/=(const field_type &k)
vector space division by scalar
Definition: densematrix.hh:324
Traits::size_type size_type
The type used for the index access and size operation.
Definition: densematrix.hh:177
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:186
std::remove_reference< row_reference >::type::Iterator ColIterator
rename the iterators for easier access
Definition: densematrix.hh:221
Traits::row_reference row_reference
The type used to represent a reference to a row (usually row_type &)
Definition: densematrix.hh:183
MAT & leftmultiply(const DenseMatrix< M2 > &M)
Multiplies M from the left to this matrix.
Definition: densematrix.hh:601
Iterator iterator
typedef for stl compliant access
Definition: densematrix.hh:217
ConstIterator ConstRowIterator
rename the iterators for easier access
Definition: densematrix.hh:254
field_type determinant() const
calculates the determinant of this matrix
DenseIterator< DenseMatrix, row_type, row_reference > Iterator
Iterator class for sequential access.
Definition: densematrix.hh:215
void umtv(const X &x, Y &y) const
y += A^T x
Definition: densematrix.hh:408
ConstIterator begin() const
begin iterator
Definition: densematrix.hh:259
@ blocklevel
The number of block levels we contain. This is 1.
Definition: densematrix.hh:191
Iterator begin()
begin iterator
Definition: densematrix.hh:224
void umhv(const X &x, Y &y) const
y += A^H x
Definition: densematrix.hh:419
FieldTraits< vt >::real_type infinity_norm() const
infinity norm (row sum norm, how to generalize for blocks?)
Definition: densematrix.hh:518
ConstIterator end() const
end iterator
Definition: densematrix.hh:265
Interface for a class of dense vectors over a given field.
Definition: densevector.hh:235
size_type size() const
size method
Definition: densevector.hh:297
size_type N() const
number of blocks in the vector (are of size 1 here)
Definition: densevector.hh:678
Error thrown if operations of a FieldMatrix fail.
Definition: densematrix.hh:140
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:68
vector space out of a tensor product of fields.
Definition: fvector.hh:93
constexpr FieldVector()
Constructor making default-initialized vector.
Definition: fvector.hh:113
Default exception class for mathematical errors.
Definition: exceptions.hh:239
FieldTraits< ft >::real_type infinity_norm_real() const
simplified infinity norm (uses Manhattan norm for complex values)
Definition: matrix.hh:1006
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:1058
A few common exception classes.
Implements a vector constructed from a given type representing a field and a compile-time given size.
#define DUNE_ASSERT_BOUNDS(cond)
If DUNE_CHECK_BOUNDS is defined: check if condition cond holds; otherwise, do nothing.
Definition: boundschecking.hh:28
#define DUNE_THROW(E, m)
Definition: exceptions.hh:216
Some useful basic math stuff.
Dune namespace.
Definition: alignment.hh:11
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.
you have to specialize this structure for any type that should be assignable to a DenseMatrix
Definition: densematrix.hh:74
Traits for type conversions and type information.
Definition of the DUNE_UNUSED macro for the case that config.h is not available.
#define DUNE_UNUSED_PARAMETER(parm)
A macro to mark intentionally unused function parameters with.
Definition: unused.hh:18
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (May 9, 22:29, 2024)