Dune Core Modules (2.6.0)

matrix.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_ISTL_MATRIX_HH
4#define DUNE_ISTL_MATRIX_HH
5
10#include <memory>
11#include <cmath>
12
14
15#include <dune/istl/istlexception.hh>
16#include <dune/istl/bvector.hh>
17
18namespace Dune {
19
20namespace MatrixImp
21{
33 template<class B, class A=std::allocator<B> >
34 class DenseMatrixBase : public Imp::block_vector_unmanaged<B,A>
35 // this derivation gives us all the blas level 1 and norms
36 // on the large array. However, access operators have to be
37 // overwritten.
38 {
39 public:
40
41 //===== type definitions and constants
42
44 typedef typename B::field_type field_type;
45
47 typedef A allocator_type;
48
50 typedef typename A::size_type size_type;
51
58
62
63 // just a shorthand
64 typedef Imp::BlockVectorWindow<B,A> window_type;
65
66 typedef window_type reference;
67
68 typedef const window_type const_reference;
69
70
71 //===== constructors and such
72
76 DenseMatrixBase () : Imp::block_vector_unmanaged<B,A>()
77 {
78 // nothing is known ...
79 rows_ = 0;
80 columns_ = 0;
81 }
82
89 DenseMatrixBase (size_type rows, size_type columns) : Imp::block_vector_unmanaged<B,A>()
90 {
91 // and we can allocate the big array in the base class
92 this->n = rows*columns;
93 columns_ = columns;
94 if (this->n>0)
95 {
96 this->p = allocator_.allocate(this->n);
97 new (this->p)B[this->n];
98 }
99 else
100 {
101 this->n = 0;
102 this->p = 0;
103 }
104
105 // we can allocate the windows now
106 rows_ = rows;
107 }
108
111 {
112 // allocate the big array in the base class
113 this->n = a.n;
114 columns_ = a.columns_;
115 if (this->n>0)
116 {
117 // allocate and construct objects
118 this->p = allocator_.allocate(this->n);
119 new (this->p)B[this->n];
120
121 // copy data
122 for (size_type i=0; i<this->n; i++)
123 this->p[i]=a.p[i];
124 }
125 else
126 {
127 this->n = 0;
128 this->p = nullptr;
129 }
130
131 // we can allocate the windows now
132 rows_ = a.rows_;
133 }
134
137 {
138 if (this->n>0) {
139 size_type i=this->n;
140 while (i)
141 this->p[--i].~B();
142 allocator_.deallocate(this->p,this->n);
143 }
144 }
145
147 void resize (size_type rows, size_type columns)
148 {
149 // deconstruct objects and deallocate memory if necessary
150 if (this->n>0) {
151 size_type i=this->n;
152 while (i)
153 this->p[--i].~B();
154 allocator_.deallocate(this->p,this->n);
155 }
156
157 // and we can allocate the big array in the base class
158 this->n = rows*columns;
159 if (this->n>0)
160 {
161 this->p = allocator_.allocate(this->n);
162 new (this->p)B[this->n];
163 }
164 else
165 {
166 this->n = 0;
167 this->p = nullptr;
168 }
169
170 // we can allocate the windows now
171 rows_ = rows;
172 columns_ = columns;
173 }
174
177 {
178 if (&a!=this) // check if this and a are different objects
179 {
180 columns_ = a.columns_;
181 // reallocate arrays if necessary
182 // Note: still the block sizes may vary !
183 if (this->n!=a.n || rows_!=a.rows_)
184 {
185 // deconstruct objects and deallocate memory if necessary
186 if (this->n>0) {
187 size_type i=this->n;
188 while (i)
189 this->p[--i].~B();
190 allocator_.deallocate(this->p,this->n);
191 }
192
193 // allocate the big array in the base class
194 this->n = a.n;
195 if (this->n>0)
196 {
197 // allocate and construct objects
198 this->p = allocator_.allocate(this->n);
199 new (this->p)B[this->n];
200 }
201 else
202 {
203 this->n = 0;
204 this->p = nullptr;
205 }
206
207 // Copy number of rows
208 rows_ = a.rows_;
209 }
210
211 // and copy the data
212 for (size_type i=0; i<this->n; i++)
213 this->p[i]=a.p[i];
214 }
215
216 return *this;
217 }
218
219
220 //===== assignment from scalar
221
224 {
225 (static_cast<Imp::block_vector_unmanaged<B,A>&>(*this)) = k;
226 return *this;
227 }
228
229
230 //===== access to components
231 // has to be overwritten from base class because it must
232 // return access to the windows
233
236 {
237#ifdef DUNE_ISTL_WITH_CHECKING
238 if (i>=rows_) DUNE_THROW(ISTLError,"index out of range");
239#endif
240 return window_type(this->p + i*columns_, columns_);
241 }
242
244 const_reference operator[] (size_type i) const
245 {
246#ifdef DUNE_ISTL_WITH_CHECKING
247 if (i<0 || i>=rows_) DUNE_THROW(ISTLError,"index out of range");
248#endif
249 return window_type(this->p + i*columns_, columns_);
250 }
251
252 // forward declaration
253 class ConstIterator;
254
257 {
258 public:
261 : window_(nullptr,0)
262 {
263 i = 0;
264 }
265
266 Iterator (Iterator& other) = default;
267 Iterator (Iterator&& other) = default;
268
270 Iterator (B* data, size_type columns, size_type _i)
271 : i(_i),
272 window_(data + _i*columns, columns)
273 {}
274
277 {
278 i = other.i;
279 // Do NOT use window_.operator=, because that copies the window content, not just the window!
280 window_.set(other.window_.getsize(),other.window_.getptr());
281 return *this;
282 }
283
286 {
287 i = other.i;
288 // Do NOT use window_.operator=, because that copies the window content, not just the window!
289 window_.set(other.window_.getsize(),other.window_.getptr());
290 return *this;
291 }
292
295 {
296 ++i;
297 window_.setptr(window_.getptr()+window_.getsize());
298 return *this;
299 }
300
303 {
304 --i;
305 window_.setptr(window_.getptr()-window_.getsize());
306 return *this;
307 }
308
310 bool operator== (const Iterator& it) const
311 {
312 return window_.getptr() == it.window_.getptr();
313 }
314
316 bool operator!= (const Iterator& it) const
317 {
318 return window_.getptr() != it.window_.getptr();
319 }
320
322 bool operator== (const ConstIterator& it) const
323 {
324 return window_.getptr() == it.window_.getptr();
325 }
326
328 bool operator!= (const ConstIterator& it) const
329 {
330 return window_.getptr() != it.window_.getptr();
331 }
332
334 window_type& operator* () const
335 {
336 return window_;
337 }
338
340 window_type* operator-> () const
341 {
342 return &window_;
343 }
344
345 // return index corresponding to pointer
346 size_type index () const
347 {
348 return i;
349 }
350
351 friend class ConstIterator;
352
353 private:
354 size_type i;
355 mutable window_type window_;
356 };
357
360 {
361 return Iterator(this->p, columns_, 0);
362 }
363
366 {
367 return Iterator(this->p, columns_, rows_);
368 }
369
373 {
374 return Iterator(this->p, columns_, rows_-1);
375 }
376
380 {
381 return Iterator(this->p, columns_, -1);
382 }
383
386 {
387 return Iterator(this->p, columns_, std::min(i,rows_));
388 }
389
392 {
393 return ConstIterator(this->p, columns_, std::min(i,rows_));
394 }
395
398 {
399 public:
402 : window_(nullptr,0)
403 {
404 i = 0;
405 }
406
408 ConstIterator (const B* data, size_type columns, size_type _i)
409 : i(_i),
410 window_(const_cast<B*>(data + _i * columns), columns)
411 {}
412
415 : i(it.i), window_(it.window_.getptr(),it.window_.getsize())
416 {}
417
418 ConstIterator& operator=(Iterator&& other)
419 {
420 i = other.i;
421 // Do NOT use window_.operator=, because that copies the window content, not just the window!
422 window_.set(other.window_.getsize(),other.window_.getptr());
423 return *this;
424 }
425
426 ConstIterator& operator=(Iterator& other)
427 {
428 i = other.i;
429 // Do NOT use window_.operator=, because that copies the window content, not just the window!
430 window_.set(other.window_.getsize(),other.window_.getptr());
431 return *this;
432 }
433
436 {
437 ++i;
438 window_.setptr(window_.getptr()+window_.getsize());
439 return *this;
440 }
441
444 {
445 --i;
446 window_.setptr(window_.getptr()-window_.getsize());
447 return *this;
448 }
449
451 bool operator== (const ConstIterator& it) const
452 {
453 return window_.getptr() == it.window_.getptr();
454 }
455
457 bool operator!= (const ConstIterator& it) const
458 {
459 return window_.getptr() != it.window_.getptr();
460 }
461
463 bool operator== (const Iterator& it) const
464 {
465 return window_.getptr() == it.window_.getptr();
466 }
467
469 bool operator!= (const Iterator& it) const
470 {
471 return window_.getptr() != it.window_.getptr();
472 }
473
475 const window_type& operator* () const
476 {
477 return window_;
478 }
479
481 const window_type* operator-> () const
482 {
483 return &window_;
484 }
485
486 // return index corresponding to pointer
487 size_type index () const
488 {
489 return i;
490 }
491
492 friend class Iterator;
493
494 private:
495 size_type i;
496 mutable window_type window_;
497 };
498
501
504
507 {
508 return ConstIterator(this->p, columns_, 0);
509 }
510
513 {
514 return ConstIterator(this->p, columns_, rows_);
515 }
516
520 {
521 return ConstIterator(this->p, columns_, rows_-1);
522 }
523
526 {
527 return ConstIterator(this->p, columns_, -1);
528 }
529
530 //===== sizes
531
533 size_type N () const
534 {
535 return rows_;
536 }
537
538
539 private:
540 size_type rows_; // number of matrix rows
541 size_type columns_; // number of matrix columns
542
543 A allocator_;
544 };
545
546} // namespace MatrixImp
547
553 template<class T, class A=std::allocator<T> >
554 class Matrix
555 {
556 public:
557
559 typedef typename T::field_type field_type;
560
562 typedef T block_type;
563
565 typedef A allocator_type;
566
568 typedef typename MatrixImp::DenseMatrixBase<T,A>::window_type row_type;
569
571 typedef typename A::size_type size_type;
572
575
577 typedef typename row_type::iterator ColIterator;
578
581
583 typedef typename row_type::const_iterator ConstColIterator;
584
585 enum {
587 blocklevel = T::blocklevel+1
588 };
589
591 Matrix() : data_(0,0), cols_(0)
592 {}
593
596 Matrix(size_type rows, size_type cols) : data_(rows,cols), cols_(cols)
597 {}
598
603 void setSize(size_type rows, size_type cols) {
604 data_.resize(rows,cols);
605 cols_ = cols;
606 }
607
610 {
611 return data_.begin();
612 }
613
616 {
617 return data_.end();
618 }
619
623 {
624 return data_.beforeEnd();
625 }
626
630 {
631 return data_.beforeBegin();
632 }
633
636 {
637 return data_.begin();
638 }
639
642 {
643 return data_.end();
644 }
645
649 {
650 return data_.beforeEnd();
651 }
652
656 {
657 return data_.beforeBegin();
658 }
659
662 {
663 data_ = t;
664 return *this;
665 }
666
669#ifdef DUNE_ISTL_WITH_CHECKING
670 if (row<0)
671 DUNE_THROW(ISTLError, "Can't access negative rows!");
672 if (row>=N())
673 DUNE_THROW(ISTLError, "Row index out of range!");
674#endif
675 return data_[row];
676 }
677
679 const row_type operator[](size_type row) const {
680#ifdef DUNE_ISTL_WITH_CHECKING
681 if (row<0)
682 DUNE_THROW(ISTLError, "Can't access negative rows!");
683 if (row>=N())
684 DUNE_THROW(ISTLError, "Row index out of range!");
685#endif
686 return data_[row];
687 }
688
690 size_type N() const {
691 return data_.N();
692 }
693
695 size_type M() const {
696 return cols_;
697 }
698
701 data_ *= scalar;
702 return (*this);
703 }
704
707 data_ /= scalar;
708 return (*this);
709 }
710
717#ifdef DUNE_ISTL_WITH_CHECKING
718 if(N()!=b.N() || M() != b.M())
719 DUNE_THROW(RangeError, "Matrix sizes do not match!");
720#endif
721 data_ += b.data_;
722 return (*this);
723 }
724
731#ifdef DUNE_ISTL_WITH_CHECKING
732 if(N()!=b.N() || M() != b.M())
733 DUNE_THROW(RangeError, "Matrix sizes do not match!");
734#endif
735 data_ -= b.data_;
736 return (*this);
737 }
738
741 Matrix out(M(), N());
742 for (size_type i=0; i<N(); i++)
743 for (size_type j=0; j<M(); j++)
744 out[j][i] = (*this)[i][j];
745
746 return out;
747 }
748
750 friend Matrix<T> operator*(const Matrix<T>& m1, const Matrix<T>& m2) {
751 Matrix<T> out(m1.N(), m2.M());
752 out = 0;
753
754 for (size_type i=0; i<out.N(); i++ ) {
755 for ( size_type j=0; j<out.M(); j++ )
756 for (size_type k=0; k<m1.M(); k++)
757 out[i][j] += m1[i][k]*m2[k][j];
758 }
759
760 return out;
761 }
762
764 template <class X, class Y>
765 friend Y operator*(const Matrix<T>& m, const X& vec) {
766#ifdef DUNE_ISTL_WITH_CHECKING
767 if (m.M()!=vec.size())
768 DUNE_THROW(ISTLError, "Vector size doesn't match the number of matrix columns!");
769#endif
770 Y out(m.N());
771 out = 0;
772
773 for (size_type i=0; i<out.size(); i++ ) {
774 for ( size_type j=0; j<vec.size(); j++ )
775 out[i] += m[i][j]*vec[j];
776 }
777
778 return out;
779 }
780
782 template <class X, class Y>
783 void mv(const X& x, Y& y) const
784 {
785#ifdef DUNE_ISTL_WITH_CHECKING
786 if (x.N()!=M()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
787 if (y.N()!=N()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
788#endif
789
790 for (size_type i=0; i<data_.N(); i++) {
791 y[i]=0;
792 for (size_type j=0; j<cols_; j++)
793 (*this)[i][j].umv(x[j], y[i]);
794
795 }
796
797 }
798
800 template<class X, class Y>
801 void mtv (const X& x, Y& y) const
802 {
803#ifdef DUNE_ISTL_WITH_CHECKING
804 if (x.N()!=N()) DUNE_THROW(ISTLError,"index out of range");
805 if (y.N()!=M()) DUNE_THROW(ISTLError,"index out of range");
806#endif
807
808 for(size_type i=0; i<y.N(); ++i)
809 y[i]=0;
810 umtv(x,y);
811 }
812
814 template <class X, class Y>
815 void umv(const X& x, Y& y) const
816 {
817#ifdef DUNE_ISTL_WITH_CHECKING
818 if (x.N()!=M()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
819 if (y.N()!=N()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
820#endif
821
822 for (size_type i=0; i<data_.N(); i++) {
823
824 for (size_type j=0; j<cols_; j++)
825 (*this)[i][j].umv(x[j], y[i]);
826
827 }
828
829 }
830
832 template<class X, class Y>
833 void mmv (const X& x, Y& y) const
834 {
835#ifdef DUNE_ISTL_WITH_CHECKING
836 if (x.N()!=M()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
837 if (y.N()!=N()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
838#endif
839 ConstRowIterator endi=end();
840 for (ConstRowIterator i=begin(); i!=endi; ++i)
841 {
842 ConstColIterator endj = (*i).end();
843 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
844 (*j).mmv(x[j.index()],y[i.index()]);
845 }
846 }
847
849 template <class X, class Y>
850 void usmv(const field_type& alpha, const X& x, Y& y) const
851 {
852#ifdef DUNE_ISTL_WITH_CHECKING
853 if (x.N()!=M()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
854 if (y.N()!=N()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
855#endif
856
857 for (size_type i=0; i<data_.N(); i++) {
858
859 for (size_type j=0; j<cols_; j++)
860 (*this)[i][j].usmv(alpha, x[j], y[i]);
861
862 }
863
864 }
865
867 template<class X, class Y>
868 void umtv (const X& x, Y& y) const
869 {
870#ifdef DUNE_ISTL_WITH_CHECKING
871 if (x.N()!=N()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
872 if (y.N()!=M()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
873#endif
874 ConstRowIterator endi=end();
875 for (ConstRowIterator i=begin(); i!=endi; ++i)
876 {
877 ConstColIterator endj = (*i).end();
878 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
879 (*j).umtv(x[i.index()],y[j.index()]);
880 }
881 }
882
884 template<class X, class Y>
885 void mmtv (const X& x, Y& y) const
886 {
887#ifdef DUNE_ISTL_WITH_CHECKING
888 if (x.N()!=N()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
889 if (y.N()!=M()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
890#endif
891 ConstRowIterator endi=end();
892 for (ConstRowIterator i=begin(); i!=endi; ++i)
893 {
894 ConstColIterator endj = (*i).end();
895 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
896 (*j).mmtv(x[i.index()],y[j.index()]);
897 }
898 }
899
901 template<class X, class Y>
902 void usmtv (const field_type& alpha, const X& x, Y& y) const
903 {
904#ifdef DUNE_ISTL_WITH_CHECKING
905 if (x.N()!=N()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
906 if (y.N()!=M()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
907#endif
908 ConstRowIterator endi=end();
909 for (ConstRowIterator i=begin(); i!=endi; ++i)
910 {
911 ConstColIterator endj = (*i).end();
912 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
913 (*j).usmtv(alpha,x[i.index()],y[j.index()]);
914 }
915 }
916
918 template<class X, class Y>
919 void umhv (const X& x, Y& y) const
920 {
921#ifdef DUNE_ISTL_WITH_CHECKING
922 if (x.N()!=N()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
923 if (y.N()!=M()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
924#endif
925 ConstRowIterator endi=end();
926 for (ConstRowIterator i=begin(); i!=endi; ++i)
927 {
928 ConstColIterator endj = (*i).end();
929 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
930 (*j).umhv(x[i.index()],y[j.index()]);
931 }
932 }
933
935 template<class X, class Y>
936 void mmhv (const X& x, Y& y) const
937 {
938#ifdef DUNE_ISTL_WITH_CHECKING
939 if (x.N()!=N()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
940 if (y.N()!=M()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
941#endif
942 ConstRowIterator endi=end();
943 for (ConstRowIterator i=begin(); i!=endi; ++i)
944 {
945 ConstColIterator endj = (*i).end();
946 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
947 (*j).mmhv(x[i.index()],y[j.index()]);
948 }
949 }
950
952 template<class X, class Y>
953 void usmhv (const field_type& alpha, const X& x, Y& y) const
954 {
955#ifdef DUNE_ISTL_WITH_CHECKING
956 if (x.N()!=N()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
957 if (y.N()!=M()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
958#endif
959 ConstRowIterator endi=end();
960 for (ConstRowIterator i=begin(); i!=endi; ++i)
961 {
962 ConstColIterator endj = (*i).end();
963 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
964 (*j).usmhv(alpha,x[i.index()],y[j.index()]);
965 }
966 }
967
968 //===== norms
969
971 typename FieldTraits<field_type>::real_type frobenius_norm () const
972 {
973 return std::sqrt(frobenius_norm2());
974 }
975
977 typename FieldTraits<field_type>::real_type frobenius_norm2 () const
978 {
979 double sum=0;
980 for (size_type i=0; i<N(); ++i)
981 for (size_type j=0; j<M(); ++j)
982 sum += data_[i][j].frobenius_norm2();
983 return sum;
984 }
985
987 template <typename ft = field_type,
988 typename std::enable_if<!has_nan<ft>::value, int>::type = 0>
989 typename FieldTraits<ft>::real_type infinity_norm() const {
990 using real_type = typename FieldTraits<ft>::real_type;
991 using std::max;
992
993 real_type norm = 0;
994 for (auto const &x : *this) {
995 real_type sum = 0;
996 for (auto const &y : x)
997 sum += y.infinity_norm();
998 norm = max(sum, norm);
999 }
1000 return norm;
1001 }
1002
1004 template <typename ft = field_type,
1005 typename std::enable_if<!has_nan<ft>::value, int>::type = 0>
1006 typename FieldTraits<ft>::real_type infinity_norm_real() const {
1007 using real_type = typename FieldTraits<ft>::real_type;
1008 using std::max;
1009
1010 real_type norm = 0;
1011 for (auto const &x : *this) {
1012 real_type sum = 0;
1013 for (auto const &y : x)
1014 sum += y.infinity_norm_real();
1015 norm = max(sum, norm);
1016 }
1017 return norm;
1018 }
1019
1021 template <typename ft = field_type,
1022 typename std::enable_if<has_nan<ft>::value, int>::type = 0>
1023 typename FieldTraits<ft>::real_type infinity_norm() const {
1024 using real_type = typename FieldTraits<ft>::real_type;
1025 using std::max;
1026
1027 real_type norm = 0;
1028 real_type isNaN = 1;
1029 for (auto const &x : *this) {
1030 real_type sum = 0;
1031 for (auto const &y : x)
1032 sum += y.infinity_norm();
1033 norm = max(sum, norm);
1034 isNaN += sum;
1035 }
1036 isNaN /= isNaN;
1037 return norm * isNaN;
1038 }
1039
1041 template <typename ft = field_type,
1042 typename std::enable_if<has_nan<ft>::value, int>::type = 0>
1043 typename FieldTraits<ft>::real_type infinity_norm_real() const {
1044 using real_type = typename FieldTraits<ft>::real_type;
1045 using std::max;
1046
1047 real_type norm = 0;
1048 real_type isNaN = 1;
1049 for (auto const &x : *this) {
1050 real_type sum = 0;
1051 for (auto const &y : x)
1052 sum += y.infinity_norm_real();
1053 norm = max(sum, norm);
1054 isNaN += sum;
1055 }
1056 isNaN /= isNaN;
1057 return norm * isNaN;
1058 }
1059
1060 //===== query
1061
1063 bool exists (size_type i, size_type j) const
1064 {
1065#ifdef DUNE_ISTL_WITH_CHECKING
1066 if (i<0 || i>=N()) DUNE_THROW(ISTLError,"row index out of range");
1067 if (j<0 || i>=M()) DUNE_THROW(ISTLError,"column index out of range");
1068#else
1070#endif
1071 return true;
1072 }
1073
1074 protected:
1075
1079
1086 };
1087
1089} // end namespace Dune
1090
1091#endif
This file implements a vector space as a tensor product of a given vector space. The number of compon...
A vector of blocks with memory management.
Definition: bvector.hh:317
derive error class from the base class in common
Definition: istlexception.hh:16
ConstIterator class for sequential access.
Definition: matrix.hh:398
const window_type * operator->() const
arrow
Definition: matrix.hh:481
const window_type & operator*() const
dereferencing
Definition: matrix.hh:475
ConstIterator & operator++()
prefix increment
Definition: matrix.hh:435
ConstIterator(const B *data, size_type columns, size_type _i)
constructor from pointer
Definition: matrix.hh:408
ConstIterator & operator--()
prefix decrement
Definition: matrix.hh:443
ConstIterator(const Iterator &it)
constructor from non_const iterator
Definition: matrix.hh:414
bool operator!=(const ConstIterator &it) const
inequality
Definition: matrix.hh:457
ConstIterator()
constructor
Definition: matrix.hh:401
bool operator==(const ConstIterator &it) const
equality
Definition: matrix.hh:451
Iterator class for sequential access.
Definition: matrix.hh:257
Iterator & operator--()
prefix decrement
Definition: matrix.hh:302
bool operator!=(const Iterator &it) const
inequality
Definition: matrix.hh:316
Iterator & operator=(Iterator &&other)
Move assignment.
Definition: matrix.hh:276
Iterator & operator++()
prefix increment
Definition: matrix.hh:294
Iterator()
constructor, no arguments
Definition: matrix.hh:260
window_type & operator*() const
dereferencing
Definition: matrix.hh:334
bool operator==(const Iterator &it) const
equality
Definition: matrix.hh:310
Iterator & operator=(Iterator &other)
Copy assignment.
Definition: matrix.hh:285
Iterator(B *data, size_type columns, size_type _i)
constructor
Definition: matrix.hh:270
window_type * operator->() const
arrow
Definition: matrix.hh:340
A Vector of blocks with different blocksizes.
Definition: matrix.hh:38
BlockVector< B, A > block_type
Same as value_type, here for historical reasons.
Definition: matrix.hh:61
DenseMatrixBase & operator=(const DenseMatrixBase &a)
assignment
Definition: matrix.hh:176
Iterator beforeBegin() const
Definition: matrix.hh:379
void resize(size_type rows, size_type columns)
same effect as constructor with same argument
Definition: matrix.hh:147
DenseMatrixBase(size_type rows, size_type columns)
Definition: matrix.hh:89
Iterator end()
end Iterator
Definition: matrix.hh:365
reference operator[](size_type i)
random access to blocks
Definition: matrix.hh:235
BlockVector< B, A > value_type
Type of the elements of the outer vector, i.e., dynamic vectors of B.
Definition: matrix.hh:57
Iterator find(size_type i)
random access returning iterator (end if not contained)
Definition: matrix.hh:385
size_type N() const
number of blocks in the vector (are of variable size here)
Definition: matrix.hh:533
ConstIterator beforeEnd() const
Definition: matrix.hh:519
ConstIterator end() const
end ConstIterator
Definition: matrix.hh:512
ConstIterator rend() const
end ConstIterator
Definition: matrix.hh:525
Iterator beforeEnd()
Definition: matrix.hh:372
ConstIterator begin() const
begin ConstIterator
Definition: matrix.hh:506
A allocator_type
export the allocator type
Definition: matrix.hh:47
DenseMatrixBase()
Definition: matrix.hh:76
ConstIterator find(size_type i) const
random access returning iterator (end if not contained)
Definition: matrix.hh:391
A::size_type size_type
The size type for the index access.
Definition: matrix.hh:50
~DenseMatrixBase()
free dynamic memory
Definition: matrix.hh:136
B::field_type field_type
export the type representing the field
Definition: matrix.hh:44
Iterator begin()
begin Iterator
Definition: matrix.hh:359
DenseMatrixBase(const DenseMatrixBase &a)
copy constructor, has copy semantics
Definition: matrix.hh:110
A generic dynamic dense matrix.
Definition: matrix.hh:555
size_type cols_
Number of columns of the matrix.
Definition: matrix.hh:1085
FieldTraits< ft >::real_type infinity_norm() const
infinity norm (row sum norm, how to generalize for blocks?)
Definition: matrix.hh:989
A allocator_type
Export the allocator.
Definition: matrix.hh:565
FieldTraits< ft >::real_type infinity_norm_real() const
simplified infinity norm (uses Manhattan norm for complex values)
Definition: matrix.hh:1006
void usmhv(const field_type &alpha, const X &x, Y &y) const
y += alpha A^H x
Definition: matrix.hh:953
void usmv(const field_type &alpha, const X &x, Y &y) const
Definition: matrix.hh:850
A::size_type size_type
Type for indices and sizes.
Definition: matrix.hh:571
MatrixImp::DenseMatrixBase< T, A > data_
Abuse DenseMatrixBase as an engine for a 2d array ISTL-style.
Definition: matrix.hh:1078
Matrix transpose() const
Return the transpose of the matrix.
Definition: matrix.hh:740
void mtv(const X &x, Y &y) const
y = A^T x
Definition: matrix.hh:801
void umv(const X &x, Y &y) const
y += A x
Definition: matrix.hh:815
void mv(const X &x, Y &y) const
y = A x
Definition: matrix.hh:783
MatrixImp::DenseMatrixBase< T, A >::ConstIterator ConstRowIterator
Const iterator over the matrix rows.
Definition: matrix.hh:580
void setSize(size_type rows, size_type cols)
Change the matrix size.
Definition: matrix.hh:603
T::field_type field_type
Export the type representing the underlying field.
Definition: matrix.hh:559
RowIterator beforeBegin()
Definition: matrix.hh:629
RowIterator beforeEnd()
Definition: matrix.hh:622
Matrix()
Create empty matrix.
Definition: matrix.hh:591
Matrix & operator-=(const Matrix &b)
Subtract the entries of another matrix from this one.
Definition: matrix.hh:730
FieldTraits< field_type >::real_type frobenius_norm2() const
square of frobenius norm, need for block recursion
Definition: matrix.hh:977
row_type::iterator ColIterator
Iterator for the entries of each row.
Definition: matrix.hh:577
ConstRowIterator beforeEnd() const
Definition: matrix.hh:648
Matrix & operator=(const field_type &t)
Assignment from scalar.
Definition: matrix.hh:661
RowIterator end()
Get iterator to one beyond last row.
Definition: matrix.hh:615
const row_type operator[](size_type row) const
The const index operator.
Definition: matrix.hh:679
ConstRowIterator end() const
Get const iterator to one beyond last row.
Definition: matrix.hh:641
friend Y operator*(const Matrix< T > &m, const X &vec)
Generic matrix-vector multiplication.
Definition: matrix.hh:765
Matrix< T > & operator*=(const field_type &scalar)
Multiplication with a scalar.
Definition: matrix.hh:700
row_type operator[](size_type row)
The index operator.
Definition: matrix.hh:668
void mmv(const X &x, Y &y) const
y -= A x
Definition: matrix.hh:833
Matrix & operator+=(const Matrix &b)
Add the entries of another matrix to this one.
Definition: matrix.hh:716
ConstRowIterator begin() const
Get const iterator to first row.
Definition: matrix.hh:635
void mmhv(const X &x, Y &y) const
y -= A^H x
Definition: matrix.hh:936
RowIterator begin()
Get iterator to first row.
Definition: matrix.hh:609
row_type::const_iterator ConstColIterator
Const iterator for the entries of each row.
Definition: matrix.hh:583
@ blocklevel
The number of nesting levels the matrix contains.
Definition: matrix.hh:587
size_type M() const
Return the number of columns.
Definition: matrix.hh:695
T block_type
Export the type representing the components.
Definition: matrix.hh:562
bool exists(size_type i, size_type j) const
return true if (i,j) is in pattern
Definition: matrix.hh:1063
Matrix< T > & operator/=(const field_type &scalar)
Division by a scalar.
Definition: matrix.hh:706
friend Matrix< T > operator*(const Matrix< T > &m1, const Matrix< T > &m2)
Generic matrix multiplication.
Definition: matrix.hh:750
void umtv(const X &x, Y &y) const
y += A^T x
Definition: matrix.hh:868
void mmtv(const X &x, Y &y) const
y -= A^T x
Definition: matrix.hh:885
MatrixImp::DenseMatrixBase< T, A >::window_type row_type
The type implementing a matrix row.
Definition: matrix.hh:568
FieldTraits< field_type >::real_type frobenius_norm() const
frobenius norm: sqrt(sum over squared values of entries)
Definition: matrix.hh:971
void usmtv(const field_type &alpha, const X &x, Y &y) const
y += alpha A^T x
Definition: matrix.hh:902
void umhv(const X &x, Y &y) const
y += A^H x
Definition: matrix.hh:919
size_type N() const
Return the number of rows.
Definition: matrix.hh:690
ConstRowIterator beforeBegin() const
Definition: matrix.hh:655
Matrix(size_type rows, size_type cols)
Create uninitialized matrix of size rows x cols.
Definition: matrix.hh:596
MatrixImp::DenseMatrixBase< T, A >::Iterator RowIterator
Iterator over the matrix rows.
Definition: matrix.hh:574
Default exception class for range errors.
Definition: exceptions.hh:252
Type traits to determine the type of reals (when working with complex numbers)
#define DUNE_UNUSED_PARAMETER(parm)
A macro to mark intentionally unused function parameters with.
Definition: unused.hh:25
#define DUNE_THROW(E, m)
Definition: exceptions.hh:216
Dune namespace.
Definition: alignedallocator.hh:10
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 24, 23:30, 2024)