Dune Core Modules (2.7.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 <utility>
11 #include <vector>
12 
14 #include <dune/common/classname.hh>
17 #include <dune/common/fvector.hh>
18 #include <dune/common/math.hh>
19 #include <dune/common/precision.hh>
20 #include <dune/common/simd/simd.hh>
22 #include <dune/common/unused.hh>
24 
25 namespace Dune
26 {
27 
28  template<typename M> class DenseMatrix;
29 
30  template<typename M>
31  struct FieldTraits< DenseMatrix<M> >
32  {
33  typedef const typename FieldTraits< typename DenseMatVecTraits<M>::value_type >::field_type field_type;
34  typedef const typename FieldTraits< typename DenseMatVecTraits<M>::value_type >::real_type real_type;
35  };
36 
37  /*
38  work around a problem of FieldMatrix/FieldVector,
39  there is no unique way to obtain the size of a class
40  */
41  template<class K, int N, int M> class FieldMatrix;
42  template<class K, int N> class FieldVector;
43  namespace {
44  template<class V>
45  struct DUNE_DEPRECATED_MSG("VectorSize is deprecated; please call the 'size()' method directly instead") VectorSize
46  {
47  static typename V::size_type size(const V & v) { return v.size(); }
48  };
49 
50  DUNE_NO_DEPRECATED_BEGIN
51  template<class K, int N>
52  struct DUNE_DEPRECATED_MSG("VectorSize is deprecated; please call the 'size()' method directly instead") VectorSize< const FieldVector<K,N> >
53  {
54  typedef FieldVector<K,N> V;
55  static typename V::size_type size(const V & v)
56  {
57  DUNE_UNUSED_PARAMETER(v);
58  return N;
59  }
60  };
61  DUNE_NO_DEPRECATED_END
62  }
63 
82  template< class DenseMatrix, class RHS >
83  struct DenseMatrixAssigner;
84 
85 #ifndef DOXYGEN
86  namespace Impl
87  {
88 
89  template< class DenseMatrix, class RHS, class = void >
90  class DenseMatrixAssigner
91  {};
92 
93  template< class DenseMatrix, class RHS >
94  class DenseMatrixAssigner< DenseMatrix, RHS, std::enable_if_t< Dune::IsNumber< RHS >::value > >
95  {
96  public:
97  static void apply ( DenseMatrix &denseMatrix, const RHS &rhs )
98  {
99  typedef typename DenseMatrix::field_type field_type;
100  std::fill( denseMatrix.begin(), denseMatrix.end(), static_cast< field_type >( rhs ) );
101  }
102  };
103 
104  template< class DenseMatrix, class RHS >
105  class DenseMatrixAssigner< DenseMatrix, RHS, std::enable_if_t< !std::is_same< typename RHS::const_iterator, void >::value
106  && std::is_convertible< typename RHS::const_iterator::value_type, typename DenseMatrix::iterator::value_type >::value > >
107  {
108  public:
109  static void apply ( DenseMatrix &denseMatrix, const RHS &rhs )
110  {
111  DUNE_ASSERT_BOUNDS(rhs.N() == denseMatrix.N());
112  DUNE_ASSERT_BOUNDS(rhs.M() == denseMatrix.M());
113  typename DenseMatrix::iterator tIt = std::begin(denseMatrix);
114  typename RHS::const_iterator sIt = std::begin(rhs);
115  for(; sIt != std::end(rhs); ++tIt, ++sIt)
116  std::copy(std::begin(*sIt), std::end(*sIt), std::begin(*tIt));
117  }
118  };
119 
120  } // namespace Impl
121 
122 
123 
124  template< class DenseMatrix, class RHS >
125  struct DenseMatrixAssigner
126  : public Impl::DenseMatrixAssigner< DenseMatrix, RHS >
127  {};
128 
129 
130  namespace Impl
131  {
132 
133  template< class DenseMatrix, class RHS >
134  std::true_type hasDenseMatrixAssigner ( DenseMatrix &, const RHS &, decltype( Dune::DenseMatrixAssigner< DenseMatrix, RHS >::apply( std::declval< DenseMatrix & >(), std::declval< const RHS & >() ) ) * = nullptr );
135 
136  std::false_type hasDenseMatrixAssigner ( ... );
137 
138  } // namespace Impl
139 
140  template< class DenseMatrix, class RHS >
141  struct HasDenseMatrixAssigner
142  : public decltype( Impl::hasDenseMatrixAssigner( std::declval< DenseMatrix & >(), std::declval< const RHS & >() ) )
143  {};
144 
145 #endif // #ifndef DOXYGEN
146 
147 
148 
150  class FMatrixError : public MathError {};
151 
162  template<typename MAT>
163  class DenseMatrix
164  {
165  typedef DenseMatVecTraits<MAT> Traits;
166 
167  // Curiously recurring template pattern
168  MAT & asImp() { return static_cast<MAT&>(*this); }
169  const MAT & asImp() const { return static_cast<const MAT&>(*this); }
170 
171  template <class>
172  friend class DenseMatrix;
173 
174  public:
175  //===== type definitions and constants
176 
178  typedef typename Traits::derived_type derived_type;
179 
181  typedef typename Traits::value_type value_type;
182 
184  typedef typename Traits::value_type field_type;
185 
187  typedef typename Traits::value_type block_type;
188 
190  typedef typename Traits::size_type size_type;
191 
193  typedef typename Traits::row_type row_type;
194 
196  typedef typename Traits::row_reference row_reference;
197 
199  typedef typename Traits::const_row_reference const_row_reference;
200 
202  enum {
204  blocklevel = 1
205  };
206 
207  private:
210  using simd_index_type = Simd::Rebind<std::size_t, value_type>;
211 
212  public:
213  //===== access to components
214 
216  row_reference operator[] ( size_type i )
217  {
218  return asImp().mat_access(i);
219  }
220 
221  const_row_reference operator[] ( size_type i ) const
222  {
223  return asImp().mat_access(i);
224  }
225 
227  size_type size() const
228  {
229  return rows();
230  }
231 
232  //===== iterator interface to rows of the matrix
234  typedef DenseIterator<DenseMatrix,row_type,row_reference> Iterator;
236  typedef Iterator iterator;
238  typedef Iterator RowIterator;
240  typedef typename std::remove_reference<row_reference>::type::Iterator ColIterator;
241 
243  Iterator begin ()
244  {
245  return Iterator(*this,0);
246  }
247 
249  Iterator end ()
250  {
251  return Iterator(*this,rows());
252  }
253 
256  Iterator beforeEnd ()
257  {
258  return Iterator(*this,rows()-1);
259  }
260 
263  Iterator beforeBegin ()
264  {
265  return Iterator(*this,-1);
266  }
267 
269  typedef DenseIterator<const DenseMatrix,const row_type,const_row_reference> ConstIterator;
271  typedef ConstIterator const_iterator;
273  typedef ConstIterator ConstRowIterator;
275  typedef typename std::remove_reference<const_row_reference>::type::ConstIterator ConstColIterator;
276 
278  ConstIterator begin () const
279  {
280  return ConstIterator(*this,0);
281  }
282 
284  ConstIterator end () const
285  {
286  return ConstIterator(*this,rows());
287  }
288 
291  ConstIterator beforeEnd () const
292  {
293  return ConstIterator(*this,rows()-1);
294  }
295 
298  ConstIterator beforeBegin () const
299  {
300  return ConstIterator(*this,-1);
301  }
302 
303  //===== assignment
304 
305  template< class RHS, class = std::enable_if_t< HasDenseMatrixAssigner< MAT, RHS >::value > >
306  derived_type &operator= ( const RHS &rhs )
307  {
308  DenseMatrixAssigner< MAT, RHS >::apply( asImp(), rhs );
309  return asImp();
310  }
311 
312  //===== vector space arithmetic
313 
315  template <class Other>
316  derived_type &operator+= (const DenseMatrix<Other>& x)
317  {
318  DUNE_ASSERT_BOUNDS(rows() == x.rows());
319  for (size_type i=0; i<rows(); i++)
320  (*this)[i] += x[i];
321  return asImp();
322  }
323 
325  derived_type operator- () const
326  {
327  MAT result;
328  typedef typename decltype(result)::size_type size_type;
329 
330  for (size_type i = 0; i < rows(); ++i)
331  for (size_type j = 0; j < cols(); ++j)
332  result[i][j] = - asImp()[i][j];
333 
334  return result;
335  }
336 
338  template <class Other>
339  derived_type &operator-= (const DenseMatrix<Other>& x)
340  {
341  DUNE_ASSERT_BOUNDS(rows() == x.rows());
342  for (size_type i=0; i<rows(); i++)
343  (*this)[i] -= x[i];
344  return asImp();
345  }
346 
348  derived_type &operator*= (const field_type& k)
349  {
350  for (size_type i=0; i<rows(); i++)
351  (*this)[i] *= k;
352  return asImp();
353  }
354 
356  derived_type &operator/= (const field_type& k)
357  {
358  for (size_type i=0; i<rows(); i++)
359  (*this)[i] /= k;
360  return asImp();
361  }
362 
364  template <class Other>
365  derived_type &axpy (const field_type &a, const DenseMatrix<Other> &x )
366  {
367  DUNE_ASSERT_BOUNDS(rows() == x.rows());
368  for( size_type i = 0; i < rows(); ++i )
369  (*this)[ i ].axpy( a, x[ i ] );
370  return asImp();
371  }
372 
374  template <class Other>
375  bool operator== (const DenseMatrix<Other>& x) const
376  {
377  DUNE_ASSERT_BOUNDS(rows() == x.rows());
378  for (size_type i=0; i<rows(); i++)
379  if ((*this)[i]!=x[i])
380  return false;
381  return true;
382  }
384  template <class Other>
385  bool operator!= (const DenseMatrix<Other>& x) const
386  {
387  return !operator==(x);
388  }
389 
390 
391  //===== linear maps
392 
394  template<class X, class Y>
395  void mv (const X& x, Y& y) const
396  {
397  auto&& xx = Impl::asVector(x);
398  auto&& yy = Impl::asVector(y);
399  DUNE_ASSERT_BOUNDS((void*)(&x) != (void*)(&y));
400  DUNE_ASSERT_BOUNDS(xx.N() == M());
401  DUNE_ASSERT_BOUNDS(yy.N() == N());
402 
403  using field_type = typename FieldTraits<Y>::field_type;
404  for (size_type i=0; i<rows(); ++i)
405  {
406  yy[i] = field_type(0);
407  for (size_type j=0; j<cols(); j++)
408  yy[i] += (*this)[i][j] * xx[j];
409  }
410  }
411 
413  template< class X, class Y >
414  void mtv ( const X &x, Y &y ) const
415  {
416  auto&& xx = Impl::asVector(x);
417  auto&& yy = Impl::asVector(y);
418  DUNE_ASSERT_BOUNDS((void*)(&x) != (void*)(&y));
419  DUNE_ASSERT_BOUNDS(xx.N() == N());
420  DUNE_ASSERT_BOUNDS(yy.N() == M());
421 
422  using field_type = typename FieldTraits<Y>::field_type;
423  for(size_type i = 0; i < cols(); ++i)
424  {
425  yy[i] = field_type(0);
426  for(size_type j = 0; j < rows(); ++j)
427  yy[i] += (*this)[j][i] * xx[j];
428  }
429  }
430 
432  template<class X, class Y>
433  void umv (const X& x, Y& y) const
434  {
435  auto&& xx = Impl::asVector(x);
436  auto&& yy = Impl::asVector(y);
437  DUNE_ASSERT_BOUNDS(xx.N() == M());
438  DUNE_ASSERT_BOUNDS(yy.N() == N());
439  for (size_type i=0; i<rows(); ++i)
440  for (size_type j=0; j<cols(); j++)
441  yy[i] += (*this)[i][j] * xx[j];
442  }
443 
445  template<class X, class Y>
446  void umtv (const X& x, Y& y) const
447  {
448  auto&& xx = Impl::asVector(x);
449  auto&& yy = Impl::asVector(y);
450  DUNE_ASSERT_BOUNDS(xx.N() == N());
451  DUNE_ASSERT_BOUNDS(yy.N() == M());
452  for(size_type i = 0; i<rows(); ++i)
453  for (size_type j=0; j<cols(); j++)
454  yy[j] += (*this)[i][j]*xx[i];
455  }
456 
458  template<class X, class Y>
459  void umhv (const X& x, Y& y) const
460  {
461  auto&& xx = Impl::asVector(x);
462  auto&& yy = Impl::asVector(y);
463  DUNE_ASSERT_BOUNDS(xx.N() == N());
464  DUNE_ASSERT_BOUNDS(yy.N() == M());
465  for (size_type i=0; i<rows(); i++)
466  for (size_type j=0; j<cols(); j++)
467  yy[j] += conjugateComplex((*this)[i][j])*xx[i];
468  }
469 
471  template<class X, class Y>
472  void mmv (const X& x, Y& y) const
473  {
474  auto&& xx = Impl::asVector(x);
475  auto&& yy = Impl::asVector(y);
476  DUNE_ASSERT_BOUNDS(xx.N() == M());
477  DUNE_ASSERT_BOUNDS(yy.N() == N());
478  for (size_type i=0; i<rows(); i++)
479  for (size_type j=0; j<cols(); j++)
480  yy[i] -= (*this)[i][j] * xx[j];
481  }
482 
484  template<class X, class Y>
485  void mmtv (const X& x, Y& y) const
486  {
487  auto&& xx = Impl::asVector(x);
488  auto&& yy = Impl::asVector(y);
489  DUNE_ASSERT_BOUNDS(xx.N() == N());
490  DUNE_ASSERT_BOUNDS(yy.N() == M());
491  for (size_type i=0; i<rows(); i++)
492  for (size_type j=0; j<cols(); j++)
493  yy[j] -= (*this)[i][j]*xx[i];
494  }
495 
497  template<class X, class Y>
498  void mmhv (const X& x, Y& y) const
499  {
500  auto&& xx = Impl::asVector(x);
501  auto&& yy = Impl::asVector(y);
502  DUNE_ASSERT_BOUNDS(xx.N() == N());
503  DUNE_ASSERT_BOUNDS(yy.N() == M());
504  for (size_type i=0; i<rows(); i++)
505  for (size_type j=0; j<cols(); j++)
506  yy[j] -= conjugateComplex((*this)[i][j])*xx[i];
507  }
508 
510  template<class X, class Y>
511  void usmv (const typename FieldTraits<Y>::field_type & alpha,
512  const X& x, Y& y) const
513  {
514  auto&& xx = Impl::asVector(x);
515  auto&& yy = Impl::asVector(y);
516  DUNE_ASSERT_BOUNDS(xx.N() == M());
517  DUNE_ASSERT_BOUNDS(yy.N() == N());
518  for (size_type i=0; i<rows(); i++)
519  for (size_type j=0; j<cols(); j++)
520  yy[i] += alpha * (*this)[i][j] * xx[j];
521  }
522 
524  template<class X, class Y>
525  void usmtv (const typename FieldTraits<Y>::field_type & alpha,
526  const X& x, Y& y) const
527  {
528  auto&& xx = Impl::asVector(x);
529  auto&& yy = Impl::asVector(y);
530  DUNE_ASSERT_BOUNDS(xx.N() == N());
531  DUNE_ASSERT_BOUNDS(yy.N() == M());
532  for (size_type i=0; i<rows(); i++)
533  for (size_type j=0; j<cols(); j++)
534  yy[j] += alpha*(*this)[i][j]*xx[i];
535  }
536 
538  template<class X, class Y>
539  void usmhv (const typename FieldTraits<Y>::field_type & alpha,
540  const X& x, Y& y) const
541  {
542  auto&& xx = Impl::asVector(x);
543  auto&& yy = Impl::asVector(y);
544  DUNE_ASSERT_BOUNDS(xx.N() == N());
545  DUNE_ASSERT_BOUNDS(yy.N() == M());
546  for (size_type i=0; i<rows(); i++)
547  for (size_type j=0; j<cols(); j++)
548  yy[j] +=
549  alpha*conjugateComplex((*this)[i][j])*xx[i];
550  }
551 
552  //===== norms
553 
555  typename FieldTraits<value_type>::real_type frobenius_norm () const
556  {
557  typename FieldTraits<value_type>::real_type sum=(0.0);
558  for (size_type i=0; i<rows(); ++i) sum += (*this)[i].two_norm2();
559  return fvmeta::sqrt(sum);
560  }
561 
563  typename FieldTraits<value_type>::real_type frobenius_norm2 () const
564  {
565  typename FieldTraits<value_type>::real_type sum=(0.0);
566  for (size_type i=0; i<rows(); ++i) sum += (*this)[i].two_norm2();
567  return sum;
568  }
569 
571  template <typename vt = value_type,
572  typename std::enable_if<!HasNaN<vt>::value, int>::type = 0>
573  typename FieldTraits<vt>::real_type infinity_norm() const {
574  using real_type = typename FieldTraits<vt>::real_type;
575  using std::max;
576 
577  real_type norm = 0;
578  for (auto const &x : *this) {
579  real_type const a = x.one_norm();
580  norm = max(a, norm);
581  }
582  return norm;
583  }
584 
586  template <typename vt = value_type,
587  typename std::enable_if<!HasNaN<vt>::value, int>::type = 0>
588  typename FieldTraits<vt>::real_type infinity_norm_real() const {
589  using real_type = typename FieldTraits<vt>::real_type;
590  using std::max;
591 
592  real_type norm = 0;
593  for (auto const &x : *this) {
594  real_type const a = x.one_norm_real();
595  norm = max(a, norm);
596  }
597  return norm;
598  }
599 
601  template <typename vt = value_type,
602  typename std::enable_if<HasNaN<vt>::value, int>::type = 0>
603  typename FieldTraits<vt>::real_type infinity_norm() const {
604  using real_type = typename FieldTraits<vt>::real_type;
605  using std::max;
606 
607  real_type norm = 0;
608  real_type isNaN = 1;
609  for (auto const &x : *this) {
610  real_type const a = x.one_norm();
611  norm = max(a, norm);
612  isNaN += a;
613  }
614  return norm * (isNaN / isNaN);
615  }
616 
618  template <typename vt = value_type,
619  typename std::enable_if<HasNaN<vt>::value, int>::type = 0>
620  typename FieldTraits<vt>::real_type infinity_norm_real() const {
621  using real_type = typename FieldTraits<vt>::real_type;
622  using std::max;
623 
624  real_type norm = 0;
625  real_type isNaN = 1;
626  for (auto const &x : *this) {
627  real_type const a = x.one_norm_real();
628  norm = max(a, norm);
629  isNaN += a;
630  }
631  return norm * (isNaN / isNaN);
632  }
633 
634  //===== solve
635 
640  template <class V1, class V2>
641  void solve (V1& x, const V2& b, bool doPivoting = true) const;
642 
647  void invert(bool doPivoting = true);
648 
650  field_type determinant (bool doPivoting = true) const;
651 
653  template<typename M2>
654  MAT& leftmultiply (const DenseMatrix<M2>& M)
655  {
656  DUNE_ASSERT_BOUNDS(M.rows() == M.cols());
657  DUNE_ASSERT_BOUNDS(M.rows() == rows());
658  AutonomousValue<MAT> C(asImp());
659 
660  for (size_type i=0; i<rows(); i++)
661  for (size_type j=0; j<cols(); j++) {
662  (*this)[i][j] = 0;
663  for (size_type k=0; k<rows(); k++)
664  (*this)[i][j] += M[i][k]*C[k][j];
665  }
666 
667  return asImp();
668  }
669 
671  template<typename M2>
672  MAT& rightmultiply (const DenseMatrix<M2>& M)
673  {
674  DUNE_ASSERT_BOUNDS(M.rows() == M.cols());
675  DUNE_ASSERT_BOUNDS(M.cols() == cols());
676  AutonomousValue<MAT> C(asImp());
677 
678  for (size_type i=0; i<rows(); i++)
679  for (size_type j=0; j<cols(); j++) {
680  (*this)[i][j] = 0;
681  for (size_type k=0; k<cols(); k++)
682  (*this)[i][j] += C[i][k]*M[k][j];
683  }
684  return asImp();
685  }
686 
687 #if 0
689  template<int l>
690  DenseMatrix<K,l,cols> leftmultiplyany (const FieldMatrix<K,l,rows>& M) const
691  {
692  FieldMatrix<K,l,cols> C;
693 
694  for (size_type i=0; i<l; i++) {
695  for (size_type j=0; j<cols(); j++) {
696  C[i][j] = 0;
697  for (size_type k=0; k<rows(); k++)
698  C[i][j] += M[i][k]*(*this)[k][j];
699  }
700  }
701  return C;
702  }
703 
705  template<int l>
706  FieldMatrix<K,rows,l> rightmultiplyany (const FieldMatrix<K,cols,l>& M) const
707  {
708  FieldMatrix<K,rows,l> C;
709 
710  for (size_type i=0; i<rows(); i++) {
711  for (size_type j=0; j<l; j++) {
712  C[i][j] = 0;
713  for (size_type k=0; k<cols(); k++)
714  C[i][j] += (*this)[i][k]*M[k][j];
715  }
716  }
717  return C;
718  }
719 #endif
720 
721  //===== sizes
722 
724  size_type N () const
725  {
726  return rows();
727  }
728 
730  size_type M () const
731  {
732  return cols();
733  }
734 
736  size_type rows() const
737  {
738  return asImp().mat_rows();
739  }
740 
742  size_type cols() const
743  {
744  return asImp().mat_cols();
745  }
746 
747  //===== query
748 
750  bool exists (size_type i, size_type j) const
751  {
752  DUNE_UNUSED_PARAMETER(i);
753  DUNE_UNUSED_PARAMETER(j);
754  DUNE_ASSERT_BOUNDS(i >= 0 && i < rows());
755  DUNE_ASSERT_BOUNDS(j >= 0 && j < cols());
756  return true;
757  }
758 
759  protected:
760 
761 #ifndef DOXYGEN
762  struct ElimPivot
763  {
764  ElimPivot(std::vector<simd_index_type> & pivot);
765 
766  void swap(std::size_t i, simd_index_type j);
767 
768  template<typename T>
769  void operator()(const T&, int, int)
770  {}
771 
772  std::vector<simd_index_type> & pivot_;
773  };
774 
775  template<typename V>
776  struct Elim
777  {
778  Elim(V& rhs);
779 
780  void swap(std::size_t i, simd_index_type j);
781 
782  void operator()(const typename V::field_type& factor, int k, int i);
783 
784  V* rhs_;
785  };
786 
787  struct ElimDet
788  {
789  ElimDet(field_type& sign) : sign_(sign)
790  { sign_ = 1; }
791 
792  void swap(std::size_t i, simd_index_type j)
793  {
794  sign_ *=
795  Simd::cond(simd_index_type(i) == j, field_type(1), field_type(-1));
796  }
797 
798  void operator()(const field_type&, int, int)
799  {}
800 
801  field_type& sign_;
802  };
803 #endif // DOXYGEN
804 
806 
844  template<class Func, class Mask>
845  static void luDecomposition(DenseMatrix<MAT>& A, Func func,
846  Mask &nonsingularLanes, bool throwEarly, bool doPivoting);
847  };
848 
849 #ifndef DOXYGEN
850  template<typename MAT>
851  DenseMatrix<MAT>::ElimPivot::ElimPivot(std::vector<simd_index_type> & pivot)
852  : pivot_(pivot)
853  {
854  typedef typename std::vector<size_type>::size_type size_type;
855  for(size_type i=0; i < pivot_.size(); ++i) pivot_[i]=i;
856  }
857 
858  template<typename MAT>
859  void DenseMatrix<MAT>::ElimPivot::swap(std::size_t i, simd_index_type j)
860  {
861  pivot_[i] =
862  Simd::cond(Simd::Scalar<simd_index_type>(i) == j, pivot_[i], j);
863  }
864 
865  template<typename MAT>
866  template<typename V>
867  DenseMatrix<MAT>::Elim<V>::Elim(V& rhs)
868  : rhs_(&rhs)
869  {}
870 
871  template<typename MAT>
872  template<typename V>
873  void DenseMatrix<MAT>::Elim<V>::swap(std::size_t i, simd_index_type j)
874  {
875  using std::swap;
876 
877  // see the comment in luDecomposition()
878  for(std::size_t l = 0; l < Simd::lanes(j); ++l)
879  swap(Simd::lane(l, (*rhs_)[ i ]),
880  Simd::lane(l, (*rhs_)[Simd::lane(l, j)]));
881  }
882 
883  template<typename MAT>
884  template<typename V>
885  void DenseMatrix<MAT>::
886  Elim<V>::operator()(const typename V::field_type& factor, int k, int i)
887  {
888  (*rhs_)[k] -= factor*(*rhs_)[i];
889  }
890 
891  template<typename MAT>
892  template<typename Func, class Mask>
893  inline void DenseMatrix<MAT>::
894  luDecomposition(DenseMatrix<MAT>& A, Func func, Mask &nonsingularLanes,
895  bool throwEarly, bool doPivoting)
896  {
897  using std::max;
898  using std::swap;
899 
900  typedef typename FieldTraits<value_type>::real_type real_type;
901 
902  // LU decomposition of A in A
903  for (size_type i=0; i<A.rows(); i++) // loop over all rows
904  {
905  real_type pivmax = fvmeta::absreal(A[i][i]);
906 
907  if (doPivoting)
908  {
909  // compute maximum of column
910  simd_index_type imax=i;
911  for (size_type k=i+1; k<A.rows(); k++)
912  {
913  auto abs = fvmeta::absreal(A[k][i]);
914  auto mask = abs > pivmax;
915  pivmax = Simd::cond(mask, abs, pivmax);
916  imax = Simd::cond(mask, simd_index_type(k), imax);
917  }
918  // swap rows
919  for (size_type j=0; j<A.rows(); j++)
920  {
921  // This is a swap operation where the second operand is scattered,
922  // and on top of that is also extracted from deep within a
923  // moderately complicated data structure (a DenseMatrix), where we
924  // can't assume much on the memory layout. On intel processors,
925  // the only instruction that might help us here is vgather, but it
926  // is unclear whether that is even faster than a software
927  // implementation, and we would also need vscatter which does not
928  // exist. So break vectorization here and do it manually.
929  for(std::size_t l = 0; l < Simd::lanes(A[i][j]); ++l)
930  swap(Simd::lane(l, A[i][j]),
931  Simd::lane(l, A[Simd::lane(l, imax)][j]));
932  }
933  func.swap(i, imax); // swap the pivot or rhs
934  }
935 
936  // singular ?
937  nonsingularLanes = nonsingularLanes && (pivmax != real_type(0));
938  if (throwEarly) {
939  if(!Simd::allTrue(nonsingularLanes))
940  DUNE_THROW(FMatrixError, "matrix is singular");
941  }
942  else { // !throwEarly
943  if(!Simd::anyTrue(nonsingularLanes))
944  return;
945  }
946 
947  // eliminate
948  for (size_type k=i+1; k<A.rows(); k++)
949  {
950  // in the simd case, A[i][i] may be close to zero in some lanes. Pray
951  // that the result is no worse than a quiet NaN.
952  field_type factor = A[k][i]/A[i][i];
953  A[k][i] = factor;
954  for (size_type j=i+1; j<A.rows(); j++)
955  A[k][j] -= factor*A[i][j];
956  func(factor, k, i);
957  }
958  }
959  }
960 
961  template<typename MAT>
962  template <class V1, class V2>
963  inline void DenseMatrix<MAT>::solve(V1& x, const V2& b, bool doPivoting) const
964  {
965  using real_type = typename FieldTraits<value_type>::real_type;
966  // never mind those ifs, because they get optimized away
967  if (rows()!=cols())
968  DUNE_THROW(FMatrixError, "Can't solve for a " << rows() << "x" << cols() << " matrix!");
969 
970  if (rows()==1) {
971 
972 #ifdef DUNE_FMatrix_WITH_CHECKING
973  if (Simd::anyTrue(fvmeta::absreal((*this)[0][0])
974  < FMatrixPrecision<>::absolute_limit()))
975  DUNE_THROW(FMatrixError,"matrix is singular");
976 #endif
977  x[0] = b[0]/(*this)[0][0];
978 
979  }
980  else if (rows()==2) {
981 
982  field_type detinv = (*this)[0][0]*(*this)[1][1]-(*this)[0][1]*(*this)[1][0];
983 #ifdef DUNE_FMatrix_WITH_CHECKING
984  if (Simd::anyTrue(fvmeta::absreal(detinv)
985  < FMatrixPrecision<>::absolute_limit()))
986  DUNE_THROW(FMatrixError,"matrix is singular");
987 #endif
988  detinv = real_type(1.0)/detinv;
989 
990  x[0] = detinv*((*this)[1][1]*b[0]-(*this)[0][1]*b[1]);
991  x[1] = detinv*((*this)[0][0]*b[1]-(*this)[1][0]*b[0]);
992 
993  }
994  else if (rows()==3) {
995 
996  field_type d = determinant(doPivoting);
997 #ifdef DUNE_FMatrix_WITH_CHECKING
998  if (Simd::anyTrue(fvmeta::absreal(d)
999  < FMatrixPrecision<>::absolute_limit()))
1000  DUNE_THROW(FMatrixError,"matrix is singular");
1001 #endif
1002 
1003  x[0] = (b[0]*(*this)[1][1]*(*this)[2][2] - b[0]*(*this)[2][1]*(*this)[1][2]
1004  - b[1] *(*this)[0][1]*(*this)[2][2] + b[1]*(*this)[2][1]*(*this)[0][2]
1005  + b[2] *(*this)[0][1]*(*this)[1][2] - b[2]*(*this)[1][1]*(*this)[0][2]) / d;
1006 
1007  x[1] = ((*this)[0][0]*b[1]*(*this)[2][2] - (*this)[0][0]*b[2]*(*this)[1][2]
1008  - (*this)[1][0] *b[0]*(*this)[2][2] + (*this)[1][0]*b[2]*(*this)[0][2]
1009  + (*this)[2][0] *b[0]*(*this)[1][2] - (*this)[2][0]*b[1]*(*this)[0][2]) / d;
1010 
1011  x[2] = ((*this)[0][0]*(*this)[1][1]*b[2] - (*this)[0][0]*(*this)[2][1]*b[1]
1012  - (*this)[1][0] *(*this)[0][1]*b[2] + (*this)[1][0]*(*this)[2][1]*b[0]
1013  + (*this)[2][0] *(*this)[0][1]*b[1] - (*this)[2][0]*(*this)[1][1]*b[0]) / d;
1014 
1015  }
1016  else {
1017 
1018  V1& rhs = x; // use x to store rhs
1019  rhs = b; // copy data
1020  Elim<V1> elim(rhs);
1021  AutonomousValue<MAT> A(asImp());
1022  Simd::Mask<typename FieldTraits<value_type>::real_type>
1023  nonsingularLanes(true);
1024 
1025  AutonomousValue<MAT>::luDecomposition(A, elim, nonsingularLanes, true, doPivoting);
1026 
1027  // backsolve
1028  for(int i=rows()-1; i>=0; i--) {
1029  for (size_type j=i+1; j<rows(); j++)
1030  rhs[i] -= A[i][j]*x[j];
1031  x[i] = rhs[i]/A[i][i];
1032  }
1033  }
1034  }
1035 
1036  template<typename MAT>
1037  inline void DenseMatrix<MAT>::invert(bool doPivoting)
1038  {
1039  using real_type = typename FieldTraits<MAT>::real_type;
1040  using std::swap;
1041 
1042  // never mind those ifs, because they get optimized away
1043  if (rows()!=cols())
1044  DUNE_THROW(FMatrixError, "Can't invert a " << rows() << "x" << cols() << " matrix!");
1045 
1046  if (rows()==1) {
1047 
1048 #ifdef DUNE_FMatrix_WITH_CHECKING
1049  if (Simd::anyTrue(fvmeta::absreal((*this)[0][0])
1050  < FMatrixPrecision<>::absolute_limit()))
1051  DUNE_THROW(FMatrixError,"matrix is singular");
1052 #endif
1053  (*this)[0][0] = real_type( 1 ) / (*this)[0][0];
1054 
1055  }
1056  else if (rows()==2) {
1057 
1058  field_type detinv = (*this)[0][0]*(*this)[1][1]-(*this)[0][1]*(*this)[1][0];
1059 #ifdef DUNE_FMatrix_WITH_CHECKING
1060  if (Simd::anyTrue(fvmeta::absreal(detinv)
1061  < FMatrixPrecision<>::absolute_limit()))
1062  DUNE_THROW(FMatrixError,"matrix is singular");
1063 #endif
1064  detinv = real_type( 1 ) / detinv;
1065 
1066  field_type temp=(*this)[0][0];
1067  (*this)[0][0] = (*this)[1][1]*detinv;
1068  (*this)[0][1] = -(*this)[0][1]*detinv;
1069  (*this)[1][0] = -(*this)[1][0]*detinv;
1070  (*this)[1][1] = temp*detinv;
1071 
1072  }
1073  else if (rows()==3)
1074  {
1075  using K = field_type;
1076  // code generated by maple
1077  K t4 = (*this)[0][0] * (*this)[1][1];
1078  K t6 = (*this)[0][0] * (*this)[1][2];
1079  K t8 = (*this)[0][1] * (*this)[1][0];
1080  K t10 = (*this)[0][2] * (*this)[1][0];
1081  K t12 = (*this)[0][1] * (*this)[2][0];
1082  K t14 = (*this)[0][2] * (*this)[2][0];
1083 
1084  K det = (t4*(*this)[2][2]-t6*(*this)[2][1]-t8*(*this)[2][2]+
1085  t10*(*this)[2][1]+t12*(*this)[1][2]-t14*(*this)[1][1]);
1086  K t17 = K(1.0)/det;
1087 
1088  K matrix01 = (*this)[0][1];
1089  K matrix00 = (*this)[0][0];
1090  K matrix10 = (*this)[1][0];
1091  K matrix11 = (*this)[1][1];
1092 
1093  (*this)[0][0] = ((*this)[1][1] * (*this)[2][2] - (*this)[1][2] * (*this)[2][1])*t17;
1094  (*this)[0][1] = -((*this)[0][1] * (*this)[2][2] - (*this)[0][2] * (*this)[2][1])*t17;
1095  (*this)[0][2] = (matrix01 * (*this)[1][2] - (*this)[0][2] * (*this)[1][1])*t17;
1096  (*this)[1][0] = -((*this)[1][0] * (*this)[2][2] - (*this)[1][2] * (*this)[2][0])*t17;
1097  (*this)[1][1] = (matrix00 * (*this)[2][2] - t14) * t17;
1098  (*this)[1][2] = -(t6-t10) * t17;
1099  (*this)[2][0] = (matrix10 * (*this)[2][1] - matrix11 * (*this)[2][0]) * t17;
1100  (*this)[2][1] = -(matrix00 * (*this)[2][1] - t12) * t17;
1101  (*this)[2][2] = (t4-t8) * t17;
1102  }
1103  else {
1104  using std::swap;
1105 
1106  AutonomousValue<MAT> A(asImp());
1107  std::vector<simd_index_type> pivot(rows());
1108  Simd::Mask<typename FieldTraits<value_type>::real_type>
1109  nonsingularLanes(true);
1110  AutonomousValue<MAT>::luDecomposition(A, ElimPivot(pivot), nonsingularLanes, true, doPivoting);
1111  auto& L=A;
1112  auto& U=A;
1113 
1114  // initialize inverse
1115  *this=field_type();
1116 
1117  for(size_type i=0; i<rows(); ++i)
1118  (*this)[i][i]=1;
1119 
1120  // L Y = I; multiple right hand sides
1121  for (size_type i=0; i<rows(); i++)
1122  for (size_type j=0; j<i; j++)
1123  for (size_type k=0; k<rows(); k++)
1124  (*this)[i][k] -= L[i][j]*(*this)[j][k];
1125 
1126  // U A^{-1} = Y
1127  for (size_type i=rows(); i>0;) {
1128  --i;
1129  for (size_type k=0; k<rows(); k++) {
1130  for (size_type j=i+1; j<rows(); j++)
1131  (*this)[i][k] -= U[i][j]*(*this)[j][k];
1132  (*this)[i][k] /= U[i][i];
1133  }
1134  }
1135 
1136  for(size_type i=rows(); i>0; ) {
1137  --i;
1138  for(std::size_t l = 0; l < Simd::lanes((*this)[0][0]); ++l)
1139  {
1140  std::size_t pi = Simd::lane(l, pivot[i]);
1141  if(i!=pi)
1142  for(size_type j=0; j<rows(); ++j)
1143  swap(Simd::lane(l, (*this)[j][pi]),
1144  Simd::lane(l, (*this)[j][ i]));
1145  }
1146  }
1147  }
1148  }
1149 
1150  // implementation of the determinant
1151  template<typename MAT>
1152  inline typename DenseMatrix<MAT>::field_type
1153  DenseMatrix<MAT>::determinant(bool doPivoting) const
1154  {
1155  // never mind those ifs, because they get optimized away
1156  if (rows()!=cols())
1157  DUNE_THROW(FMatrixError, "There is no determinant for a " << rows() << "x" << cols() << " matrix!");
1158 
1159  if (rows()==1)
1160  return (*this)[0][0];
1161 
1162  if (rows()==2)
1163  return (*this)[0][0]*(*this)[1][1] - (*this)[0][1]*(*this)[1][0];
1164 
1165  if (rows()==3) {
1166  // code generated by maple
1167  field_type t4 = (*this)[0][0] * (*this)[1][1];
1168  field_type t6 = (*this)[0][0] * (*this)[1][2];
1169  field_type t8 = (*this)[0][1] * (*this)[1][0];
1170  field_type t10 = (*this)[0][2] * (*this)[1][0];
1171  field_type t12 = (*this)[0][1] * (*this)[2][0];
1172  field_type t14 = (*this)[0][2] * (*this)[2][0];
1173 
1174  return (t4*(*this)[2][2]-t6*(*this)[2][1]-t8*(*this)[2][2]+
1175  t10*(*this)[2][1]+t12*(*this)[1][2]-t14*(*this)[1][1]);
1176 
1177  }
1178 
1179  AutonomousValue<MAT> A(asImp());
1180  field_type det;
1181  Simd::Mask<typename FieldTraits<value_type>::real_type>
1182  nonsingularLanes(true);
1183 
1184  AutonomousValue<MAT>::luDecomposition(A, ElimDet(det), nonsingularLanes, false, doPivoting);
1185  det = Simd::cond(nonsingularLanes, det, field_type(0));
1186 
1187  for (size_type i = 0; i < rows(); ++i)
1188  det *= A[i][i];
1189  return det;
1190  }
1191 
1192 #endif // DOXYGEN
1193 
1194  namespace DenseMatrixHelp {
1195 
1197  template <typename MAT, typename V1, typename V2>
1198  static inline void multAssign(const DenseMatrix<MAT> &matrix, const DenseVector<V1> & x, DenseVector<V2> & ret)
1199  {
1200  DUNE_ASSERT_BOUNDS(x.size() == matrix.cols());
1201  DUNE_ASSERT_BOUNDS(ret.size() == matrix.rows());
1202  typedef typename DenseMatrix<MAT>::size_type size_type;
1203 
1204  for(size_type i=0; i<matrix.rows(); ++i)
1205  {
1206  ret[i] = 0.0;
1207  for(size_type j=0; j<matrix.cols(); ++j)
1208  {
1209  ret[i] += matrix[i][j]*x[j];
1210  }
1211  }
1212  }
1213 
1214 #if 0
1216  template <typename K, int rows, int cols>
1217  static inline void multAssignTransposed( const FieldMatrix<K,rows,cols> &matrix, const FieldVector<K,rows> & x, FieldVector<K,cols> & ret)
1218  {
1219  typedef typename FieldMatrix<K,rows,cols>::size_type size_type;
1220 
1221  for(size_type i=0; i<cols(); ++i)
1222  {
1223  ret[i] = 0.0;
1224  for(size_type j=0; j<rows(); ++j)
1225  ret[i] += matrix[j][i]*x[j];
1226  }
1227  }
1228 
1230  template <typename K, int rows, int cols>
1231  static inline FieldVector<K,rows> mult(const FieldMatrix<K,rows,cols> &matrix, const FieldVector<K,cols> & x)
1232  {
1233  FieldVector<K,rows> ret;
1234  multAssign(matrix,x,ret);
1235  return ret;
1236  }
1237 
1239  template <typename K, int rows, int cols>
1240  static inline FieldVector<K,cols> multTransposed(const FieldMatrix<K,rows,cols> &matrix, const FieldVector<K,rows> & x)
1241  {
1242  FieldVector<K,cols> ret;
1243  multAssignTransposed( matrix, x, ret );
1244  return ret;
1245  }
1246 #endif
1247 
1248  } // end namespace DenseMatrixHelp
1249 
1251  template<typename MAT>
1252  std::ostream& operator<< (std::ostream& s, const DenseMatrix<MAT>& a)
1253  {
1254  for (typename DenseMatrix<MAT>::size_type i=0; i<a.rows(); i++)
1255  s << a[i] << std::endl;
1256  return s;
1257  }
1258 
1261 } // end namespace Dune
1262 
1263 #endif
Macro for wrapping boundary checks.
A free function to provide the demangled class name of a given object or type as a string.
Definition of the DUNE_DEPRECATED macro for the case that config.h is not available.
A few common exception classes.
Implements a vector constructed from a given type representing a field and a compile-time given size.
#define DUNE_DEPRECATED_MSG(text)
Mark some entity as deprecated.
Definition: deprecated.hh:169
Some useful basic math stuff.
Dune namespace.
Definition: alignedallocator.hh:14
Various precision settings for calculations with FieldMatrix and FieldVector.
Implements a scalar vector view wrapper around an existing scalar.
Include file for users of the SIMD abstraction layer.
Traits for type conversions and type information.
Definition of the DUNE_UNUSED macro for the case that config.h is not available.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (May 9, 22:29, 2024)