Dune Core Modules (2.8.0)

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