dune-common 2.11
Loading...
Searching...
No Matches
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// SPDX-FileCopyrightInfo: Copyright © DUNE Project contributors, see file LICENSE.md in module root
4// SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
5#ifndef DUNE_DENSEMATRIX_HH
6#define DUNE_DENSEMATRIX_HH
7
8#include <cmath>
9#include <cstddef>
10#include <iostream>
11#include <iterator>
12#include <type_traits>
13#include <utility>
14#include <vector>
15
20#include <dune/common/math.hh>
26
27namespace Dune
28{
29
30 template<typename M> class DenseMatrix;
31
32 template<typename M>
38
39 template<class K, int N, int M> class FieldMatrix;
40 template<class K, int N> class FieldVector;
41
46
52
53
54
60 template< class DenseMatrix, class RHS >
62
63#ifndef DOXYGEN
64 namespace Impl
65 {
66
67 template< class DenseMatrix, class RHS >
69 {};
70
71 template< class DenseMatrix, class RHS >
74 {
75 public:
76 constexpr static void apply ( DenseMatrix &denseMatrix, const RHS &rhs )
77 {
78 typedef typename DenseMatrix::field_type field_type;
79 std::fill( denseMatrix.begin(), denseMatrix.end(), static_cast< field_type >( rhs ) );
80 }
81 };
82
83 template< class DenseMatrix, class RHS >
84 requires Std::indirectly_copyable<
85 decltype(std::begin(*std::declval<typename RHS::const_iterator>())),
86 decltype(std::begin(*std::declval<typename DenseMatrix::iterator>()))>
87 class DenseMatrixAssigner<DenseMatrix, RHS>
88 {
89 public:
90 constexpr static void apply ( DenseMatrix &denseMatrix, const RHS &rhs )
91 {
92 DUNE_ASSERT_BOUNDS(rhs.N() == denseMatrix.N());
93 DUNE_ASSERT_BOUNDS(rhs.M() == denseMatrix.M());
94 typename DenseMatrix::iterator tIt = std::begin(denseMatrix);
95 typename RHS::const_iterator sIt = std::begin(rhs);
96 for(; sIt != std::end(rhs); ++tIt, ++sIt)
97 std::copy(std::begin(*sIt), std::end(*sIt), std::begin(*tIt));
98 }
99 };
100
101 } // namespace Impl
102
103
104
105 template< class DenseMatrix, class RHS >
107 : public Impl::DenseMatrixAssigner< DenseMatrix, RHS >
108 {};
109
110
111 namespace Impl
112 {
113
114 template< class DenseMatrix, class RHS >
115 std::true_type hasDenseMatrixAssigner ( DenseMatrix &, const RHS &, decltype( Dune::DenseMatrixAssigner< DenseMatrix, RHS >::apply( std::declval< DenseMatrix & >(), std::declval< const RHS & >() ) ) * = nullptr );
116
117 std::false_type hasDenseMatrixAssigner ( ... );
118
119 } // namespace Impl
120
121 template< class DenseMatrix, class RHS >
122 struct HasDenseMatrixAssigner
123 : public decltype( Impl::hasDenseMatrixAssigner( std::declval< DenseMatrix & >(), std::declval< const RHS & >() ) )
124 {};
125
126#endif // #ifndef DOXYGEN
127
128
129
131 class FMatrixError : public MathError {};
132
143 template<typename MAT>
145 {
146 typedef DenseMatVecTraits<MAT> Traits;
147
148 // Curiously recurring template pattern
149 constexpr MAT & asImp() { return static_cast<MAT&>(*this); }
150 constexpr const MAT & asImp() const { return static_cast<const MAT&>(*this); }
151
152 template <class>
153 friend class DenseMatrix;
154
155 public:
156 //===== type definitions and constants
157
159 typedef typename Traits::derived_type derived_type;
160
162 typedef typename Traits::value_type value_type;
163
165 typedef typename Traits::value_type field_type;
166
168 typedef typename Traits::value_type block_type;
169
171 typedef typename Traits::size_type size_type;
172
174 typedef typename Traits::row_type row_type;
175
177 typedef typename Traits::row_reference row_reference;
178
180 typedef typename Traits::const_row_reference const_row_reference;
181
183 constexpr static int blocklevel = 1;
184
185 private:
188 using simd_index_type = Simd::Rebind<std::size_t, value_type>;
189
190 public:
191 //===== access to components
192
195 {
196 return asImp().mat_access(i);
197 }
198
200 {
201 return asImp().mat_access(i);
202 }
203
205 constexpr size_type size() const
206 {
207 return rows();
208 }
209
210 //===== iterator interface to rows of the matrix
218 typedef typename std::remove_reference<row_reference>::type::Iterator ColIterator;
219
221 constexpr Iterator begin ()
222 {
223 return Iterator(*this,0);
224 }
225
227 constexpr Iterator end ()
228 {
229 return Iterator(*this,rows());
230 }
231
234 constexpr Iterator beforeEnd ()
235 {
236 return Iterator(*this,rows()-1);
237 }
238
242 {
243 return Iterator(*this,-1);
244 }
245
253 typedef typename std::remove_reference<const_row_reference>::type::ConstIterator ConstColIterator;
254
256 constexpr ConstIterator begin () const
257 {
258 return ConstIterator(*this,0);
259 }
260
262 constexpr ConstIterator end () const
263 {
264 return ConstIterator(*this,rows());
265 }
266
269 constexpr ConstIterator beforeEnd () const
270 {
271 return ConstIterator(*this,rows()-1);
272 }
273
276 constexpr ConstIterator beforeBegin () const
277 {
278 return ConstIterator(*this,-1);
279 }
280
281 //===== assignment
282
283 template< class RHS, class = std::enable_if_t< HasDenseMatrixAssigner< MAT, RHS >::value > >
284 constexpr derived_type &operator= ( const RHS &rhs )
285 {
287 return asImp();
288 }
289
290 //===== vector space arithmetic
291
293 template <class Other>
295 {
296 DUNE_ASSERT_BOUNDS(rows() == x.rows());
297 for (size_type i=0; i<rows(); i++)
298 (*this)[i] += x[i];
299 return asImp();
300 }
301
303 constexpr derived_type operator- () const
304 {
305 MAT result;
306 using idx_type = typename decltype(result)::size_type;
307
308 for (idx_type i = 0; i < rows(); ++i)
309 for (idx_type j = 0; j < cols(); ++j)
310 result[i][j] = - asImp()[i][j];
311
312 return result;
313 }
314
316 template <class Other>
318 {
319 DUNE_ASSERT_BOUNDS(rows() == x.rows());
320 for (size_type i=0; i<rows(); i++)
321 (*this)[i] -= x[i];
322 return asImp();
323 }
324
327 {
328 for (size_type i=0; i<rows(); i++)
329 (*this)[i] *= k;
330 return asImp();
331 }
332
335 {
336 for (size_type i=0; i<rows(); i++)
337 (*this)[i] /= k;
338 return asImp();
339 }
340
342 template <class Other>
343 constexpr derived_type &axpy (const field_type &a, const DenseMatrix<Other> &x )
344 {
345 DUNE_ASSERT_BOUNDS(rows() == x.rows());
346 for( size_type i = 0; i < rows(); ++i )
347 (*this)[ i ].axpy( a, x[ i ] );
348 return asImp();
349 }
350
352 template <class Other>
353 constexpr bool operator== (const DenseMatrix<Other>& x) const
354 {
355 DUNE_ASSERT_BOUNDS(rows() == x.rows());
356 for (size_type i=0; i<rows(); i++)
357 if ((*this)[i]!=x[i])
358 return false;
359 return true;
360 }
361
362 template <class Other>
363 constexpr bool operator!= (const DenseMatrix<Other>& x) const
364 {
365 return !operator==(x);
366 }
367
368
369 //===== linear maps
370
372 template<class X, class Y>
373 constexpr void mv (const X& x, Y& y) const
374 {
375 auto&& xx = Impl::asVector(x);
376 auto&& yy = Impl::asVector(y);
377 DUNE_ASSERT_BOUNDS((void*)(&x) != (void*)(&y));
378 DUNE_ASSERT_BOUNDS(xx.N() == M());
379 DUNE_ASSERT_BOUNDS(yy.N() == N());
380
381 using y_field_type = typename FieldTraits<Y>::field_type;
382 for (size_type i=0; i<rows(); ++i)
383 {
384 yy[i] = y_field_type(0);
385 for (size_type j=0; j<cols(); j++)
386 yy[i] += (*this)[i][j] * xx[j];
387 }
388 }
389
391 template< class X, class Y >
392 constexpr void mtv ( const X &x, Y &y ) const
393 {
394 auto&& xx = Impl::asVector(x);
395 auto&& yy = Impl::asVector(y);
396 DUNE_ASSERT_BOUNDS((void*)(&x) != (void*)(&y));
397 DUNE_ASSERT_BOUNDS(xx.N() == N());
398 DUNE_ASSERT_BOUNDS(yy.N() == M());
399
400 using y_field_type = typename FieldTraits<Y>::field_type;
401 for(size_type i = 0; i < cols(); ++i)
402 {
403 yy[i] = y_field_type(0);
404 for(size_type j = 0; j < rows(); ++j)
405 yy[i] += (*this)[j][i] * xx[j];
406 }
407 }
408
410 template<class X, class Y>
411 constexpr void umv (const X& x, Y& y) const
412 {
413 auto&& xx = Impl::asVector(x);
414 auto&& yy = Impl::asVector(y);
415 DUNE_ASSERT_BOUNDS(xx.N() == M());
416 DUNE_ASSERT_BOUNDS(yy.N() == N());
417 for (size_type i=0; i<rows(); ++i)
418 for (size_type j=0; j<cols(); j++)
419 yy[i] += (*this)[i][j] * xx[j];
420 }
421
423 template<class X, class Y>
424 constexpr void umtv (const X& x, Y& y) const
425 {
426 auto&& xx = Impl::asVector(x);
427 auto&& yy = Impl::asVector(y);
428 DUNE_ASSERT_BOUNDS(xx.N() == N());
429 DUNE_ASSERT_BOUNDS(yy.N() == M());
430 for(size_type i = 0; i<rows(); ++i)
431 for (size_type j=0; j<cols(); j++)
432 yy[j] += (*this)[i][j]*xx[i];
433 }
434
436 template<class X, class Y>
437 constexpr void umhv (const X& x, Y& y) const
438 {
439 auto&& xx = Impl::asVector(x);
440 auto&& yy = Impl::asVector(y);
441 DUNE_ASSERT_BOUNDS(xx.N() == N());
442 DUNE_ASSERT_BOUNDS(yy.N() == M());
443 for (size_type i=0; i<rows(); i++)
444 for (size_type j=0; j<cols(); j++)
445 yy[j] += conjugateComplex((*this)[i][j])*xx[i];
446 }
447
449 template<class X, class Y>
450 constexpr void mmv (const X& x, Y& y) const
451 {
452 auto&& xx = Impl::asVector(x);
453 auto&& yy = Impl::asVector(y);
454 DUNE_ASSERT_BOUNDS(xx.N() == M());
455 DUNE_ASSERT_BOUNDS(yy.N() == N());
456 for (size_type i=0; i<rows(); i++)
457 for (size_type j=0; j<cols(); j++)
458 yy[i] -= (*this)[i][j] * xx[j];
459 }
460
462 template<class X, class Y>
463 constexpr void mmtv (const X& x, Y& y) const
464 {
465 auto&& xx = Impl::asVector(x);
466 auto&& yy = Impl::asVector(y);
467 DUNE_ASSERT_BOUNDS(xx.N() == N());
468 DUNE_ASSERT_BOUNDS(yy.N() == M());
469 for (size_type i=0; i<rows(); i++)
470 for (size_type j=0; j<cols(); j++)
471 yy[j] -= (*this)[i][j]*xx[i];
472 }
473
475 template<class X, class Y>
476 constexpr void mmhv (const X& x, Y& y) const
477 {
478 auto&& xx = Impl::asVector(x);
479 auto&& yy = Impl::asVector(y);
480 DUNE_ASSERT_BOUNDS(xx.N() == N());
481 DUNE_ASSERT_BOUNDS(yy.N() == M());
482 for (size_type i=0; i<rows(); i++)
483 for (size_type j=0; j<cols(); j++)
484 yy[j] -= conjugateComplex((*this)[i][j])*xx[i];
485 }
486
488 template<class X, class Y>
489 constexpr void usmv (const typename FieldTraits<Y>::field_type & alpha,
490 const X& x, Y& y) const
491 {
492 auto&& xx = Impl::asVector(x);
493 auto&& yy = Impl::asVector(y);
494 DUNE_ASSERT_BOUNDS(xx.N() == M());
495 DUNE_ASSERT_BOUNDS(yy.N() == N());
496 for (size_type i=0; i<rows(); i++)
497 for (size_type j=0; j<cols(); j++)
498 yy[i] += alpha * (*this)[i][j] * xx[j];
499 }
500
502 template<class X, class Y>
503 constexpr void usmtv (const typename FieldTraits<Y>::field_type & alpha,
504 const X& x, Y& y) const
505 {
506 auto&& xx = Impl::asVector(x);
507 auto&& yy = Impl::asVector(y);
508 DUNE_ASSERT_BOUNDS(xx.N() == N());
509 DUNE_ASSERT_BOUNDS(yy.N() == M());
510 for (size_type i=0; i<rows(); i++)
511 for (size_type j=0; j<cols(); j++)
512 yy[j] += alpha*(*this)[i][j]*xx[i];
513 }
514
516 template<class X, class Y>
517 constexpr void usmhv (const typename FieldTraits<Y>::field_type & alpha,
518 const X& x, Y& y) const
519 {
520 auto&& xx = Impl::asVector(x);
521 auto&& yy = Impl::asVector(y);
522 DUNE_ASSERT_BOUNDS(xx.N() == N());
523 DUNE_ASSERT_BOUNDS(yy.N() == M());
524 for (size_type i=0; i<rows(); i++)
525 for (size_type j=0; j<cols(); j++)
526 yy[j] +=
527 alpha*conjugateComplex((*this)[i][j])*xx[i];
528 }
529
530 //===== norms
531
534 {
535 typename FieldTraits<value_type>::real_type sum=(0.0);
536 for (size_type i=0; i<rows(); ++i) sum += (*this)[i].two_norm2();
537 return fvmeta::sqrt(sum);
538 }
539
542 {
543 typename FieldTraits<value_type>::real_type sum=(0.0);
544 for (size_type i=0; i<rows(); ++i) sum += (*this)[i].two_norm2();
545 return sum;
546 }
547
549 template <typename vt = value_type,
550 typename std::enable_if<!HasNaN<vt>::value, int>::type = 0>
551 constexpr typename FieldTraits<vt>::real_type infinity_norm() const {
552 using real_type = typename FieldTraits<vt>::real_type;
553 using std::max;
554
555 real_type norm = 0;
556 for (auto const &x : *this) {
557 real_type const a = x.one_norm();
558 norm = max(a, norm);
559 }
560 return norm;
561 }
562
564 template <typename vt = value_type,
565 typename std::enable_if<!HasNaN<vt>::value, int>::type = 0>
567 using real_type = typename FieldTraits<vt>::real_type;
568 using std::max;
569
570 real_type norm = 0;
571 for (auto const &x : *this) {
572 real_type const a = x.one_norm_real();
573 norm = max(a, norm);
574 }
575 return norm;
576 }
577
579 template <typename vt = value_type,
580 typename std::enable_if<HasNaN<vt>::value, int>::type = 0>
581 constexpr typename FieldTraits<vt>::real_type infinity_norm() const {
582 using real_type = typename FieldTraits<vt>::real_type;
583 using std::max;
584
585 real_type norm = 0;
586 real_type is_nan = 1;
587 for (auto const &x : *this) {
588 real_type const a = x.one_norm();
589 norm = max(a, norm);
590 is_nan += a;
591 }
592 return norm * (is_nan / is_nan);
593 }
594
596 template <typename vt = value_type,
597 typename std::enable_if<HasNaN<vt>::value, int>::type = 0>
599 using real_type = typename FieldTraits<vt>::real_type;
600 using std::max;
601
602 real_type norm = 0;
603 real_type is_nan = 1;
604 for (auto const &x : *this) {
605 real_type const a = x.one_norm_real();
606 norm = max(a, norm);
607 is_nan += a;
608 }
609 return norm * (is_nan / is_nan);
610 }
611
612 //===== solve
613
618 template <class V1, class V2>
619 void solve (V1& x, const V2& b, bool doPivoting = true) const;
620
625 void invert(bool doPivoting = true);
626
628 field_type determinant (bool doPivoting = true) const;
629
631 template<typename M2>
633 {
634 DUNE_ASSERT_BOUNDS(M.rows() == M.cols());
635 DUNE_ASSERT_BOUNDS(M.rows() == rows());
636 AutonomousValue<MAT> C(asImp());
637
638 for (size_type i=0; i<rows(); i++)
639 for (size_type j=0; j<cols(); j++) {
640 (*this)[i][j] = 0;
641 for (size_type k=0; k<rows(); k++)
642 (*this)[i][j] += M[i][k]*C[k][j];
643 }
644
645 return asImp();
646 }
647
649 template<typename M2>
651 {
652 DUNE_ASSERT_BOUNDS(M.rows() == M.cols());
653 DUNE_ASSERT_BOUNDS(M.cols() == cols());
654 AutonomousValue<MAT> C(asImp());
655
656 for (size_type i=0; i<rows(); i++)
657 for (size_type j=0; j<cols(); j++) {
658 (*this)[i][j] = 0;
659 for (size_type k=0; k<cols(); k++)
660 (*this)[i][j] += C[i][k]*M[k][j];
661 }
662 return asImp();
663 }
664
665#if 0
667 template<int l>
668 DenseMatrix<K,l,cols> leftmultiplyany (const FieldMatrix<K,l,rows>& M) const
669 {
671
672 for (size_type i=0; i<l; i++) {
673 for (size_type j=0; j<cols(); j++) {
674 C[i][j] = 0;
675 for (size_type k=0; k<rows(); k++)
676 C[i][j] += M[i][k]*(*this)[k][j];
677 }
678 }
679 return C;
680 }
681
683 template<int l>
684 FieldMatrix<K,rows,l> rightmultiplyany (const FieldMatrix<K,cols,l>& M) const
685 {
686 FieldMatrix<K,rows,l> C;
687
688 for (size_type i=0; i<rows(); i++) {
689 for (size_type j=0; j<l; j++) {
690 C[i][j] = 0;
691 for (size_type k=0; k<cols(); k++)
692 C[i][j] += (*this)[i][k]*M[k][j];
693 }
694 }
695 return C;
696 }
697#endif
698
699 //===== sizes
700
702 constexpr size_type N () const
703 {
704 return rows();
705 }
706
708 constexpr size_type M () const
709 {
710 return cols();
711 }
712
714 constexpr size_type rows() const
715 {
716 return asImp().mat_rows();
717 }
718
720 constexpr size_type cols() const
721 {
722 return asImp().mat_cols();
723 }
724
725 //===== query
726
728 constexpr bool exists ([[maybe_unused]] size_type i, [[maybe_unused]] size_type j) const
729 {
730 DUNE_ASSERT_BOUNDS(i >= 0 && i < rows());
731 DUNE_ASSERT_BOUNDS(j >= 0 && j < cols());
732 return true;
733 }
734
735 protected:
736
737#ifndef DOXYGEN
738 struct ElimPivot
739 {
740 ElimPivot(std::vector<simd_index_type> & pivot);
741
742 void swap(std::size_t i, simd_index_type j);
743
744 template<typename T>
745 void operator()(const T&, int, int)
746 {}
747
748 std::vector<simd_index_type> & pivot_;
749 };
750
751 template<typename V>
752 struct Elim
753 {
754 Elim(V& rhs);
755
756 void swap(std::size_t i, simd_index_type j);
757
758 void operator()(const typename V::field_type& factor, int k, int i);
759
760 V* rhs_;
761 };
762
763 struct ElimDet
764 {
765 ElimDet(field_type& sign) : sign_(sign)
766 { sign_ = 1; }
767
768 void swap(std::size_t i, simd_index_type j)
769 {
770 sign_ *=
771 Simd::cond(simd_index_type(i) == j, field_type(1), field_type(-1));
772 }
773
774 void operator()(const field_type&, int, int)
775 {}
776
777 field_type& sign_;
778 };
779#endif // DOXYGEN
780
782
820 template<class Func, class Mask>
821 static void luDecomposition(DenseMatrix<MAT>& A, Func func,
822 Mask &nonsingularLanes, bool throwEarly, bool doPivoting);
823 };
824
825#ifndef DOXYGEN
826 template<typename MAT>
827 DenseMatrix<MAT>::ElimPivot::ElimPivot(std::vector<simd_index_type> & pivot)
828 : pivot_(pivot)
829 {
830 for(typename std::vector<size_type>::size_type i=0; i < pivot_.size(); ++i)
831 pivot_[i]=i;
832 }
833
834 template<typename MAT>
835 void DenseMatrix<MAT>::ElimPivot::swap(std::size_t i, simd_index_type j)
836 {
837 pivot_[i] =
838 Simd::cond(Simd::Scalar<simd_index_type>(i) == j, pivot_[i], j);
839 }
840
841 template<typename MAT>
842 template<typename V>
844 : rhs_(&rhs)
845 {}
846
847 template<typename MAT>
848 template<typename V>
849 void DenseMatrix<MAT>::Elim<V>::swap(std::size_t i, simd_index_type j)
850 {
851 using std::swap;
852
853 // see the comment in luDecomposition()
854 for(std::size_t l = 0; l < Simd::lanes(j); ++l)
855 swap(Simd::lane(l, (*rhs_)[ i ]),
856 Simd::lane(l, (*rhs_)[Simd::lane(l, j)]));
857 }
858
859 template<typename MAT>
860 template<typename V>
862 Elim<V>::operator()(const typename V::field_type& factor, int k, int i)
863 {
864 (*rhs_)[k] -= factor*(*rhs_)[i];
865 }
866
867 template<typename MAT>
868 template<typename Func, class Mask>
869 inline void DenseMatrix<MAT>::
870 luDecomposition(DenseMatrix<MAT>& A, Func func, Mask &nonsingularLanes,
871 bool throwEarly, bool doPivoting)
872 {
873 using std::max;
874 using std::swap;
875
877
878 // LU decomposition of A in A
879 for (size_type i=0; i<A.rows(); i++) // loop over all rows
880 {
881 real_type pivmax = fvmeta::absreal(A[i][i]);
882
883 if (doPivoting)
884 {
885 // compute maximum of column
886 simd_index_type imax=i;
887 for (size_type k=i+1; k<A.rows(); k++)
888 {
889 auto abs = fvmeta::absreal(A[k][i]);
890 auto mask = abs > pivmax;
891 pivmax = Simd::cond(mask, abs, pivmax);
892 imax = Simd::cond(mask, simd_index_type(k), imax);
893 }
894 // swap rows
895 for (size_type j=0; j<A.rows(); j++)
896 {
897 // This is a swap operation where the second operand is scattered,
898 // and on top of that is also extracted from deep within a
899 // moderately complicated data structure (a DenseMatrix), where we
900 // can't assume much on the memory layout. On intel processors,
901 // the only instruction that might help us here is vgather, but it
902 // is unclear whether that is even faster than a software
903 // implementation, and we would also need vscatter which does not
904 // exist. So break vectorization here and do it manually.
905 for(std::size_t l = 0; l < Simd::lanes(A[i][j]); ++l)
906 swap(Simd::lane(l, A[i][j]),
907 Simd::lane(l, A[Simd::lane(l, imax)][j]));
908 }
909 func.swap(i, imax); // swap the pivot or rhs
910 }
911
912 // singular ?
913 nonsingularLanes = nonsingularLanes && (pivmax != real_type(0));
914 if (throwEarly) {
915 if(!Simd::allTrue(nonsingularLanes))
916 DUNE_THROW(FMatrixError, "matrix is singular");
917 }
918 else { // !throwEarly
919 if(!Simd::anyTrue(nonsingularLanes))
920 return;
921 }
922
923 // eliminate
924 for (size_type k=i+1; k<A.rows(); k++)
925 {
926 // in the simd case, A[i][i] may be close to zero in some lanes. Pray
927 // that the result is no worse than a quiet NaN.
928 field_type factor = A[k][i]/A[i][i];
929 A[k][i] = factor;
930 for (size_type j=i+1; j<A.rows(); j++)
931 A[k][j] -= factor*A[i][j];
932 func(factor, k, i);
933 }
934 }
935 }
936
937 template<typename MAT>
938 template <class V1, class V2>
939 inline void DenseMatrix<MAT>::solve(V1& x, const V2& b, bool doPivoting) const
940 {
942 // never mind those ifs, because they get optimized away
943 if (rows()!=cols())
944 DUNE_THROW(FMatrixError, "Can't solve for a " << rows() << "x" << cols() << " matrix!");
945
946 if (rows()==1) {
947
948#ifdef DUNE_FMatrix_WITH_CHECKING
949 if (Simd::anyTrue(fvmeta::absreal((*this)[0][0])
951 DUNE_THROW(FMatrixError,"matrix is singular");
952#endif
953 x[0] = b[0]/(*this)[0][0];
954
955 }
956 else if (rows()==2) {
957
958 field_type detinv = (*this)[0][0]*(*this)[1][1]-(*this)[0][1]*(*this)[1][0];
959#ifdef DUNE_FMatrix_WITH_CHECKING
960 if (Simd::anyTrue(fvmeta::absreal(detinv)
962 DUNE_THROW(FMatrixError,"matrix is singular");
963#endif
964 detinv = real_type(1.0)/detinv;
965
966 x[0] = detinv*((*this)[1][1]*b[0]-(*this)[0][1]*b[1]);
967 x[1] = detinv*((*this)[0][0]*b[1]-(*this)[1][0]*b[0]);
968
969 }
970 else if (rows()==3) {
971
972 field_type d = determinant(doPivoting);
973#ifdef DUNE_FMatrix_WITH_CHECKING
974 if (Simd::anyTrue(fvmeta::absreal(d)
976 DUNE_THROW(FMatrixError,"matrix is singular");
977#endif
978
979 x[0] = (b[0]*(*this)[1][1]*(*this)[2][2] - b[0]*(*this)[2][1]*(*this)[1][2]
980 - b[1] *(*this)[0][1]*(*this)[2][2] + b[1]*(*this)[2][1]*(*this)[0][2]
981 + b[2] *(*this)[0][1]*(*this)[1][2] - b[2]*(*this)[1][1]*(*this)[0][2]) / d;
982
983 x[1] = ((*this)[0][0]*b[1]*(*this)[2][2] - (*this)[0][0]*b[2]*(*this)[1][2]
984 - (*this)[1][0] *b[0]*(*this)[2][2] + (*this)[1][0]*b[2]*(*this)[0][2]
985 + (*this)[2][0] *b[0]*(*this)[1][2] - (*this)[2][0]*b[1]*(*this)[0][2]) / d;
986
987 x[2] = ((*this)[0][0]*(*this)[1][1]*b[2] - (*this)[0][0]*(*this)[2][1]*b[1]
988 - (*this)[1][0] *(*this)[0][1]*b[2] + (*this)[1][0]*(*this)[2][1]*b[0]
989 + (*this)[2][0] *(*this)[0][1]*b[1] - (*this)[2][0]*(*this)[1][1]*b[0]) / d;
990
991 }
992 else {
993
994 V1& rhs = x; // use x to store rhs
995 rhs = b; // copy data
996 Elim<V1> elim(rhs);
997 AutonomousValue<MAT> A(asImp());
999 nonsingularLanes(true);
1000
1001 AutonomousValue<MAT>::luDecomposition(A, elim, nonsingularLanes, true, doPivoting);
1002
1003 // backsolve
1004 for(int i=rows()-1; i>=0; i--) {
1005 for (size_type j=i+1; j<rows(); j++)
1006 rhs[i] -= A[i][j]*x[j];
1007 x[i] = rhs[i]/A[i][i];
1008 }
1009 }
1010 }
1011
1012 template<typename MAT>
1013 inline void DenseMatrix<MAT>::invert(bool doPivoting)
1014 {
1015 using real_type = typename FieldTraits<MAT>::real_type;
1016 using std::swap;
1017
1018 // never mind those ifs, because they get optimized away
1019 if (rows()!=cols())
1020 DUNE_THROW(FMatrixError, "Can't invert a " << rows() << "x" << cols() << " matrix!");
1021
1022 if (rows()==1) {
1023
1024#ifdef DUNE_FMatrix_WITH_CHECKING
1025 if (Simd::anyTrue(fvmeta::absreal((*this)[0][0])
1027 DUNE_THROW(FMatrixError,"matrix is singular");
1028#endif
1029 (*this)[0][0] = real_type( 1 ) / (*this)[0][0];
1030
1031 }
1032 else if (rows()==2) {
1033
1034 field_type detinv = (*this)[0][0]*(*this)[1][1]-(*this)[0][1]*(*this)[1][0];
1035#ifdef DUNE_FMatrix_WITH_CHECKING
1036 if (Simd::anyTrue(fvmeta::absreal(detinv)
1038 DUNE_THROW(FMatrixError,"matrix is singular");
1039#endif
1040 detinv = real_type( 1 ) / detinv;
1041
1042 field_type temp=(*this)[0][0];
1043 (*this)[0][0] = (*this)[1][1]*detinv;
1044 (*this)[0][1] = -(*this)[0][1]*detinv;
1045 (*this)[1][0] = -(*this)[1][0]*detinv;
1046 (*this)[1][1] = temp*detinv;
1047
1048 }
1049 else if (rows()==3)
1050 {
1051 using K = field_type;
1052 // code generated by maple
1053 K t4 = (*this)[0][0] * (*this)[1][1];
1054 K t6 = (*this)[0][0] * (*this)[1][2];
1055 K t8 = (*this)[0][1] * (*this)[1][0];
1056 K t10 = (*this)[0][2] * (*this)[1][0];
1057 K t12 = (*this)[0][1] * (*this)[2][0];
1058 K t14 = (*this)[0][2] * (*this)[2][0];
1059
1060 K det = (t4*(*this)[2][2]-t6*(*this)[2][1]-t8*(*this)[2][2]+
1061 t10*(*this)[2][1]+t12*(*this)[1][2]-t14*(*this)[1][1]);
1062 K t17 = K(1.0)/det;
1063
1064 K matrix01 = (*this)[0][1];
1065 K matrix00 = (*this)[0][0];
1066 K matrix10 = (*this)[1][0];
1067 K matrix11 = (*this)[1][1];
1068
1069 (*this)[0][0] = ((*this)[1][1] * (*this)[2][2] - (*this)[1][2] * (*this)[2][1])*t17;
1070 (*this)[0][1] = -((*this)[0][1] * (*this)[2][2] - (*this)[0][2] * (*this)[2][1])*t17;
1071 (*this)[0][2] = (matrix01 * (*this)[1][2] - (*this)[0][2] * (*this)[1][1])*t17;
1072 (*this)[1][0] = -((*this)[1][0] * (*this)[2][2] - (*this)[1][2] * (*this)[2][0])*t17;
1073 (*this)[1][1] = (matrix00 * (*this)[2][2] - t14) * t17;
1074 (*this)[1][2] = -(t6-t10) * t17;
1075 (*this)[2][0] = (matrix10 * (*this)[2][1] - matrix11 * (*this)[2][0]) * t17;
1076 (*this)[2][1] = -(matrix00 * (*this)[2][1] - t12) * t17;
1077 (*this)[2][2] = (t4-t8) * t17;
1078 }
1079 else {
1080 using std::swap;
1081
1082 AutonomousValue<MAT> A(asImp());
1083 std::vector<simd_index_type> pivot(rows());
1085 nonsingularLanes(true);
1086 AutonomousValue<MAT>::luDecomposition(A, ElimPivot(pivot), nonsingularLanes, true, doPivoting);
1087 auto& L=A;
1088 auto& U=A;
1089
1090 // initialize inverse
1091 *this=field_type();
1092
1093 for(size_type i=0; i<rows(); ++i)
1094 (*this)[i][i]=1;
1095
1096 // L Y = I; multiple right hand sides
1097 for (size_type i=0; i<rows(); i++)
1098 for (size_type j=0; j<i; j++)
1099 for (size_type k=0; k<rows(); k++)
1100 (*this)[i][k] -= L[i][j]*(*this)[j][k];
1101
1102 // U A^{-1} = Y
1103 for (size_type i=rows(); i>0;) {
1104 --i;
1105 for (size_type k=0; k<rows(); k++) {
1106 for (size_type j=i+1; j<rows(); j++)
1107 (*this)[i][k] -= U[i][j]*(*this)[j][k];
1108 (*this)[i][k] /= U[i][i];
1109 }
1110 }
1111
1112 for(size_type i=rows(); i>0; ) {
1113 --i;
1114 for(std::size_t l = 0; l < Simd::lanes((*this)[0][0]); ++l)
1115 {
1116 std::size_t pi = Simd::lane(l, pivot[i]);
1117 if(i!=pi)
1118 for(size_type j=0; j<rows(); ++j)
1119 swap(Simd::lane(l, (*this)[j][pi]),
1120 Simd::lane(l, (*this)[j][ i]));
1121 }
1122 }
1123 }
1124 }
1125
1126 // implementation of the determinant
1127 template<typename MAT>
1128 inline typename DenseMatrix<MAT>::field_type
1129 DenseMatrix<MAT>::determinant(bool doPivoting) const
1130 {
1131 // never mind those ifs, because they get optimized away
1132 if (rows()!=cols())
1133 DUNE_THROW(FMatrixError, "There is no determinant for a " << rows() << "x" << cols() << " matrix!");
1134
1135 if (rows()==1)
1136 return (*this)[0][0];
1137
1138 if (rows()==2)
1139 return (*this)[0][0]*(*this)[1][1] - (*this)[0][1]*(*this)[1][0];
1140
1141 if (rows()==3) {
1142 // code generated by maple
1143 field_type t4 = (*this)[0][0] * (*this)[1][1];
1144 field_type t6 = (*this)[0][0] * (*this)[1][2];
1145 field_type t8 = (*this)[0][1] * (*this)[1][0];
1146 field_type t10 = (*this)[0][2] * (*this)[1][0];
1147 field_type t12 = (*this)[0][1] * (*this)[2][0];
1148 field_type t14 = (*this)[0][2] * (*this)[2][0];
1149
1150 return (t4*(*this)[2][2]-t6*(*this)[2][1]-t8*(*this)[2][2]+
1151 t10*(*this)[2][1]+t12*(*this)[1][2]-t14*(*this)[1][1]);
1152
1153 }
1154
1155 AutonomousValue<MAT> A(asImp());
1156 field_type det;
1158 nonsingularLanes(true);
1159
1160 AutonomousValue<MAT>::luDecomposition(A, ElimDet(det), nonsingularLanes, false, doPivoting);
1161 det = Simd::cond(nonsingularLanes, det, field_type(0));
1162
1163 for (size_type i = 0; i < rows(); ++i)
1164 det *= A[i][i];
1165 return det;
1166 }
1167
1168#endif // DOXYGEN
1169
1171
1173 template <typename MAT, typename V1, typename V2>
1174 static inline void multAssign(const DenseMatrix<MAT> &matrix, const DenseVector<V1> & x, DenseVector<V2> & ret)
1175 {
1176 DUNE_ASSERT_BOUNDS(x.size() == matrix.cols());
1177 DUNE_ASSERT_BOUNDS(ret.size() == matrix.rows());
1178 typedef typename DenseMatrix<MAT>::size_type size_type;
1179
1180 for(size_type i=0; i<matrix.rows(); ++i)
1181 {
1182 ret[i] = 0.0;
1183 for(size_type j=0; j<matrix.cols(); ++j)
1184 {
1185 ret[i] += matrix[i][j]*x[j];
1186 }
1187 }
1188 }
1189
1190#if 0
1192 template <typename K, int rows, int cols>
1193 static inline void multAssignTransposed( const FieldMatrix<K,rows,cols> &matrix, const FieldVector<K,rows> & x, FieldVector<K,cols> & ret)
1194 {
1195 typedef typename FieldMatrix<K,rows,cols>::size_type size_type;
1196
1197 for(size_type i=0; i<cols(); ++i)
1198 {
1199 ret[i] = 0.0;
1200 for(size_type j=0; j<rows(); ++j)
1201 ret[i] += matrix[j][i]*x[j];
1202 }
1203 }
1204
1206 template <typename K, int rows, int cols>
1207 static inline FieldVector<K,rows> mult(const FieldMatrix<K,rows,cols> &matrix, const FieldVector<K,cols> & x)
1208 {
1209 FieldVector<K,rows> ret;
1210 multAssign(matrix,x,ret);
1211 return ret;
1212 }
1213
1215 template <typename K, int rows, int cols>
1216 static inline FieldVector<K,cols> multTransposed(const FieldMatrix<K,rows,cols> &matrix, const FieldVector<K,rows> & x)
1217 {
1218 FieldVector<K,cols> ret;
1219 multAssignTransposed( matrix, x, ret );
1220 return ret;
1221 }
1222#endif
1223
1224 } // end namespace DenseMatrixHelp
1225
1227 template<typename MAT>
1228 std::ostream& operator<< (std::ostream& s, const DenseMatrix<MAT>& a)
1229 {
1230 for (typename DenseMatrix<MAT>::size_type i=0; i<a.rows(); i++)
1231 s << a[i] << std::endl;
1232 return s;
1233 }
1234
1236
1237} // end namespace Dune
1238
1239#endif
Traits for type conversions and type information.
Implements a scalar vector view wrapper around an existing scalar.
Various precision settings for calculations with FieldMatrix and FieldVector.
Some useful basic math stuff.
Implements a vector constructed from a given type representing a field and a compile-time given size.
A few common exception classes.
A free function to provide the demangled class name of a given object or type as a string.
Macro for wrapping boundary checks.
#define DUNE_ASSERT_BOUNDS(cond)
If DUNE_CHECK_BOUNDS is defined: check if condition cond holds; otherwise, do nothing.
Definition boundschecking.hh:30
typename AutonomousValueType< T >::type AutonomousValue
Type free of internal references that T can be converted to.
Definition typetraits.hh:566
std::ostream & operator<<(std::ostream &s, const bigunsignedint< k > &x)
Definition bigunsignedint.hh:301
#define DUNE_THROW(E,...)
Definition exceptions.hh:314
bool anyTrue(const Mask &mask)
Whether any entry is true.
Definition simd/interface.hh:429
V cond(M &&mask, const V &ifTrue, const V &ifFalse)
Like the ?: operator.
Definition simd/interface.hh:386
bool allTrue(const Mask &mask)
Whether all entries are true.
Definition simd/interface.hh:439
typename Overloads::ScalarType< std::decay_t< V > >::type Scalar
Element type of some SIMD type.
Definition simd/interface.hh:235
Rebind< bool, V > Mask
Mask type type of some SIMD type.
Definition simd/interface.hh:289
typename Overloads::RebindType< std::decay_t< S >, std::decay_t< V > >::type Rebind
Construct SIMD type with different scalar type.
Definition simd/interface.hh:252
constexpr std::size_t lanes()
Number of lanes in a SIMD type.
Definition simd/interface.hh:305
decltype(auto) lane(std::size_t l, V &&v)
Extract an element of a SIMD type.
Definition simd/interface.hh:324
Mask< V > mask(ADLTag< 0, std::is_same< V, Mask< V > >::value >, const V &v)
implements Simd::mask()
Definition defaults.hh:153
Dune namespace
Definition alignedallocator.hh:13
constexpr int sign(const T &val)
Return the sign of the value.
Definition math.hh:162
void swap(T &v1, T &v2, bool mask)
Definition simd.hh:472
constexpr K conjugateComplex(const K &x)
compute conjugate complex of x
Definition math.hh:146
Definition densematrix.hh:1170
static void multAssign(const DenseMatrix< MAT > &matrix, const DenseVector< V1 > &x, DenseVector< V2 > &ret)
calculates ret = matrix * x
Definition densematrix.hh:1174
constexpr T abs(T t)
Definition cmath.hh:27
Whether this type acts as a scalar in the context of (hierarchically blocked) containers.
Definition typetraits.hh:194
A dense n x m matrix.
Definition densematrix.hh:145
ConstIterator const_iterator
typedef for stl compliant access
Definition densematrix.hh:249
constexpr derived_type & operator+=(const DenseMatrix< Other > &x)
vector space addition
Definition densematrix.hh:294
constexpr Iterator beforeBegin()
Definition densematrix.hh:241
constexpr void umhv(const X &x, Y &y) const
y += A^H x
Definition densematrix.hh:437
constexpr derived_type & operator*=(const field_type &k)
vector space multiplication with scalar
Definition densematrix.hh:326
constexpr void usmhv(const typename FieldTraits< Y >::field_type &alpha, const X &x, Y &y) const
y += alpha A^H x
Definition densematrix.hh:517
constexpr derived_type & operator-=(const DenseMatrix< Other > &x)
vector space subtraction
Definition densematrix.hh:317
void solve(V1 &x, const V2 &b, bool doPivoting=true) const
Solve system A x = b.
Traits::value_type field_type
export the type representing the field
Definition densematrix.hh:165
constexpr ConstIterator beforeEnd() const
Definition densematrix.hh:269
constexpr bool operator!=(const DenseMatrix< Other > &x) const
Binary matrix incomparison.
Definition densematrix.hh:363
constexpr derived_type & axpy(const field_type &a, const DenseMatrix< Other > &x)
vector space axpy operation (*this += a x)
Definition densematrix.hh:343
std::remove_reference< const_row_reference >::type::ConstIterator ConstColIterator
rename the iterators for easier access
Definition densematrix.hh:253
void invert(bool doPivoting=true)
Compute inverse.
constexpr void usmtv(const typename FieldTraits< Y >::field_type &alpha, const X &x, Y &y) const
y += alpha A^T x
Definition densematrix.hh:503
static void luDecomposition(DenseMatrix< MAT > &A, Func func, Mask &nonsingularLanes, bool throwEarly, bool doPivoting)
do an LU-Decomposition on matrix A
Traits::value_type block_type
export the type representing the components
Definition densematrix.hh:168
constexpr size_type cols() const
number of columns
Definition densematrix.hh:720
constexpr size_type size() const
size method (number of rows)
Definition densematrix.hh:205
constexpr size_type M() const
Definition densematrix.hh:708
constexpr void mmhv(const X &x, Y &y) const
y -= A^H x
Definition densematrix.hh:476
MAT & rightmultiply(const DenseMatrix< M2 > &M)
Multiplies M from the right to this matrix.
Definition densematrix.hh:650
constexpr Iterator beforeEnd()
Definition densematrix.hh:234
constexpr ConstIterator beforeBegin() const
Definition densematrix.hh:276
constexpr FieldTraits< value_type >::real_type frobenius_norm() const
frobenius norm: sqrt(sum over squared values of entries)
Definition densematrix.hh:533
Traits::value_type value_type
export the type representing the field
Definition densematrix.hh:162
constexpr Iterator end()
end iterator
Definition densematrix.hh:227
constexpr derived_type & operator/=(const field_type &k)
vector space division by scalar
Definition densematrix.hh:334
constexpr derived_type operator-() const
Matrix negation.
Definition densematrix.hh:303
constexpr size_type rows() const
number of rows
Definition densematrix.hh:714
constexpr void mv(const X &x, Y &y) const
y = A x
Definition densematrix.hh:373
MAT & leftmultiply(const DenseMatrix< M2 > &M)
Multiplies M from the left to this matrix.
Definition densematrix.hh:632
constexpr void mtv(const X &x, Y &y) const
y = A^T x
Definition densematrix.hh:392
constexpr FieldTraits< vt >::real_type infinity_norm() const
infinity norm (row sum norm, how to generalize for blocks?)
Definition densematrix.hh:551
constexpr row_reference operator[](size_type i)
random access
Definition densematrix.hh:194
constexpr void umv(const X &x, Y &y) const
y += A x
Definition densematrix.hh:411
Traits::derived_type derived_type
type of derived matrix class
Definition densematrix.hh:159
constexpr bool exists(size_type i, size_type j) const
return true when (i,j) is in pattern
Definition densematrix.hh:728
constexpr void usmv(const typename FieldTraits< Y >::field_type &alpha, const X &x, Y &y) const
y += alpha A x
Definition densematrix.hh:489
Iterator RowIterator
rename the iterators for easier access
Definition densematrix.hh:216
static constexpr int blocklevel
Definition densematrix.hh:183
constexpr ConstIterator begin() const
begin iterator
Definition densematrix.hh:256
constexpr void umtv(const X &x, Y &y) const
y += A^T x
Definition densematrix.hh:424
DenseIterator< const DenseMatrix, const row_type, const_row_reference > ConstIterator
Iterator class for sequential access.
Definition densematrix.hh:247
constexpr void mmtv(const X &x, Y &y) const
y -= A^T x
Definition densematrix.hh:463
constexpr ConstIterator end() const
end iterator
Definition densematrix.hh:262
constexpr FieldTraits< value_type >::real_type frobenius_norm2() const
square of frobenius norm, need for block recursion
Definition densematrix.hh:541
Traits::row_type row_type
The type used to represent a row (must fulfill the Dune::DenseVector interface).
Definition densematrix.hh:174
constexpr size_type N() const
Definition densematrix.hh:702
constexpr bool operator==(const DenseMatrix< Other > &x) const
Binary matrix comparison.
Definition densematrix.hh:353
Traits::size_type size_type
The type used for the index access and size operation.
Definition densematrix.hh:171
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:180
std::remove_reference< row_reference >::type::Iterator ColIterator
rename the iterators for easier access
Definition densematrix.hh:218
Traits::row_reference row_reference
The type used to represent a reference to a row (usually row_type &).
Definition densematrix.hh:177
friend class DenseMatrix
Definition densematrix.hh:153
constexpr void mmv(const X &x, Y &y) const
y -= A x
Definition densematrix.hh:450
Iterator iterator
typedef for stl compliant access
Definition densematrix.hh:214
constexpr FieldTraits< vt >::real_type infinity_norm_real() const
simplified infinity norm (uses Manhattan norm for complex values)
Definition densematrix.hh:566
ConstIterator ConstRowIterator
rename the iterators for easier access
Definition densematrix.hh:251
constexpr Iterator begin()
begin iterator
Definition densematrix.hh:221
constexpr derived_type & operator=(const RHS &rhs)
Definition densematrix.hh:284
DenseIterator< DenseMatrix, row_type, row_reference > Iterator
Iterator class for sequential access.
Definition densematrix.hh:212
field_type determinant(bool doPivoting=true) const
calculates the determinant of this matrix
const FieldTraits< typenameDenseMatVecTraits< M >::value_type >::real_type real_type
Definition densematrix.hh:36
const FieldTraits< typenameDenseMatVecTraits< M >::value_type >::field_type field_type
Definition densematrix.hh:35
A dense n x m matrix.
Definition fmatrix.hh:117
Base::size_type size_type
Definition fmatrix.hh:128
vector space out of a tensor product of fields.
Definition fvector.hh:97
you have to specialize this structure for any type that should be assignable to a DenseMatrix
Definition densematrix.hh:61
Error thrown if operations of a FieldMatrix fail.
Definition densematrix.hh:131
Interface for a class of dense vectors over a given field.
Definition densevector.hh:244
constexpr size_type size() const
size method
Definition densevector.hh:351
Generic iterator class for dense vector and matrix implementations.
Definition densevector.hh:132
Default exception class for mathematical errors.
Definition exceptions.hh:335
Definition ftraits.hh:26
T field_type
export the type representing the field
Definition ftraits.hh:28
T real_type
export the type representing the real type of the field
Definition ftraits.hh:30
Definition matvectraits.hh:31
static ctype absolute_limit()
return threshold to declare matrix singular
Definition precision.hh:28
Include file for users of the SIMD abstraction layer.