Dune Core Modules (unstable)

diagonalmatrix.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_DIAGONAL_MATRIX_HH
6#define DUNE_DIAGONAL_MATRIX_HH
7
12#include <algorithm>
13#include <cassert>
14#include <cmath>
15#include <complex>
16#include <cstddef>
17#include <initializer_list>
18#include <iostream>
19#include <memory>
20
28
29
30namespace Dune {
31
32 template< class K, int n > class DiagonalRowVectorConst;
33 template< class K, int n > class DiagonalRowVector;
34 template< class DiagonalMatrixType > class DiagonalMatrixWrapper;
35 template< class C, class T, class R> class ContainerWrapperIterator;
36
51 template<class K, int n>
53 {
54 typedef DiagonalMatrixWrapper< DiagonalMatrix<K,n> > WrapperType;
55
56 public:
57 //===== type definitions and constants
58
60 typedef K value_type;
61 typedef value_type field_type;
62
64 typedef K block_type;
65
67 typedef std::size_t size_type;
68
70 constexpr static int blocklevel = 1;
71
73 typedef DiagonalRowVector<K,n> row_type;
74 typedef row_type reference;
75 typedef row_type row_reference;
76 typedef DiagonalRowVectorConst<K,n> const_row_type;
77 typedef const_row_type const_reference;
78 typedef const_row_type const_row_reference;
79
81 constexpr static int rows = n;
83 constexpr static int cols = n;
84
85 //==== size
86
87 static constexpr size_type size ()
88 {
89 return rows;
90 }
91
92 //===== constructors
93
95 constexpr DiagonalMatrix() = default;
96
98 DiagonalMatrix (const K& k)
99 : diag_(k)
100 {}
101
104 : diag_(diag)
105 {}
106
115 DiagonalMatrix (std::initializer_list<K> const &l)
116 {
117 std::copy_n(l.begin(), std::min(static_cast<std::size_t>(rows),
118 l.size()),
119 diag_.begin());
120 }
121
124 {
125 diag_ = k;
126 return *this;
127 }
128
130 bool identical(const DiagonalMatrix<K,n>& other) const
131 {
132 return (this==&other);
133 }
134
137 {
138 return *this;
139 }
140
141 //===== iterator interface to rows of the matrix
150
153 {
154 return Iterator(WrapperType(this),0);
155 }
156
159 {
160 return Iterator(WrapperType(this),n);
161 }
162
166 {
167 return Iterator(WrapperType(this),n-1);
168 }
169
173 {
174 return Iterator(WrapperType(this),-1);
175 }
176
177
186
189 {
190 return ConstIterator(WrapperType(this),0);
191 }
192
195 {
196 return ConstIterator(WrapperType(this),n);
197 }
198
202 {
203 return ConstIterator(WrapperType(this),n-1);
204 }
205
209 {
210 return ConstIterator(WrapperType(this),-1);
211 }
212
213
214
215 //===== vector space arithmetic
216
219 {
220 diag_ += y.diag_;
221 return *this;
222 }
223
226 {
227 diag_ -= y.diag_;
228 return *this;
229 }
230
233 {
234 diag_ += k;
235 return *this;
236 }
237
240 {
241 diag_ -= k;
242 return *this;
243 }
244
247 {
248 diag_ *= k;
249 return *this;
250 }
251
254 {
255 diag_ /= k;
256 return *this;
257 }
258
259 //===== comparison ops
260
262 bool operator==(const DiagonalMatrix& other) const
263 {
264 return diag_==other.diagonal();
265 }
266
268 bool operator!=(const DiagonalMatrix& other) const
269 {
270 return diag_!=other.diagonal();
271 }
272
273
274 //===== linear maps
275
277 template<class X, class Y>
278 void mv (const X& x, Y& y) const
279 {
280#ifdef DUNE_FMatrix_WITH_CHECKING
281 if (x.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
282 if (y.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
283#endif
284 for (size_type i=0; i<n; ++i)
285 y[i] = diag_[i] * x[i];
286 }
287
289 template<class X, class Y>
290 void mtv (const X& x, Y& y) const
291 {
292 mv(x, y);
293 }
294
296 template<class X, class Y>
297 void umv (const X& x, Y& y) const
298 {
299#ifdef DUNE_FMatrix_WITH_CHECKING
300 if (x.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
301 if (y.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
302#endif
303 for (size_type i=0; i<n; ++i)
304 y[i] += diag_[i] * x[i];
305 }
306
308 template<class X, class Y>
309 void umtv (const X& x, Y& y) const
310 {
311#ifdef DUNE_FMatrix_WITH_CHECKING
312 if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
313 if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
314#endif
315 for (size_type i=0; i<n; ++i)
316 y[i] += diag_[i] * x[i];
317 }
318
320 template<class X, class Y>
321 void umhv (const X& x, Y& y) const
322 {
323#ifdef DUNE_FMatrix_WITH_CHECKING
324 if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
325 if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
326#endif
327 for (size_type i=0; i<n; i++)
328 y[i] += conjugateComplex(diag_[i])*x[i];
329 }
330
332 template<class X, class Y>
333 void mmv (const X& x, Y& y) const
334 {
335#ifdef DUNE_FMatrix_WITH_CHECKING
336 if (x.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
337 if (y.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
338#endif
339 for (size_type i=0; i<n; ++i)
340 y[i] -= diag_[i] * x[i];
341 }
342
344 template<class X, class Y>
345 void mmtv (const X& x, Y& y) const
346 {
347#ifdef DUNE_FMatrix_WITH_CHECKING
348 if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
349 if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
350#endif
351 for (size_type i=0; i<n; ++i)
352 y[i] -= diag_[i] * x[i];
353 }
354
356 template<class X, class Y>
357 void mmhv (const X& x, Y& y) const
358 {
359#ifdef DUNE_FMatrix_WITH_CHECKING
360 if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
361 if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
362#endif
363 for (size_type i=0; i<n; i++)
364 y[i] -= conjugateComplex(diag_[i])*x[i];
365 }
366
368 template<class X, class Y>
369 void usmv (const typename FieldTraits<Y>::field_type & alpha,
370 const X& x, Y& y) const
371 {
372#ifdef DUNE_FMatrix_WITH_CHECKING
373 if (x.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
374 if (y.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
375#endif
376 for (size_type i=0; i<n; i++)
377 y[i] += alpha * diag_[i] * x[i];
378 }
379
381 template<class X, class Y>
382 void usmtv (const typename FieldTraits<Y>::field_type & alpha,
383 const X& x, Y& y) const
384 {
385#ifdef DUNE_FMatrix_WITH_CHECKING
386 if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
387 if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
388#endif
389 for (size_type i=0; i<n; i++)
390 y[i] += alpha * diag_[i] * x[i];
391 }
392
394 template<class X, class Y>
395 void usmhv (const typename FieldTraits<Y>::field_type & alpha,
396 const X& x, Y& y) const
397 {
398#ifdef DUNE_FMatrix_WITH_CHECKING
399 if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
400 if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
401#endif
402 for (size_type i=0; i<n; i++)
403 y[i] += alpha * conjugateComplex(diag_[i]) * x[i];
404 }
405
406 //===== norms
407
409 double frobenius_norm () const
410 {
411 return diag_.two_norm();
412 }
413
415 double frobenius_norm2 () const
416 {
417 return diag_.two_norm2();
418 }
419
421 double infinity_norm () const
422 {
423 return diag_.infinity_norm();
424 }
425
427 double infinity_norm_real () const
428 {
429 return diag_.infinity_norm_real();
430 }
431
432
433
434 //===== solve
435
437 template<class V>
438 void solve (V& x, const V& b) const
439 {
440 for (int i=0; i<n; i++)
441 x[i] = b[i]/diag_[i];
442 }
443
445 void invert()
446 {
447 using real_type = typename FieldTraits<K>::real_type;
448 for (int i=0; i<n; i++)
449 diag_[i] = real_type(1.0)/diag_[i];
450 }
451
453 K determinant () const
454 {
455 K det = diag_[0];
456 for (int i=1; i<n; i++)
457 det *= diag_[i];
458 return det;
459 }
460
461
462
463 //===== matrix-matrix multiplication
464
467 template <class OtherScalar>
468 friend auto operator* ( const DiagonalMatrix& matrixA,
469 const DiagonalMatrix<OtherScalar, n>& matrixB)
470 {
472 for(int i=0; i<n; ++i)
473 result.diagonal(i) = matrixA.diagonal(i)*matrixB.diagonal(i);
474 return result;
475 }
476
477
478
479 //===== sizes
480
482 static constexpr size_type N ()
483 {
484 return n;
485 }
486
488 static constexpr size_type M ()
489 {
490 return n;
491 }
492
493
494
495 //===== query
496
498 bool exists (size_type i, size_type j) const
499 {
500 DUNE_ASSERT_BOUNDS(i >= 0 && i < n);
501 DUNE_ASSERT_BOUNDS(j >= 0 && j < n);
502 return i==j;
503 }
504
505
506
508 friend std::ostream& operator<< (std::ostream& s, const DiagonalMatrix<K,n>& a)
509 {
510 for (size_type i=0; i<n; i++) {
511 for (size_type j=0; j<n; j++)
512 s << ((i==j) ? a.diag_[i] : 0) << " ";
513 s << std::endl;
514 }
515 return s;
516 }
517
520 {
521 return reference(const_cast<K*>(&diag_[i]), i);
522 }
523
525 const_reference operator[](size_type i) const
526 {
527 return const_reference(const_cast<K*>(&diag_[i]), i);
528 }
529
531 const K& diagonal(size_type i) const
532 {
533 return diag_[i];
534 }
535
538 {
539 return diag_[i];
540 }
541
544 {
545 return diag_;
546 }
547
550 {
551 return diag_;
552 }
553
554 private:
555
556 // the data, a FieldVector storing the diagonal
557 FieldVector<K,n> diag_;
558 };
559
560 template< class K, int n >
561 struct FieldTraits< DiagonalMatrix<K,n> >
562 {
563 typedef typename FieldTraits<K>::field_type field_type;
564 typedef typename FieldTraits<K>::real_type real_type;
565 };
566
567
568#ifndef DOXYGEN // hide specialization
571 template< class K >
572 class DiagonalMatrix<K, 1> : public FieldMatrix<K, 1, 1>
573 {
574 typedef FieldMatrix<K,1,1> Base;
575 public:
577 typedef typename Base::size_type size_type;
578
581 constexpr static int blocklevel = 1;
582
583 typedef typename Base::row_type row_type;
584
585 typedef typename Base::row_reference row_reference;
586 typedef typename Base::const_row_reference const_row_reference;
587
590 constexpr static int rows = 1;
593 constexpr static int cols = 1;
594
595
597 constexpr DiagonalMatrix() = default;
598
600 DiagonalMatrix(const K& scalar)
601 {
602 (*this)[0][0] = scalar;
603 }
604
606 const K& diagonal(size_type) const
607 {
608 return (*this)[0][0];
609 }
610
613 {
614 return (*this)[0][0];
615 }
616
618 const FieldVector<K,1>& diagonal() const
619 {
620 return (*this)[0];
621 }
622
624 FieldVector<K,1>& diagonal()
625 {
626 return (*this)[0];
627 }
628
630 DiagonalMatrix<K, 1> transposed() const
631 {
632 return *this;
633 }
634
637 template <class OtherScalar>
638 friend auto operator* ( const DiagonalMatrix& matrixA,
639 const DiagonalMatrix<OtherScalar, 1>& matrixB)
640 {
641 return DiagonalMatrix<typename PromotionTraits<K,OtherScalar>::PromotedType, 1>{matrixA.diagonal(0)*matrixB.diagonal(0)};
642 }
643
644 };
645#endif
646
647
648 template<class DiagonalMatrixType>
649 class DiagonalMatrixWrapper
650 {
651 typedef typename DiagonalMatrixType::reference reference;
652 typedef typename DiagonalMatrixType::const_reference const_reference;
653 typedef typename DiagonalMatrixType::field_type K;
654 typedef DiagonalRowVector<K, DiagonalMatrixType::rows> row_type;
655 typedef std::size_t size_type;
656 typedef DiagonalMatrixWrapper< DiagonalMatrixType> MyType;
657
658 friend class ContainerWrapperIterator<const MyType, reference, reference>;
659 friend class ContainerWrapperIterator<const MyType, const_reference, const_reference>;
660
661 public:
662
663 DiagonalMatrixWrapper() :
664 mat_(0)
665 {}
666
667 DiagonalMatrixWrapper(const DiagonalMatrixType* mat) :
668 mat_(const_cast<DiagonalMatrixType*>(mat))
669 {}
670
671 size_type realIndex(int i) const
672 {
673 return i;
674 }
675
676 row_type* pointer(int i) const
677 {
678 row_ = row_type(&(mat_->diagonal(i)), i);
679 return &row_;
680 }
681
682 bool identical(const DiagonalMatrixWrapper& other) const
683 {
684 return mat_==other.mat_;
685 }
686
687 private:
688
689 mutable DiagonalMatrixType* mat_;
690 mutable row_type row_;
691 };
692
696 template< class K, int n >
697 class DiagonalRowVectorConst
698 {
699 template<class DiagonalMatrixType>
700 friend class DiagonalMatrixWrapper;
701 friend class ContainerWrapperIterator<DiagonalRowVectorConst<K,n>, const K, const K&>;
702
703 public:
704 // remember size of vector
705 constexpr static int dimension = n;
706
707 // standard constructor and everything is sufficient ...
708
709 //===== type definitions and constants
710
712 typedef K field_type;
713
715 typedef K block_type;
716
718 typedef std::size_t size_type;
719
721 constexpr static int blocklevel = 1;
722
724 constexpr static int size = n;
725
728 p_(0),
729 row_(0)
730 {}
731
733 explicit DiagonalRowVectorConst (K* p, int col) :
734 p_(p),
735 row_(col)
736 {}
737
738 //===== access to components
739
741 const K& operator[] ([[maybe_unused]] size_type i) const
742 {
743 DUNE_ASSERT_BOUNDS(i == row_);
744 return *p_;
745 }
746
747 // check if row is identical to other row (not only identical values)
748 // since this is a proxy class we need to check equality of the stored pointer
749 bool identical(const DiagonalRowVectorConst<K,n>& other) const
750 {
751 return ((p_ == other.p_)and (row_ == other.row_));
752 }
753
758
761 {
762 return ConstIterator(*this,0);
763 }
764
767 {
768 return ConstIterator(*this,1);
769 }
770
774 {
775 return ConstIterator(*this,0);
776 }
777
781 {
782 return ConstIterator(*this,-1);
783 }
784
786 bool operator== (const DiagonalRowVectorConst& y) const
787 {
788 return ((p_==y.p_)and (row_==y.row_));
789 }
790
791 //===== sizes
792
794 size_type N () const
795 {
796 return n;
797 }
798
800 size_type dim () const
801 {
802 return n;
803 }
804
807 {
808 return row_;
809 }
810
812 const K& diagonal() const
813 {
814 return *p_;
815 }
816
817 protected:
818
819 size_type realIndex([[maybe_unused]] int i) const
820 {
821 return rowIndex();
822 }
823
824 K* pointer([[maybe_unused]] size_type i) const
825 {
826 return const_cast<K*>(p_);
827 }
828
829 DiagonalRowVectorConst* operator&()
830 {
831 return this;
832 }
833
834 // the data, very simply a pointer to the diagonal value and the row number
835 K* p_;
836 size_type row_;
837 };
838
839 template< class K, int n >
840 class DiagonalRowVector : public DiagonalRowVectorConst<K,n>
841 {
842 template<class DiagonalMatrixType>
843 friend class DiagonalMatrixWrapper;
844 friend class ContainerWrapperIterator<DiagonalRowVector<K,n>, K, K&>;
845
846 public:
847 // standard constructor and everything is sufficient ...
848
849 //===== type definitions and constants
850
852 typedef K field_type;
853
855 typedef K block_type;
856
858 typedef std::size_t size_type;
859
861 DiagonalRowVector() : DiagonalRowVectorConst<K,n>()
862 {}
863
865 explicit DiagonalRowVector (K* p, int col) : DiagonalRowVectorConst<K,n>(p, col)
866 {}
867
868 //===== assignment from scalar
870 DiagonalRowVector& operator= (const K& k)
871 {
872 *p_ = k;
873 return *this;
874 }
875
876 //===== access to components
877
879 K& operator[] ([[maybe_unused]] size_type i)
880 {
881 DUNE_ASSERT_BOUNDS(i == row_);
882 return *p_;
883 }
884
889
892 {
893 return Iterator(*this, 0);
894 }
895
898 {
899 return Iterator(*this, 1);
900 }
901
905 {
906 return Iterator(*this, 0);
907 }
908
912 {
913 return Iterator(*this, -1);
914 }
915
920
921 using DiagonalRowVectorConst<K,n>::identical;
922 using DiagonalRowVectorConst<K,n>::operator[];
923 using DiagonalRowVectorConst<K,n>::operator==;
924 using DiagonalRowVectorConst<K,n>::begin;
925 using DiagonalRowVectorConst<K,n>::end;
926 using DiagonalRowVectorConst<K,n>::beforeEnd;
927 using DiagonalRowVectorConst<K,n>::beforeBegin;
928 using DiagonalRowVectorConst<K,n>::N;
929 using DiagonalRowVectorConst<K,n>::dim;
930 using DiagonalRowVectorConst<K,n>::rowIndex;
931 using DiagonalRowVectorConst<K,n>::diagonal;
932
933 protected:
934
935 DiagonalRowVector* operator&()
936 {
937 return this;
938 }
939
940 private:
941
942 using DiagonalRowVectorConst<K,n>::p_;
943 using DiagonalRowVectorConst<K,n>::row_;
944 };
945
946
947 // implement type traits
948 template<class K, int n>
949 struct const_reference< DiagonalRowVector<K,n> >
950 {
951 typedef DiagonalRowVectorConst<K,n> type;
952 };
953
954 template<class K, int n>
955 struct const_reference< DiagonalRowVectorConst<K,n> >
956 {
957 typedef DiagonalRowVectorConst<K,n> type;
958 };
959
960 template<class K, int n>
961 struct mutable_reference< DiagonalRowVector<K,n> >
962 {
963 typedef DiagonalRowVector<K,n> type;
964 };
965
966 template<class K, int n>
967 struct mutable_reference< DiagonalRowVectorConst<K,n> >
968 {
969 typedef DiagonalRowVector<K,n> type;
970 };
971
972
973
996 template<class CW, class T, class R>
997 class ContainerWrapperIterator : public BidirectionalIteratorFacade<ContainerWrapperIterator<CW,T,R>,T, R, int>
998 {
999 typedef typename std::remove_const<CW>::type NonConstCW;
1000
1001 friend class ContainerWrapperIterator<CW, typename mutable_reference<T>::type, typename mutable_reference<R>::type>;
1002 friend class ContainerWrapperIterator<CW, typename const_reference<T>::type, typename const_reference<R>::type>;
1003
1004 typedef ContainerWrapperIterator<CW, typename mutable_reference<T>::type, typename mutable_reference<R>::type> MyType;
1005 typedef ContainerWrapperIterator<CW, typename const_reference<T>::type, typename const_reference<R>::type> MyConstType;
1006
1007 public:
1008
1009 // Constructors needed by the facade iterators.
1011 containerWrapper_(),
1012 position_(0)
1013 {}
1014
1015 ContainerWrapperIterator(CW containerWrapper, int position) :
1016 containerWrapper_(containerWrapper),
1017 position_(position)
1018 {}
1019
1020 template<class OtherContainerWrapperIteratorType>
1021 ContainerWrapperIterator(OtherContainerWrapperIteratorType& other) :
1022 containerWrapper_(other.containerWrapper_),
1023 position_(other.position_)
1024 {}
1025
1026 ContainerWrapperIterator(const MyType& other) :
1027 containerWrapper_(other.containerWrapper_),
1028 position_(other.position_)
1029 {}
1030
1032 containerWrapper_(other.containerWrapper_),
1033 position_(other.position_)
1034 {}
1035
1036 template<class OtherContainerWrapperIteratorType>
1037 ContainerWrapperIterator& operator=(OtherContainerWrapperIteratorType& other)
1038 {
1039 containerWrapper_ = other.containerWrapper_;
1040 position_ = other.position_;
1041 return *this;
1042 }
1043
1044 // This operator is needed since we can not get the address of the
1045 // temporary object returned by dereference
1046 T* operator->() const
1047 {
1048 return containerWrapper_.pointer(position_);
1049 }
1050
1051 // Methods needed by the forward iterator
1052 bool equals(const MyType& other) const
1053 {
1054 return position_ == other.position_ && containerWrapper_.identical(other.containerWrapper_);
1055 }
1056
1057 bool equals(const MyConstType& other) const
1058 {
1059 return position_ == other.position_ && containerWrapper_.identical(other.containerWrapper_);
1060 }
1061
1062 R dereference() const
1063 {
1064 return *containerWrapper_.pointer(position_);
1065 }
1066
1067 void increment()
1068 {
1069 ++position_;
1070 }
1071
1072 // Additional function needed by BidirectionalIterator
1073 void decrement()
1074 {
1075 --position_;
1076 }
1077
1078 // Additional function needed by RandomAccessIterator
1079 R elementAt(int i) const
1080 {
1081 return *containerWrapper_.pointer(position_+i);
1082 }
1083
1084 void advance(int n)
1085 {
1086 position_=position_+n;
1087 }
1088
1089 template<class OtherContainerWrapperIteratorType>
1090 std::ptrdiff_t distanceTo(OtherContainerWrapperIteratorType& other) const
1091 {
1092 assert(containerWrapper_.identical(other));
1093 return other.position_ - position_;
1094 }
1095
1096 std::ptrdiff_t index() const
1097 {
1098 return containerWrapper_.realIndex(position_);
1099 }
1100
1101 private:
1102 NonConstCW containerWrapper_;
1103 size_t position_;
1104 };
1105
1106 template <class DenseMatrix, class field, int N>
1108 static void apply(DenseMatrix& denseMatrix,
1109 DiagonalMatrix<field, N> const& rhs) {
1110 DUNE_ASSERT_BOUNDS(denseMatrix.M() == N);
1111 DUNE_ASSERT_BOUNDS(denseMatrix.N() == N);
1112 denseMatrix = field(0);
1113 for (int i = 0; i < N; ++i)
1114 denseMatrix[i][i] = rhs.diagonal()[i];
1115 }
1116 };
1117 /* @} */
1118} // end namespace
1119#endif
Macro for wrapping boundary checks.
Facade class for stl conformant bidirectional iterators.
Definition: iteratorfacades.hh:275
Iterator class for sparse vector-like containers.
Definition: diagonalmatrix.hh:998
A dense n x m matrix.
Definition: densematrix.hh:140
constexpr size_type M() const
number of columns
Definition: densematrix.hh:703
constexpr size_type N() const
number of rows
Definition: densematrix.hh:697
FieldTraits< value_type >::real_type two_norm2() const
square of two norm (sum over squared values of entries), need for block recursion
Definition: densevector.hh:650
FieldTraits< value_type >::real_type two_norm() const
two norm sqrt(sum over squared values of entries)
Definition: densevector.hh:641
Iterator begin()
begin iterator
Definition: densevector.hh:347
FieldTraits< vt >::real_type infinity_norm() const
infinity norm (maximum of absolute values of entries)
Definition: densevector.hh:661
FieldTraits< vt >::real_type infinity_norm_real() const
simplified infinity norm (uses Manhattan norm for complex values)
Definition: densevector.hh:677
A diagonal matrix of static size.
Definition: diagonalmatrix.hh:53
Error thrown if operations of a FieldMatrix fail.
Definition: densematrix.hh:126
Implements a matrix constructed from a given type representing a field and a compile-time given numbe...
A few common exception classes.
Implements a matrix constructed from a given type representing a field and compile-time given number ...
Implements a vector constructed from a given type representing a field and a compile-time given size.
Implements a generic iterator class for writing stl conformant iterators.
#define DUNE_ASSERT_BOUNDS(cond)
If DUNE_CHECK_BOUNDS is defined: check if condition cond holds; otherwise, do nothing.
Definition: boundschecking.hh:30
ConstIterator beforeBegin() const
Definition: diagonalmatrix.hh:208
DiagonalMatrix< K, n > transposed() const
Return transposed of the matrix as DiagonalMatrix.
Definition: diagonalmatrix.hh:136
void mmhv(const X &x, Y &y) const
y -= A^H x
Definition: diagonalmatrix.hh:357
DiagonalMatrix & operator*=(const K &k)
vector space multiplication with scalar
Definition: diagonalmatrix.hh:246
std::size_t size_type
The type used for the index access and size operations.
Definition: diagonalmatrix.hh:67
size_type dim() const
dimension of the vector space
Definition: diagonalmatrix.hh:800
ConstIterator ConstRowIterator
rename the iterators for easier access
Definition: diagonalmatrix.hh:183
static constexpr int rows
The number of rows.
Definition: diagonalmatrix.hh:81
static constexpr size_type M()
number of blocks in column direction
Definition: diagonalmatrix.hh:488
FieldVector< K, n > & diagonal()
Get reference to diagonal vector.
Definition: diagonalmatrix.hh:549
void usmhv(const typename FieldTraits< Y >::field_type &alpha, const X &x, Y &y) const
y += alpha A^H x
Definition: diagonalmatrix.hh:395
Iterator iterator
typedef for stl compliant access
Definition: diagonalmatrix.hh:888
static constexpr int blocklevel
The number of block levels we contain. This is the leaf, that is, 1.
Definition: diagonalmatrix.hh:70
K field_type
export the type representing the field
Definition: diagonalmatrix.hh:712
ConstIterator beforeEnd() const
Definition: diagonalmatrix.hh:773
const_row_type::ConstIterator ConstColIterator
rename the iterators for easier access
Definition: diagonalmatrix.hh:185
bool exists(size_type i, size_type j) const
return true when (i,j) is in pattern
Definition: diagonalmatrix.hh:498
ContainerWrapperIterator< const WrapperType, const_reference, const_reference > ConstIterator
Iterator class for sequential access.
Definition: diagonalmatrix.hh:179
DiagonalRowVector(K *p, int col)
Constructor making vector with identical coordinates.
Definition: diagonalmatrix.hh:865
K & diagonal(size_type i)
Get reference to diagonal entry.
Definition: diagonalmatrix.hh:537
void solve(V &x, const V &b) const
Solve system A x = b.
Definition: diagonalmatrix.hh:438
static constexpr size_type N()
number of blocks in row direction
Definition: diagonalmatrix.hh:482
Iterator beforeBegin()
Definition: diagonalmatrix.hh:911
const_reference operator[](size_type i) const
Return const_reference object as row replacement.
Definition: diagonalmatrix.hh:525
Iterator iterator
typedef for stl compliant access
Definition: diagonalmatrix.hh:145
ConstIterator begin() const
begin ConstIterator
Definition: diagonalmatrix.hh:760
ConstIterator const_iterator
typedef for stl compliant access
Definition: diagonalmatrix.hh:919
DiagonalMatrix & operator-=(const DiagonalMatrix &y)
vector space subtraction
Definition: diagonalmatrix.hh:225
DiagonalRowVectorConst(K *p, int col)
Constructor making vector with identical coordinates.
Definition: diagonalmatrix.hh:733
void mmtv(const X &x, Y &y) const
y -= A^T x
Definition: diagonalmatrix.hh:345
DiagonalMatrix(const K &k)
Constructor initializing the whole matrix with a scalar.
Definition: diagonalmatrix.hh:98
ContainerWrapperIterator< DiagonalRowVector< K, n >, K, K & > Iterator
Iterator class for sequential access.
Definition: diagonalmatrix.hh:886
Iterator beforeEnd()
Definition: diagonalmatrix.hh:904
void umtv(const X &x, Y &y) const
y += A^T x
Definition: diagonalmatrix.hh:309
ConstIterator const_iterator
typedef for stl compliant access
Definition: diagonalmatrix.hh:757
double infinity_norm_real() const
simplified infinity norm (uses Manhattan norm for complex values)
Definition: diagonalmatrix.hh:427
void umv(const X &x, Y &y) const
y += A x
Definition: diagonalmatrix.hh:297
ContainerWrapperIterator< const WrapperType, reference, reference > Iterator
Iterator class for sequential access.
Definition: diagonalmatrix.hh:143
void mv(const X &x, Y &y) const
y = A x
Definition: diagonalmatrix.hh:278
double frobenius_norm() const
frobenius norm: sqrt(sum over squared values of entries)
Definition: diagonalmatrix.hh:409
ConstIterator end() const
end iterator
Definition: diagonalmatrix.hh:194
static constexpr int cols
The number of columns.
Definition: diagonalmatrix.hh:83
size_type rowIndex() const
index of this row in surrounding matrix
Definition: diagonalmatrix.hh:806
bool operator!=(const DiagonalMatrix &other) const
incomparison operator
Definition: diagonalmatrix.hh:268
ConstIterator begin() const
begin iterator
Definition: diagonalmatrix.hh:188
Iterator begin()
begin iterator
Definition: diagonalmatrix.hh:891
void mmv(const X &x, Y &y) const
y -= A x
Definition: diagonalmatrix.hh:333
DiagonalMatrix & operator/=(const K &k)
vector space division by scalar
Definition: diagonalmatrix.hh:253
ConstIterator const_iterator
typedef for stl compliant access
Definition: diagonalmatrix.hh:181
bool identical(const DiagonalMatrix< K, n > &other) const
Check if matrix is the same object as the other matrix.
Definition: diagonalmatrix.hh:130
friend auto operator*(const DiagonalMatrix &matrixA, const DiagonalMatrix< OtherScalar, n > &matrixB)
Matrix-matrix multiplication.
Definition: diagonalmatrix.hh:468
Iterator end()
end iterator
Definition: diagonalmatrix.hh:158
DiagonalRowVector()
Constructor making uninitialized vector.
Definition: diagonalmatrix.hh:861
double frobenius_norm2() const
square of frobenius norm, need for block recursion
Definition: diagonalmatrix.hh:415
constexpr DiagonalMatrix()=default
Default constructor.
void mtv(const X &x, Y &y) const
y = A^T x
Definition: diagonalmatrix.hh:290
ContainerWrapperIterator< DiagonalRowVectorConst< K, n >, const K, const K & > ConstIterator
ConstIterator class for sequential access.
Definition: diagonalmatrix.hh:917
void usmtv(const typename FieldTraits< Y >::field_type &alpha, const X &x, Y &y) const
y += alpha A^T x
Definition: diagonalmatrix.hh:382
const FieldVector< K, n > & diagonal() const
Get const reference to diagonal vector.
Definition: diagonalmatrix.hh:543
void invert()
Compute inverse.
Definition: diagonalmatrix.hh:445
const K & diagonal(size_type i) const
Get const reference to diagonal entry.
Definition: diagonalmatrix.hh:531
size_type N() const
number of blocks in the vector (are of size 1 here)
Definition: diagonalmatrix.hh:794
DiagonalMatrix(std::initializer_list< K > const &l)
Construct diagonal matrix from an initializer list.
Definition: diagonalmatrix.hh:115
row_type::Iterator ColIterator
rename the iterators for easier access
Definition: diagonalmatrix.hh:149
K field_type
export the type representing the field
Definition: diagonalmatrix.hh:852
reference operator[](size_type i)
Return reference object as row replacement.
Definition: diagonalmatrix.hh:519
Iterator end()
end iterator
Definition: diagonalmatrix.hh:897
DiagonalMatrix(const FieldVector< K, n > &diag)
Constructor initializing the diagonal with a vector.
Definition: diagonalmatrix.hh:103
std::size_t size_type
The type used for the index access and size operation.
Definition: diagonalmatrix.hh:858
void usmv(const typename FieldTraits< Y >::field_type &alpha, const X &x, Y &y) const
y += alpha A x
Definition: diagonalmatrix.hh:369
double infinity_norm() const
infinity norm (row sum norm, how to generalize for blocks?)
Definition: diagonalmatrix.hh:421
DiagonalMatrix & operator+=(const DiagonalMatrix &y)
vector space addition
Definition: diagonalmatrix.hh:218
K block_type
export the type representing the components
Definition: diagonalmatrix.hh:64
void umhv(const X &x, Y &y) const
y += A^H x
Definition: diagonalmatrix.hh:321
Iterator beforeBegin()
Definition: diagonalmatrix.hh:172
ConstIterator end() const
end ConstIterator
Definition: diagonalmatrix.hh:766
std::size_t size_type
The type used for the index access and size operation.
Definition: diagonalmatrix.hh:718
DiagonalMatrix & operator=(const K &k)
Assignment from a scalar.
Definition: diagonalmatrix.hh:123
DiagonalRowVectorConst()
Constructor making uninitialized vector.
Definition: diagonalmatrix.hh:727
Iterator beforeEnd()
Definition: diagonalmatrix.hh:165
friend std::ostream & operator<<(std::ostream &s, const DiagonalMatrix< K, n > &a)
Sends the matrix to an output stream.
Definition: diagonalmatrix.hh:508
K value_type
export the type representing the field
Definition: diagonalmatrix.hh:60
K block_type
export the type representing the components
Definition: diagonalmatrix.hh:715
K determinant() const
calculates the determinant of this matrix
Definition: diagonalmatrix.hh:453
ConstIterator beforeBegin() const
Definition: diagonalmatrix.hh:780
ContainerWrapperIterator< DiagonalRowVectorConst< K, n >, const K, const K & > ConstIterator
ConstIterator class for sequential access.
Definition: diagonalmatrix.hh:755
Iterator begin()
begin iterator
Definition: diagonalmatrix.hh:152
bool operator==(const DiagonalMatrix &other) const
comparison operator
Definition: diagonalmatrix.hh:262
DiagonalRowVector< K, n > row_type
Each row is implemented by a field vector.
Definition: diagonalmatrix.hh:73
Iterator RowIterator
rename the iterators for easier access
Definition: diagonalmatrix.hh:147
K block_type
export the type representing the components
Definition: diagonalmatrix.hh:855
ConstIterator beforeEnd() const
Definition: diagonalmatrix.hh:201
const K & diagonal() const
the diagonal value
Definition: diagonalmatrix.hh:812
#define DUNE_THROW(E, m)
Definition: exceptions.hh:218
constexpr auto min
Function object that returns the smaller of the given values.
Definition: hybridutilities.hh:506
EnableIfInterOperable< T1, T2, bool >::type operator==(const ForwardIteratorFacade< T1, V1, R1, D > &lhs, const ForwardIteratorFacade< T2, V2, R2, D > &rhs)
Checks for equality.
Definition: iteratorfacades.hh:238
Dune namespace.
Definition: alignedallocator.hh:13
constexpr std::integral_constant< std::size_t, sizeof...(II)> size(std::integer_sequence< T, II... >)
Return the size of the sequence.
Definition: integersequence.hh:75
K conjugateComplex(const K &x)
compute conjugate complex of x
Definition: math.hh:164
you have to specialize this structure for any type that should be assignable to a DenseMatrix
Definition: densematrix.hh:59
get the 'mutable' version of a reference to a const object
Definition: genericiterator.hh:116
Traits for type conversions and type information.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 15, 22:36, 2024)