Dune Core Modules (2.7.0)

bcrsmatrix.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
4#ifndef DUNE_ISTL_BCRSMATRIX_HH
5#define DUNE_ISTL_BCRSMATRIX_HH
6
7#include <cmath>
8#include <complex>
9#include <set>
10#include <iostream>
11#include <algorithm>
12#include <numeric>
13#include <vector>
14#include <map>
15
16#include "istlexception.hh"
17#include "bvector.hh"
18#include "matrixutils.hh"
25
30namespace Dune {
31
71 template<typename M>
72 struct MatrixDimension;
73
75
81 template<typename size_type>
83 {
85 double avg;
87 size_type maximum;
89 size_type overflow_total;
91
94 double mem_ratio;
95 };
96
98
110 template<class M_>
112 {
113
114 public:
115
117 typedef M_ Matrix;
118
121
123 typedef typename Matrix::size_type size_type;
124
126
132 {
133
134 public:
135
138 {
139 return _m.entry(_i,j);
140 }
141
142#ifndef DOXYGEN
143
145 : _m(m)
146 , _i(i)
147 {}
148
149#endif
150
151 private:
152
153 Matrix& _m;
154 size_type _i;
155
156 };
157
159
166 : _m(m)
167 {
168 if (m.buildMode() != Matrix::implicit)
169 DUNE_THROW(BCRSMatrixError,"You can only create an ImplicitBuilder for a matrix in implicit build mode");
170 if (m.buildStage() != Matrix::building)
171 DUNE_THROW(BCRSMatrixError,"You can only create an ImplicitBuilder for a matrix with set size that has not been compressed() yet");
172 }
173
175
189 ImplicitMatrixBuilder(Matrix& m, size_type rows, size_type cols, size_type avg_cols_per_row, double overflow_fraction)
190 : _m(m)
191 {
192 if (m.buildStage() != Matrix::notAllocated)
193 DUNE_THROW(BCRSMatrixError,"You can only set up a matrix for this ImplicitBuilder if it has no memory allocated yet");
194 m.setBuildMode(Matrix::implicit);
195 m.setImplicitBuildModeParameters(avg_cols_per_row,overflow_fraction);
196 m.setSize(rows,cols);
197 }
198
201 {
202 return row_object(_m,i);
203 }
204
206 size_type N() const
207 {
208 return _m.N();
209 }
210
212 size_type M() const
213 {
214 return _m.M();
215 }
216
217 private:
218
219 Matrix& _m;
220
221 };
222
423 template<class B, class A=std::allocator<B> >
425 {
426 friend struct MatrixDimension<BCRSMatrix>;
427 public:
441 built=3
442 };
443
444 //===== type definitions and constants
445
447 using field_type = typename Imp::BlockTraits<B>::field_type;
448
450 typedef B block_type;
451
453 typedef A allocator_type;
454
456 typedef Imp::CompressedBlockVectorWindow<B,A> row_type;
457
459 typedef typename A::size_type size_type;
460
462 typedef ::Dune::CompressionStatistics<size_type> CompressionStatistics;
463
465 static constexpr unsigned int blocklevel = Imp::BlockTraits<B>::blockLevel()+1;
466
501 unknown
502 };
503
504 //===== random access interface to rows of the matrix
505
508 {
509#ifdef DUNE_ISTL_WITH_CHECKING
510 if (build_mode == implicit && ready != built)
511 DUNE_THROW(BCRSMatrixError,"You cannot use operator[] in implicit build mode before calling compress()");
512 if (r==0) DUNE_THROW(BCRSMatrixError,"row not initialized yet");
513 if (i>=n) DUNE_THROW(BCRSMatrixError,"index out of range");
514#endif
515 return r[i];
516 }
517
520 {
521#ifdef DUNE_ISTL_WITH_CHECKING
522 if (build_mode == implicit && ready != built)
523 DUNE_THROW(BCRSMatrixError,"You cannot use operator[] in implicit build mode before calling compress()");
524 if (built!=ready) DUNE_THROW(BCRSMatrixError,"row not initialized yet");
525 if (i>=n) DUNE_THROW(BCRSMatrixError,"index out of range");
526#endif
527 return r[i];
528 }
529
530
531 //===== iterator interface to rows of the matrix
532
534 template<class T>
536 : public RandomAccessIteratorFacade<RealRowIterator<T>, T>
537 {
538
539 public:
541 typedef typename std::remove_const<T>::type ValueType;
542
545 friend class RealRowIterator<const ValueType>;
546 friend class RealRowIterator<ValueType>;
547
550 : p(_p), i(_i)
551 {}
552
555 : p(0), i(0)
556 {}
557
559 : p(it.p), i(it.i)
560 {}
561
562
565 {
566 return i;
567 }
568
569 std::ptrdiff_t distanceTo(const RealRowIterator<ValueType>& other) const
570 {
571 assert(other.p==p);
572 return (other.i-i);
573 }
574
575 std::ptrdiff_t distanceTo(const RealRowIterator<const ValueType>& other) const
576 {
577 assert(other.p==p);
578 return (other.i-i);
579 }
580
582 bool equals (const RealRowIterator<ValueType>& other) const
583 {
584 assert(other.p==p);
585 return i==other.i;
586 }
587
589 bool equals (const RealRowIterator<const ValueType>& other) const
590 {
591 assert(other.p==p);
592 return i==other.i;
593 }
594
595 private:
597 void increment()
598 {
599 ++i;
600 }
601
603 void decrement()
604 {
605 --i;
606 }
607
608 void advance(std::ptrdiff_t diff)
609 {
610 i+=diff;
611 }
612
613 T& elementAt(std::ptrdiff_t diff) const
614 {
615 return p[i+diff];
616 }
617
619 row_type& dereference () const
620 {
621 return p[i];
622 }
623
624 row_type* p;
625 size_type i;
626 };
627
631
634 {
635 return Iterator(r,0);
636 }
637
640 {
641 return Iterator(r,n);
642 }
643
647 {
648 return Iterator(r,n-1);
649 }
650
654 {
655 return Iterator(r,-1);
656 }
657
660
662 typedef typename row_type::Iterator ColIterator;
663
667
668
671 {
672 return ConstIterator(r,0);
673 }
674
677 {
678 return ConstIterator(r,n);
679 }
680
684 {
685 return ConstIterator(r,n-1);
686 }
687
691 {
692 return ConstIterator(r,-1);
693 }
694
697
699 typedef typename row_type::ConstIterator ConstColIterator;
700
701 //===== constructors & resizers
702
703 // we use a negative overflowsize to indicate that the implicit
704 // mode parameters have not been set yet
705
708 : build_mode(unknown), ready(notAllocated), n(0), m(0), nnz_(0),
709 allocationSize_(0), r(0), a(0),
710 avg(0), overflowsize(-1.0)
711 {}
712
715 : build_mode(bm), ready(notAllocated), n(0), m(0), nnz_(0),
716 allocationSize_(0), r(0), a(0),
717 avg(0), overflowsize(-1.0)
718 {
719 allocate(_n, _m, _nnz,true,false);
720 }
721
724 : build_mode(bm), ready(notAllocated), n(0), m(0), nnz_(0),
725 allocationSize_(0), r(0), a(0),
726 avg(0), overflowsize(-1.0)
727 {
728 allocate(_n, _m,0,true,false);
729 }
730
732
742 BCRSMatrix (size_type _n, size_type _m, size_type _avg, double _overflowsize, BuildMode bm)
743 : build_mode(bm), ready(notAllocated), n(0), m(0), nnz_(0),
744 allocationSize_(0), r(0), a(0),
745 avg(_avg), overflowsize(_overflowsize)
746 {
747 if (bm != implicit)
748 DUNE_THROW(BCRSMatrixError,"Only call this constructor when using the implicit build mode");
749 // Prevent user from setting a negative overflowsize:
750 // 1) It doesn't make sense
751 // 2) We use a negative overflow value to indicate that the parameters
752 // have not been set yet
753 if (_overflowsize < 0.0)
754 DUNE_THROW(BCRSMatrixError,"You cannot set a negative overflow fraction");
755 implicit_allocate(_n,_m);
756 }
757
764 : build_mode(Mat.build_mode), ready(notAllocated), n(0), m(0), nnz_(0),
765 allocationSize_(0), r(0), a(0),
766 avg(Mat.avg), overflowsize(Mat.overflowsize)
767 {
768 if (!(Mat.ready == notAllocated || Mat.ready == built))
769 DUNE_THROW(InvalidStateException,"BCRSMatrix can only be copy-constructed when source matrix is completely empty (size not set) or fully built)");
770
771 // deep copy in global array
772 size_type _nnz = Mat.nonzeroes();
773
774 j_ = Mat.j_; // enable column index sharing, release array in case of row-wise allocation
775 allocate(Mat.n, Mat.m, _nnz, true, true);
776
777 // build window structure
779 }
780
783 {
784 deallocate();
785 }
786
792 {
793 if (ready == notAllocated)
794 {
795 build_mode = bm;
796 return;
797 }
798 if (ready == building && (build_mode == unknown || build_mode == random || build_mode == row_wise) && (bm == row_wise || bm == random))
799 build_mode = bm;
800 else
801 DUNE_THROW(InvalidStateException, "Matrix structure cannot be changed at this stage anymore (ready == "<<ready<<").");
802 }
803
819 void setSize(size_type rows, size_type columns, size_type nnz=0)
820 {
821 // deallocate already setup memory
822 deallocate();
823
824 if (build_mode == implicit)
825 {
826 if (nnz>0)
827 DUNE_THROW(Dune::BCRSMatrixError,"number of non-zeroes may not be set in implicit mode, use setImplicitBuildModeParameters() instead");
828
829 // implicit allocates differently
830 implicit_allocate(rows,columns);
831 }
832 else
833 {
834 // allocate matrix memory
835 allocate(rows, columns, nnz, true, false);
836 }
837 }
838
847 void setImplicitBuildModeParameters(size_type _avg, double _overflow)
848 {
849 // Prevent user from setting a negative overflowsize:
850 // 1) It doesn't make sense
851 // 2) We use a negative overflow value to indicate that the parameters
852 // have not been set yet
853 if (_overflow < 0.0)
854 DUNE_THROW(BCRSMatrixError,"You cannot set a negative overflow fraction");
855
856 // make sure the parameters aren't changed after memory has been allocated
857 if (ready != notAllocated)
858 DUNE_THROW(InvalidStateException,"You cannot modify build mode parameters at this stage anymore");
859 avg = _avg;
860 overflowsize = _overflow;
861 }
862
870 {
871 // return immediately when self-assignment
872 if (&Mat==this) return *this;
873
874 if (!((ready == notAllocated || ready == built) && (Mat.ready == notAllocated || Mat.ready == built)))
875 DUNE_THROW(InvalidStateException,"BCRSMatrix can only be copied when both target and source are empty or fully built)");
876
877 // make it simple: ALWAYS throw away memory for a and j_
878 // and deallocate rows only if n != Mat.n
879 deallocate(n!=Mat.n);
880
881 // reallocate the rows if required
882 if (n>0 && n!=Mat.n) {
883 // free rows
884 for(row_type *riter=r+(n-1), *rend=r-1; riter!=rend; --riter)
885 rowAllocator_.destroy(riter);
886 rowAllocator_.deallocate(r,n);
887 }
888
889 nnz_ = Mat.nonzeroes();
890
891 // allocate a, share j_
892 j_ = Mat.j_;
893 allocate(Mat.n, Mat.m, nnz_, n!=Mat.n, true);
894
895 // build window structure
897 return *this;
898 }
899
902 {
903
904 if (!(ready == notAllocated || ready == built))
905 DUNE_THROW(InvalidStateException,"Scalar assignment only works on fully built BCRSMatrix)");
906
907 for (size_type i=0; i<n; i++) r[i] = k;
908 return *this;
909 }
910
911 //===== row-wise creation interface
912
915 {
916 public:
919 : Mat(_Mat), i(_i), nnz(0), current_row(nullptr, Mat.j_.get(), 0)
920 {
921 if (Mat.build_mode == unknown && Mat.ready == building)
922 {
923 Mat.build_mode = row_wise;
924 }
925 if (i==0 && Mat.ready != building)
926 DUNE_THROW(BCRSMatrixError,"creation only allowed for uninitialized matrix");
927 if(Mat.build_mode!=row_wise)
928 DUNE_THROW(BCRSMatrixError,"creation only allowed if row wise allocation was requested in the constructor");
929 if(i==0 && _Mat.N()==0)
930 // empty Matrix is always built.
931 Mat.ready = built;
932 }
933
936 {
937 // this should only be called if matrix is in creation
938 if (Mat.ready != building)
939 DUNE_THROW(BCRSMatrixError,"matrix already built up");
940
941 // row i is defined through the pattern
942 // get memory for the row and initialize the j_ array
943 // this depends on the allocation mode
944
945 // compute size of the row
946 size_type s = pattern.size();
947
948 if(s>0) {
949 // update number of nonzeroes including this row
950 nnz += s;
951
952 // alloc memory / set window
953 if (Mat.nnz_ > 0)
954 {
955 // memory is allocated in one long array
956
957 // check if that memory is sufficient
958 if (nnz > Mat.nnz_)
959 DUNE_THROW(BCRSMatrixError,"allocated nnz too small");
960
961 // set row i
962 Mat.r[i].set(s,nullptr,current_row.getindexptr());
963 current_row.setindexptr(current_row.getindexptr()+s);
964 }else{
965 // memory is allocated individually per row
966 // allocate and set row i
967 B* b = Mat.allocator_.allocate(s);
968 // use placement new to call constructor that allocates
969 // additional memory.
970 new (b) B[s];
971 size_type* j = Mat.sizeAllocator_.allocate(s);
972 Mat.r[i].set(s,b,j);
973 }
974 }else
975 // setup empty row
976 Mat.r[i].set(0,nullptr,nullptr);
977
978 // initialize the j array for row i from pattern
979 std::copy(pattern.cbegin(), pattern.cend(), Mat.r[i].getindexptr());
980
981 // now go to next row
982 i++;
983 pattern.clear();
984
985 // check if this was last row
986 if (i==Mat.n)
987 {
988 Mat.ready = built;
989 if(Mat.nnz_ > 0)
990 {
991 // Set nnz to the exact number of nonzero blocks inserted
992 // as some methods rely on it
993 Mat.nnz_ = nnz;
994 // allocate data array
995 Mat.allocateData();
996 Mat.setDataPointers();
997 }
998 }
999 // done
1000 return *this;
1001 }
1002
1004 bool operator!= (const CreateIterator& it) const
1005 {
1006 return (i!=it.i) || (&Mat!=&it.Mat);
1007 }
1008
1010 bool operator== (const CreateIterator& it) const
1011 {
1012 return (i==it.i) && (&Mat==&it.Mat);
1013 }
1014
1017 {
1018 return i;
1019 }
1020
1023 {
1024 pattern.insert(j);
1025 }
1026
1029 {
1030 return pattern.find(j) != pattern.end();
1031 }
1038 {
1039 return pattern.size();
1040 }
1041
1042 private:
1043 BCRSMatrix& Mat; // the matrix we are defining
1044 size_type i; // current row to be defined
1045 size_type nnz; // count total number of nonzeros
1046 typedef std::set<size_type,std::less<size_type> > PatternType;
1047 PatternType pattern; // used to compile entries in a row
1048 row_type current_row; // row pointing to the current row to setup
1049 };
1050
1052 friend class CreateIterator;
1053
1056 {
1057 return CreateIterator(*this,0);
1058 }
1059
1062 {
1063 return CreateIterator(*this,n);
1064 }
1065
1066
1067 //===== random creation interface
1068
1076 {
1077 if (build_mode!=random)
1078 DUNE_THROW(BCRSMatrixError,"requires random build mode");
1079 if (ready != building)
1080 DUNE_THROW(BCRSMatrixError,"matrix row sizes already built up");
1081
1082 r[i].setsize(s);
1083 }
1084
1087 {
1088#ifdef DUNE_ISTL_WITH_CHECKING
1089 if (r==0) DUNE_THROW(BCRSMatrixError,"row not initialized yet");
1090 if (i>=n) DUNE_THROW(BCRSMatrixError,"index out of range");
1091#endif
1092 return r[i].getsize();
1093 }
1094
1097 {
1098 if (build_mode!=random)
1099 DUNE_THROW(BCRSMatrixError,"requires random build mode");
1100 if (ready != building)
1101 DUNE_THROW(BCRSMatrixError,"matrix row sizes already built up");
1102
1103 r[i].setsize(r[i].getsize()+s);
1104 }
1105
1108 {
1109 if (build_mode!=random)
1110 DUNE_THROW(BCRSMatrixError,"requires random build mode");
1111 if (ready != building)
1112 DUNE_THROW(BCRSMatrixError,"matrix row sizes already built up");
1113
1114 // compute total size, check positivity
1115 size_type total=0;
1116 for (size_type i=0; i<n; i++)
1117 {
1118 total += r[i].getsize();
1119 }
1120
1121 if(nnz_ == 0)
1122 // allocate/check memory
1123 allocate(n,m,total,false,false);
1124 else if(nnz_ < total)
1125 DUNE_THROW(BCRSMatrixError,"Specified number of nonzeros ("<<nnz_<<") not "
1126 <<"sufficient for calculated nonzeros ("<<total<<"! ");
1127
1128 // set the window pointers correctly
1130
1131 // initialize j_ array with m (an invalid column index)
1132 // this indicates an unused entry
1133 for (size_type k=0; k<nnz_; k++)
1134 j_.get()[k] = m;
1135 ready = rowSizesBuilt;
1136 }
1137
1139
1150 {
1151 if (build_mode!=random)
1152 DUNE_THROW(BCRSMatrixError,"requires random build mode");
1153 if (ready==built)
1154 DUNE_THROW(BCRSMatrixError,"matrix already built up");
1155 if (ready==building)
1156 DUNE_THROW(BCRSMatrixError,"matrix row sizes not built up yet");
1157 if (ready==notAllocated)
1158 DUNE_THROW(BCRSMatrixError,"matrix size not set and no memory allocated yet");
1159
1160 if (col >= m)
1161 DUNE_THROW(BCRSMatrixError,"column index exceeds matrix size");
1162
1163 // get row range
1164 size_type* const first = r[row].getindexptr();
1165 size_type* const last = first + r[row].getsize();
1166
1167 // find correct insertion position for new column index
1168 size_type* pos = std::lower_bound(first,last,col);
1169
1170 // check if index is already in row
1171 if (pos!=last && *pos == col) return;
1172
1173 // find end of already inserted column indices
1174 size_type* end = std::lower_bound(pos,last,m);
1175 if (end==last)
1176 DUNE_THROW(BCRSMatrixError,"row is too small");
1177
1178 // insert new column index at correct position
1179 std::copy_backward(pos,end,end+1);
1180 *pos = col;
1181 }
1182
1184
1191 template<typename It>
1192 void setIndices(size_type row, It begin, It end)
1193 {
1194 size_type row_size = r[row].size();
1195 size_type* col_begin = r[row].getindexptr();
1196 size_type* col_end;
1197 // consistency check between allocated row size and number of passed column indices
1198 if ((col_end = std::copy(begin,end,r[row].getindexptr())) != col_begin + row_size)
1199 DUNE_THROW(BCRSMatrixError,"Given size of row " << row
1200 << " (" << row_size
1201 << ") does not match number of passed entries (" << (col_end - col_begin) << ")");
1202 std::sort(col_begin,col_end);
1203 }
1204
1207 {
1208 if (build_mode!=random)
1209 DUNE_THROW(BCRSMatrixError,"requires random build mode");
1210 if (ready==built)
1211 DUNE_THROW(BCRSMatrixError,"matrix already built up");
1212 if (ready==building)
1213 DUNE_THROW(BCRSMatrixError,"row sizes are not built up yet");
1214 if (ready==notAllocated)
1215 DUNE_THROW(BCRSMatrixError,"matrix size not set and no memory allocated yet");
1216
1217 // check if there are undefined indices
1218 RowIterator endi=end();
1219 for (RowIterator i=begin(); i!=endi; ++i)
1220 {
1221 ColIterator endj = (*i).end();
1222 for (ColIterator j=(*i).begin(); j!=endj; ++j) {
1223 if (j.index() >= m) {
1224 dwarn << "WARNING: size of row "<< i.index()<<" is "<<j.offset()<<". But was specified as being "<< (*i).end().offset()
1225 <<". This means you are wasting valuable space and creating additional cache misses!"<<std::endl;
1226 nnz_ -= ((*i).end().offset() - j.offset());
1227 r[i.index()].setsize(j.offset());
1228 break;
1229 }
1230 }
1231 }
1232
1233 allocateData();
1235
1236 // if not, set matrix to built
1237 ready = built;
1238 }
1239
1240 //===== implicit creation interface
1241
1243
1255 {
1256#ifdef DUNE_ISTL_WITH_CHECKING
1257 if (build_mode!=implicit)
1258 DUNE_THROW(BCRSMatrixError,"requires implicit build mode");
1259 if (ready==built)
1260 DUNE_THROW(BCRSMatrixError,"matrix already built up, use operator[] for entry access now");
1261 if (ready==notAllocated)
1262 DUNE_THROW(BCRSMatrixError,"matrix size not set and no memory allocated yet");
1263 if (ready!=building)
1264 DUNE_THROW(InvalidStateException,"You may only use entry() during the 'building' stage");
1265
1266 if (row >= n)
1267 DUNE_THROW(BCRSMatrixError,"row index exceeds matrix size");
1268 if (col >= m)
1269 DUNE_THROW(BCRSMatrixError,"column index exceeds matrix size");
1270#endif
1271
1272 size_type* begin = r[row].getindexptr();
1273 size_type* end = begin + r[row].getsize();
1274
1275 size_type* pos = std::find(begin, end, col);
1276
1277 //treat the case that there was a match in the array
1278 if (pos != end)
1279 if (*pos == col)
1280 {
1281 std::ptrdiff_t offset = pos - r[row].getindexptr();
1282 B* aptr = r[row].getptr() + offset;
1283
1284 return *aptr;
1285 }
1286
1287 //determine whether overflow has to be taken into account or not
1288 if (r[row].getsize() == avg)
1289 return overflow[std::make_pair(row,col)];
1290 else
1291 {
1292 //modify index array
1293 *end = col;
1294
1295 //do simulatenous operations on data array a
1296 std::ptrdiff_t offset = end - r[row].getindexptr();
1297 B* apos = r[row].getptr() + offset;
1298
1299 //increase rowsize
1300 r[row].setsize(r[row].getsize()+1);
1301
1302 //return reference to the newly created entry
1303 return *apos;
1304 }
1305 }
1306
1308
1319 {
1320 if (build_mode!=implicit)
1321 DUNE_THROW(BCRSMatrixError,"requires implicit build mode");
1322 if (ready==built)
1323 DUNE_THROW(BCRSMatrixError,"matrix already built up, no more need for compression");
1324 if (ready==notAllocated)
1325 DUNE_THROW(BCRSMatrixError,"matrix size not set and no memory allocated yet");
1326 if (ready!=building)
1327 DUNE_THROW(InvalidStateException,"You may only call compress() at the end of the 'building' stage");
1328
1329 //calculate statistics
1331 stats.overflow_total = overflow.size();
1332 stats.maximum = 0;
1333
1334 //get insertion iterators pointing to one before start (for later use of ++it)
1335 size_type* jiit = j_.get();
1336 B* aiit = a;
1337
1338 //get iterator to the smallest overflow element
1339 typename OverflowType::iterator oit = overflow.begin();
1340
1341 //store a copy of index pointers on which to perform sortation
1342 std::vector<size_type*> perm;
1343
1344 //iterate over all rows and copy elements into their position in the compressed array
1345 for (size_type i=0; i<n; i++)
1346 {
1347 //get old pointers into a and j and size without overflow changes
1348 size_type* begin = r[i].getindexptr();
1349 //B* apos = r[i].getptr();
1350 size_type size = r[i].getsize();
1351
1352 perm.resize(size);
1353
1354 typename std::vector<size_type*>::iterator it = perm.begin();
1355 for (size_type* iit = begin; iit < begin + size; ++iit, ++it)
1356 *it = iit;
1357
1358 //sort permutation array
1359 std::sort(perm.begin(),perm.end(),PointerCompare<size_type>());
1360
1361 //change row window pointer to their new positions
1362 r[i].setindexptr(jiit);
1363 r[i].setptr(aiit);
1364
1365 for (it = perm.begin(); it != perm.end(); ++it)
1366 {
1367 //check whether there are elements in the overflow area which take precedence
1368 while ((oit!=overflow.end()) && (oit->first < std::make_pair(i,**it)))
1369 {
1370 //check whether there is enough memory to write to
1371 if (jiit > begin)
1373 "Allocated memory for BCRSMatrix exhausted during compress()!"
1374 "Please increase either the average number of entries per row or the overflow fraction."
1375 );
1376 //copy an element from the overflow area to the insertion position in a and j
1377 *jiit = oit->first.second;
1378 ++jiit;
1379 *aiit = oit->second;
1380 ++aiit;
1381 ++oit;
1382 r[i].setsize(r[i].getsize()+1);
1383 }
1384
1385 //check whether there is enough memory to write to
1386 if (jiit > begin)
1388 "Allocated memory for BCRSMatrix exhausted during compress()!"
1389 "Please increase either the average number of entries per row or the overflow fraction."
1390 );
1391
1392 //copy element from array
1393 *jiit = **it;
1394 ++jiit;
1395 B* apos = *it - j_.get() + a;
1396 *aiit = *apos;
1397 ++aiit;
1398 }
1399
1400 //copy remaining elements from the overflow area
1401 while ((oit!=overflow.end()) && (oit->first.first == i))
1402 {
1403 //check whether there is enough memory to write to
1404 if (jiit > begin)
1406 "Allocated memory for BCRSMatrix exhausted during compress()!"
1407 "Please increase either the average number of entries per row or the overflow fraction."
1408 );
1409
1410 //copy and element from the overflow area to the insertion position in a and j
1411 *jiit = oit->first.second;
1412 ++jiit;
1413 *aiit = oit->second;
1414 ++aiit;
1415 ++oit;
1416 r[i].setsize(r[i].getsize()+1);
1417 }
1418
1419 // update maximum row size
1420 if (r[i].getsize()>stats.maximum)
1421 stats.maximum = r[i].getsize();
1422 }
1423
1424 // overflow area may be cleared
1425 overflow.clear();
1426
1427 //determine average number of entries and memory usage
1428 std::ptrdiff_t diff = (r[n-1].getindexptr() + r[n-1].getsize() - j_.get());
1429 nnz_ = diff;
1430 stats.avg = (double) (nnz_) / (double) n;
1431 stats.mem_ratio = (double) (nnz_) / (double) allocationSize_;
1432
1433 //matrix is now built
1434 ready = built;
1435
1436 return stats;
1437 }
1438
1439 //===== vector space arithmetic
1440
1443 {
1444#ifdef DUNE_ISTL_WITH_CHECKING
1445 if (ready != built)
1446 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1447#endif
1448
1449 if (nnz_ > 0)
1450 {
1451 // process 1D array
1452 for (size_type i=0; i<nnz_; i++)
1453 a[i] *= k;
1454 }
1455 else
1456 {
1457 RowIterator endi=end();
1458 for (RowIterator i=begin(); i!=endi; ++i)
1459 {
1460 ColIterator endj = (*i).end();
1461 for (ColIterator j=(*i).begin(); j!=endj; ++j)
1462 (*j) *= k;
1463 }
1464 }
1465
1466 return *this;
1467 }
1468
1471 {
1472#ifdef DUNE_ISTL_WITH_CHECKING
1473 if (ready != built)
1474 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1475#endif
1476
1477 if (nnz_ > 0)
1478 {
1479 // process 1D array
1480 for (size_type i=0; i<nnz_; i++)
1481 a[i] /= k;
1482 }
1483 else
1484 {
1485 RowIterator endi=end();
1486 for (RowIterator i=begin(); i!=endi; ++i)
1487 {
1488 ColIterator endj = (*i).end();
1489 for (ColIterator j=(*i).begin(); j!=endj; ++j)
1490 (*j) /= k;
1491 }
1492 }
1493
1494 return *this;
1495 }
1496
1497
1504 {
1505#ifdef DUNE_ISTL_WITH_CHECKING
1506 if (ready != built || b.ready != built)
1507 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1508 if(N()!=b.N() || M() != b.M())
1509 DUNE_THROW(RangeError, "Matrix sizes do not match!");
1510#endif
1511 RowIterator endi=end();
1512 ConstRowIterator j=b.begin();
1513 for (RowIterator i=begin(); i!=endi; ++i, ++j) {
1514 i->operator+=(*j);
1515 }
1516
1517 return *this;
1518 }
1519
1526 {
1527#ifdef DUNE_ISTL_WITH_CHECKING
1528 if (ready != built || b.ready != built)
1529 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1530 if(N()!=b.N() || M() != b.M())
1531 DUNE_THROW(RangeError, "Matrix sizes do not match!");
1532#endif
1533 RowIterator endi=end();
1534 ConstRowIterator j=b.begin();
1535 for (RowIterator i=begin(); i!=endi; ++i, ++j) {
1536 i->operator-=(*j);
1537 }
1538
1539 return *this;
1540 }
1541
1551 {
1552#ifdef DUNE_ISTL_WITH_CHECKING
1553 if (ready != built || b.ready != built)
1554 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1555 if(N()!=b.N() || M() != b.M())
1556 DUNE_THROW(RangeError, "Matrix sizes do not match!");
1557#endif
1558 RowIterator endi=end();
1559 ConstRowIterator j=b.begin();
1560 for(RowIterator i=begin(); i!=endi; ++i, ++j)
1561 i->axpy(alpha, *j);
1562
1563 return *this;
1564 }
1565
1566 //===== linear maps
1567
1569 template<class X, class Y>
1570 void mv (const X& x, Y& y) const
1571 {
1572#ifdef DUNE_ISTL_WITH_CHECKING
1573 if (ready != built)
1574 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1575 if (x.N()!=M()) DUNE_THROW(BCRSMatrixError,
1576 "Size mismatch: M: " << N() << "x" << M() << " x: " << x.N());
1577 if (y.N()!=N()) DUNE_THROW(BCRSMatrixError,
1578 "Size mismatch: M: " << N() << "x" << M() << " y: " << y.N());
1579#endif
1580 ConstRowIterator endi=end();
1581 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1582 {
1583 y[i.index()]=0;
1584 ConstColIterator endj = (*i).end();
1585 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1586 {
1587 auto&& xj = Impl::asVector(x[j.index()]);
1588 auto&& yi = Impl::asVector(y[i.index()]);
1589 Impl::asMatrix(*j).umv(xj, yi);
1590 }
1591 }
1592 }
1593
1595 template<class X, class Y>
1596 void umv (const X& x, Y& y) const
1597 {
1598#ifdef DUNE_ISTL_WITH_CHECKING
1599 if (ready != built)
1600 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1601 if (x.N()!=M()) DUNE_THROW(BCRSMatrixError,"index out of range");
1602 if (y.N()!=N()) DUNE_THROW(BCRSMatrixError,"index out of range");
1603#endif
1604 ConstRowIterator endi=end();
1605 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1606 {
1607 ConstColIterator endj = (*i).end();
1608 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1609 {
1610 auto&& xj = Impl::asVector(x[j.index()]);
1611 auto&& yi = Impl::asVector(y[i.index()]);
1612 Impl::asMatrix(*j).umv(xj,yi);
1613 }
1614 }
1615 }
1616
1618 template<class X, class Y>
1619 void mmv (const X& x, Y& y) const
1620 {
1621#ifdef DUNE_ISTL_WITH_CHECKING
1622 if (ready != built)
1623 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1624 if (x.N()!=M()) DUNE_THROW(BCRSMatrixError,"index out of range");
1625 if (y.N()!=N()) DUNE_THROW(BCRSMatrixError,"index out of range");
1626#endif
1627 ConstRowIterator endi=end();
1628 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1629 {
1630 ConstColIterator endj = (*i).end();
1631 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1632 {
1633 auto&& xj = Impl::asVector(x[j.index()]);
1634 auto&& yi = Impl::asVector(y[i.index()]);
1635 Impl::asMatrix(*j).mmv(xj,yi);
1636 }
1637 }
1638 }
1639
1641 template<class X, class Y, class F>
1642 void usmv (F&& alpha, const X& x, Y& y) const
1643 {
1644#ifdef DUNE_ISTL_WITH_CHECKING
1645 if (ready != built)
1646 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1647 if (x.N()!=M()) DUNE_THROW(BCRSMatrixError,"index out of range");
1648 if (y.N()!=N()) DUNE_THROW(BCRSMatrixError,"index out of range");
1649#endif
1650 ConstRowIterator endi=end();
1651 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1652 {
1653 ConstColIterator endj = (*i).end();
1654 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1655 {
1656 auto&& xj = Impl::asVector(x[j.index()]);
1657 auto&& yi = Impl::asVector(y[i.index()]);
1658 Impl::asMatrix(*j).usmv(alpha,xj,yi);
1659 }
1660 }
1661 }
1662
1664 template<class X, class Y>
1665 void mtv (const X& x, Y& y) const
1666 {
1667#ifdef DUNE_ISTL_WITH_CHECKING
1668 if (ready != built)
1669 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1670 if (x.N()!=N()) DUNE_THROW(BCRSMatrixError,"index out of range");
1671 if (y.N()!=M()) DUNE_THROW(BCRSMatrixError,"index out of range");
1672#endif
1673 for(size_type i=0; i<y.N(); ++i)
1674 y[i]=0;
1675 umtv(x,y);
1676 }
1677
1679 template<class X, class Y>
1680 void umtv (const X& x, Y& y) const
1681 {
1682#ifdef DUNE_ISTL_WITH_CHECKING
1683 if (ready != built)
1684 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1685 if (x.N()!=N()) DUNE_THROW(BCRSMatrixError,"index out of range");
1686 if (y.N()!=M()) DUNE_THROW(BCRSMatrixError,"index out of range");
1687#endif
1688 ConstRowIterator endi=end();
1689 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1690 {
1691 ConstColIterator endj = (*i).end();
1692 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1693 {
1694 auto&& xi = Impl::asVector(x[i.index()]);
1695 auto&& yj = Impl::asVector(y[j.index()]);
1696 Impl::asMatrix(*j).umtv(xi,yj);
1697 }
1698 }
1699 }
1700
1702 template<class X, class Y>
1703 void mmtv (const X& x, Y& y) const
1704 {
1705#ifdef DUNE_ISTL_WITH_CHECKING
1706 if (x.N()!=N()) DUNE_THROW(BCRSMatrixError,"index out of range");
1707 if (y.N()!=M()) DUNE_THROW(BCRSMatrixError,"index out of range");
1708#endif
1709 ConstRowIterator endi=end();
1710 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1711 {
1712 ConstColIterator endj = (*i).end();
1713 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1714 {
1715 auto&& xi = Impl::asVector(x[i.index()]);
1716 auto&& yj = Impl::asVector(y[j.index()]);
1717 Impl::asMatrix(*j).mmtv(xi,yj);
1718 }
1719 }
1720 }
1721
1723 template<class X, class Y>
1724 void usmtv (const field_type& alpha, const X& x, Y& y) const
1725 {
1726#ifdef DUNE_ISTL_WITH_CHECKING
1727 if (ready != built)
1728 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1729 if (x.N()!=N()) DUNE_THROW(BCRSMatrixError,"index out of range");
1730 if (y.N()!=M()) DUNE_THROW(BCRSMatrixError,"index out of range");
1731#endif
1732 ConstRowIterator endi=end();
1733 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1734 {
1735 ConstColIterator endj = (*i).end();
1736 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1737 {
1738 auto&& xi = Impl::asVector(x[i.index()]);
1739 auto&& yj = Impl::asVector(y[j.index()]);
1740 Impl::asMatrix(*j).usmtv(alpha,xi,yj);
1741 }
1742 }
1743 }
1744
1746 template<class X, class Y>
1747 void umhv (const X& x, Y& y) const
1748 {
1749#ifdef DUNE_ISTL_WITH_CHECKING
1750 if (ready != built)
1751 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1752 if (x.N()!=N()) DUNE_THROW(BCRSMatrixError,"index out of range");
1753 if (y.N()!=M()) DUNE_THROW(BCRSMatrixError,"index out of range");
1754#endif
1755 ConstRowIterator endi=end();
1756 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1757 {
1758 ConstColIterator endj = (*i).end();
1759 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1760 {
1761 auto&& xi = Impl::asVector(x[i.index()]);
1762 auto&& yj = Impl::asVector(y[j.index()]);
1763 Impl::asMatrix(*j).umhv(xi,yj);
1764 }
1765 }
1766 }
1767
1769 template<class X, class Y>
1770 void mmhv (const X& x, Y& y) const
1771 {
1772#ifdef DUNE_ISTL_WITH_CHECKING
1773 if (ready != built)
1774 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1775 if (x.N()!=N()) DUNE_THROW(BCRSMatrixError,"index out of range");
1776 if (y.N()!=M()) DUNE_THROW(BCRSMatrixError,"index out of range");
1777#endif
1778 ConstRowIterator endi=end();
1779 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1780 {
1781 ConstColIterator endj = (*i).end();
1782 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1783 {
1784 auto&& xi = Impl::asVector(x[i.index()]);
1785 auto&& yj = Impl::asVector(y[j.index()]);
1786 Impl::asMatrix(*j).mmhv(xi,yj);
1787 }
1788 }
1789 }
1790
1792 template<class X, class Y>
1793 void usmhv (const field_type& alpha, const X& x, Y& y) const
1794 {
1795#ifdef DUNE_ISTL_WITH_CHECKING
1796 if (ready != built)
1797 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1798 if (x.N()!=N()) DUNE_THROW(BCRSMatrixError,"index out of range");
1799 if (y.N()!=M()) DUNE_THROW(BCRSMatrixError,"index out of range");
1800#endif
1801 ConstRowIterator endi=end();
1802 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1803 {
1804 ConstColIterator endj = (*i).end();
1805 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1806 {
1807 auto&& xi = Impl::asVector(x[i.index()]);
1808 auto&& yj = Impl::asVector(y[j.index()]);
1809 Impl::asMatrix(*j).usmhv(alpha,xi,yj);
1810 }
1811 }
1812 }
1813
1814
1815 //===== norms
1816
1818 typename FieldTraits<field_type>::real_type frobenius_norm2 () const
1819 {
1820#ifdef DUNE_ISTL_WITH_CHECKING
1821 if (ready != built)
1822 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1823#endif
1824
1825 typename FieldTraits<field_type>::real_type sum=0;
1826
1827 for (auto&& row : (*this))
1828 for (auto&& entry : row)
1829 sum += Impl::asMatrix(entry).frobenius_norm2();
1830
1831 return sum;
1832 }
1833
1835 typename FieldTraits<field_type>::real_type frobenius_norm () const
1836 {
1837 return sqrt(frobenius_norm2());
1838 }
1839
1841 template <typename ft = field_type,
1842 typename std::enable_if<!HasNaN<ft>::value, int>::type = 0>
1843 typename FieldTraits<ft>::real_type infinity_norm() const {
1844 if (ready != built)
1845 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1846
1847 using real_type = typename FieldTraits<ft>::real_type;
1848 using std::max;
1849
1850 real_type norm = 0;
1851 for (auto const &x : *this) {
1852 real_type sum = 0;
1853 for (auto const &y : x)
1854 sum += Impl::asMatrix(y).infinity_norm();
1855 norm = max(sum, norm);
1856 }
1857 return norm;
1858 }
1859
1861 template <typename ft = field_type,
1862 typename std::enable_if<!HasNaN<ft>::value, int>::type = 0>
1863 typename FieldTraits<ft>::real_type infinity_norm_real() const {
1864 if (ready != built)
1865 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1866
1867 using real_type = typename FieldTraits<ft>::real_type;
1868 using std::max;
1869
1870 real_type norm = 0;
1871 for (auto const &x : *this) {
1872 real_type sum = 0;
1873 for (auto const &y : x)
1874 sum += Impl::asMatrix(y).infinity_norm_real();
1875 norm = max(sum, norm);
1876 }
1877 return norm;
1878 }
1879
1881 template <typename ft = field_type,
1882 typename std::enable_if<HasNaN<ft>::value, int>::type = 0>
1883 typename FieldTraits<ft>::real_type infinity_norm() const {
1884 if (ready != built)
1885 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1886
1887 using real_type = typename FieldTraits<ft>::real_type;
1888 using std::max;
1889
1890 real_type norm = 0;
1891 real_type isNaN = 1;
1892 for (auto const &x : *this) {
1893 real_type sum = 0;
1894 for (auto const &y : x)
1895 sum += Impl::asMatrix(y).infinity_norm();
1896 norm = max(sum, norm);
1897 isNaN += sum;
1898 }
1899
1900 return norm * (isNaN / isNaN);
1901 }
1902
1904 template <typename ft = field_type,
1905 typename std::enable_if<HasNaN<ft>::value, int>::type = 0>
1906 typename FieldTraits<ft>::real_type infinity_norm_real() const {
1907 if (ready != built)
1908 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1909
1910 using real_type = typename FieldTraits<ft>::real_type;
1911 using std::max;
1912
1913 real_type norm = 0;
1914 real_type isNaN = 1;
1915
1916 for (auto const &x : *this) {
1917 real_type sum = 0;
1918 for (auto const &y : x)
1919 sum += Impl::asMatrix(y).infinity_norm_real();
1920 norm = max(sum, norm);
1921 isNaN += sum;
1922 }
1923
1924 return norm * (isNaN / isNaN);
1925 }
1926
1927 //===== sizes
1928
1930 size_type N () const
1931 {
1932 return n;
1933 }
1934
1936 size_type M () const
1937 {
1938 return m;
1939 }
1940
1943 {
1944 // in case of row-wise allocation
1945 if( nnz_ <= 0 )
1946 nnz_ = std::accumulate( begin(), end(), size_type( 0 ), [] ( size_type s, const row_type &row ) { return s+row.getsize(); } );
1947 return nnz_;
1948 }
1949
1952 {
1953 return ready;
1954 }
1955
1958 {
1959 return build_mode;
1960 }
1961
1962 //===== query
1963
1965 bool exists (size_type i, size_type j) const
1966 {
1967#ifdef DUNE_ISTL_WITH_CHECKING
1968 if (i<0 || i>=n) DUNE_THROW(BCRSMatrixError,"row index out of range");
1969 if (j<0 || j>=m) DUNE_THROW(BCRSMatrixError,"column index out of range");
1970#endif
1971 return (r[i].size() && r[i].find(j) != r[i].end());
1972 }
1973
1974
1975 protected:
1976 // state information
1977 BuildMode build_mode; // row wise or whole matrix
1978 BuildStage ready; // indicate the stage the matrix building is in
1979
1980 // The allocator used for memory management
1981 typename A::template rebind<B>::other allocator_;
1982
1983 typename A::template rebind<row_type>::other rowAllocator_;
1984
1985 typename A::template rebind<size_type>::other sizeAllocator_;
1986
1987 // size of the matrix
1988 size_type n; // number of rows
1989 size_type m; // number of columns
1990 mutable size_type nnz_; // number of nonzeroes contained in the matrix
1991 size_type allocationSize_; //allocated size of a and j arrays, except in implicit mode: nnz_==allocationsSize_
1992 // zero means that memory is allocated separately for each row.
1993
1994 // the rows are dynamically allocated
1995 row_type* r; // [n] the individual rows having pointers into a,j arrays
1996
1997 // dynamically allocated memory
1998 B* a; // [allocationSize] non-zero entries of the matrix in row-wise ordering
1999 // If a single array of column indices is used, it can be shared
2000 // between different matrices with the same sparsity pattern
2001 std::shared_ptr<size_type> j_; // [allocationSize] column indices of entries
2002
2003 // additional data is needed in implicit buildmode
2004 size_type avg;
2005 double overflowsize;
2006
2007 typedef std::map<std::pair<size_type,size_type>, B> OverflowType;
2008 OverflowType overflow;
2009
2010 void setWindowPointers(ConstRowIterator row)
2011 {
2012 row_type current_row(a,j_.get(),0); // Pointers to current row data
2013 for (size_type i=0; i<n; i++, ++row) {
2014 // set row i
2015 size_type s = row->getsize();
2016
2017 if (s>0) {
2018 // setup pointers and size
2019 r[i].set(s,current_row.getptr(), current_row.getindexptr());
2020 // update pointer for next row
2021 current_row.setptr(current_row.getptr()+s);
2022 current_row.setindexptr(current_row.getindexptr()+s);
2023 } else{
2024 // empty row
2025 r[i].set(0,nullptr,nullptr);
2026 }
2027 }
2028 }
2029
2031
2036 {
2037 size_type* jptr = j_.get();
2038 for (size_type i=0; i<n; ++i, ++row) {
2039 // set row i
2040 size_type s = row->getsize();
2041
2042 if (s>0) {
2043 // setup pointers and size
2044 r[i].setsize(s);
2045 r[i].setindexptr(jptr);
2046 } else{
2047 // empty row
2048 r[i].set(0,nullptr,nullptr);
2049 }
2050
2051 // advance position in global array
2052 jptr += s;
2053 }
2054 }
2055
2057
2062 {
2063 B* aptr = a;
2064 for (size_type i=0; i<n; ++i) {
2065 // set row i
2066 if (r[i].getsize() > 0) {
2067 // setup pointers and size
2068 r[i].setptr(aptr);
2069 } else{
2070 // empty row
2071 r[i].set(0,nullptr,nullptr);
2072 }
2073
2074 // advance position in global array
2075 aptr += r[i].getsize();
2076 }
2077 }
2078
2081 {
2082 setWindowPointers(Mat.begin());
2083
2084 // copy data
2085 for (size_type i=0; i<n; i++) r[i] = Mat.r[i];
2086
2087 // finish off
2088 build_mode = row_wise; // dummy
2089 ready = built;
2090 }
2091
2097 void deallocate(bool deallocateRows=true)
2098 {
2099
2100 if (notAllocated)
2101 return;
2102
2103 if (allocationSize_>0)
2104 {
2105 // a,j_ have been allocated as one long vector
2106 j_.reset();
2107 if (a)
2108 {
2109 for(B *aiter=a+(allocationSize_-1), *aend=a-1; aiter!=aend; --aiter)
2110 allocator_.destroy(aiter);
2111 allocator_.deallocate(a,allocationSize_);
2112 a = nullptr;
2113 }
2114 }
2115 else if (r)
2116 {
2117 // check if memory for rows have been allocated individually
2118 for (size_type i=0; i<n; i++)
2119 if (r[i].getsize()>0)
2120 {
2121 for (B *col=r[i].getptr()+(r[i].getsize()-1),
2122 *colend = r[i].getptr()-1; col!=colend; --col) {
2123 allocator_.destroy(col);
2124 }
2125 sizeAllocator_.deallocate(r[i].getindexptr(),1);
2126 allocator_.deallocate(r[i].getptr(),1);
2127 // clear out row data in case we don't want to deallocate the rows
2128 // otherwise we might run into a double free problem here later
2129 r[i].set(0,nullptr,nullptr);
2130 }
2131 }
2132
2133 // deallocate the rows
2134 if (n>0 && deallocateRows && r) {
2135 for(row_type *riter=r+(n-1), *rend=r-1; riter!=rend; --riter)
2136 rowAllocator_.destroy(riter);
2137 rowAllocator_.deallocate(r,n);
2138 r = nullptr;
2139 }
2140
2141 // Mark matrix as not built at all.
2142 ready=notAllocated;
2143
2144 }
2145
2148 {
2149 typename A::template rebind<size_type>::other& sizeAllocator_;
2150
2151 public:
2152 Deallocator(typename A::template rebind<size_type>::other& sizeAllocator)
2153 : sizeAllocator_(sizeAllocator)
2154 {}
2155
2156 void operator()(size_type* p) { sizeAllocator_.deallocate(p,1); }
2157 };
2158
2159
2177 void allocate(size_type rows, size_type columns, size_type allocationSize, bool allocateRows, bool allocate_data)
2178 {
2179 // Store size
2180 n = rows;
2181 m = columns;
2182 nnz_ = allocationSize;
2183 allocationSize_ = allocationSize;
2184
2185 // allocate rows
2186 if(allocateRows) {
2187 if (n>0) {
2188 if (r)
2189 DUNE_THROW(InvalidStateException,"Rows have already been allocated, cannot allocate a second time");
2190 r = rowAllocator_.allocate(rows);
2191 // initialize row entries
2192 for(row_type* ri=r; ri!=r+rows; ++ri)
2193 rowAllocator_.construct(ri, row_type());
2194 }else{
2195 r = 0;
2196 }
2197 }
2198
2199 // allocate a and j_ array
2200 if (allocate_data)
2201 allocateData();
2202 if (allocationSize_>0) {
2203 // allocate column indices only if not yet present (enable sharing)
2204 if (!j_.get())
2205 j_.reset(sizeAllocator_.allocate(allocationSize_),Deallocator(sizeAllocator_));
2206 }else{
2207 j_.reset();
2208 }
2209
2210 // Mark the matrix as not built.
2211 ready = building;
2212 }
2213
2214 void allocateData()
2215 {
2216 if (a)
2217 DUNE_THROW(InvalidStateException,"Cannot allocate data array (already allocated)");
2218 if (allocationSize_>0) {
2219 a = allocator_.allocate(allocationSize_);
2220 // use placement new to call constructor that allocates
2221 // additional memory.
2222 new (a) B[allocationSize_];
2223 } else {
2224 a = nullptr;
2225 }
2226 }
2227
2234 {
2235 if (build_mode != implicit)
2236 DUNE_THROW(InvalidStateException,"implicit_allocate() may only be called in implicit build mode");
2237 if (ready != notAllocated)
2238 DUNE_THROW(InvalidStateException,"memory has already been allocated");
2239
2240 // check to make sure the user has actually set the parameters
2241 if (overflowsize < 0)
2242 DUNE_THROW(InvalidStateException,"You have to set the implicit build mode parameters before starting to build the matrix");
2243 //calculate size of overflow region, add buffer for row sort!
2244 size_type osize = (size_type) (_n*avg)*overflowsize + 4*avg;
2245 allocationSize_ = _n*avg + osize;
2246
2247 allocate(_n, _m, allocationSize_,true,true);
2248
2249 //set row pointers correctly
2250 size_type* jptr = j_.get() + osize;
2251 B* aptr = a + osize;
2252 for (size_type i=0; i<n; i++)
2253 {
2254 r[i].set(0,aptr,jptr);
2255 jptr = jptr + avg;
2256 aptr = aptr + avg;
2257 }
2258
2259 ready = building;
2260 }
2261 };
2262
2263
2266} // end namespace
2267
2268#endif
This file implements a vector space as a tensor product of a given vector space. The number of compon...
Error specific to BCRSMatrix.
Definition: istlexception.hh:21
Iterator class for sequential creation of blocks
Definition: bcrsmatrix.hh:915
bool operator==(const CreateIterator &it) const
equality
Definition: bcrsmatrix.hh:1010
CreateIterator & operator++()
prefix increment
Definition: bcrsmatrix.hh:935
size_type index() const
The number of the row that the iterator currently points to.
Definition: bcrsmatrix.hh:1016
bool operator!=(const CreateIterator &it) const
inequality
Definition: bcrsmatrix.hh:1004
CreateIterator(BCRSMatrix &_Mat, size_type _i)
constructor
Definition: bcrsmatrix.hh:918
void insert(size_type j)
put column index in row
Definition: bcrsmatrix.hh:1022
bool contains(size_type j)
return true if column index is in row
Definition: bcrsmatrix.hh:1028
size_type size() const
Get the current row size.
Definition: bcrsmatrix.hh:1037
Class used by shared_ptr to deallocate memory using the proper allocator.
Definition: bcrsmatrix.hh:2148
Iterator access to matrix rows
Definition: bcrsmatrix.hh:537
RealRowIterator()
empty constructor, use with care!
Definition: bcrsmatrix.hh:554
bool equals(const RealRowIterator< ValueType > &other) const
equality
Definition: bcrsmatrix.hh:582
std::remove_const< T >::type ValueType
The unqualified value type.
Definition: bcrsmatrix.hh:541
bool equals(const RealRowIterator< const ValueType > &other) const
equality
Definition: bcrsmatrix.hh:589
RealRowIterator(row_type *_p, size_type _i)
constructor
Definition: bcrsmatrix.hh:549
size_type index() const
return index
Definition: bcrsmatrix.hh:564
A sparse block matrix with compressed row storage.
Definition: bcrsmatrix.hh:425
typename Imp::BlockTraits< B >::field_type field_type
export the type representing the field
Definition: bcrsmatrix.hh:447
bool exists(size_type i, size_type j) const
return true if (i,j) is in pattern
Definition: bcrsmatrix.hh:1965
BuildStage buildStage() const
The current build stage of the matrix.
Definition: bcrsmatrix.hh:1951
Iterator begin()
Get iterator to first row.
Definition: bcrsmatrix.hh:633
friend class CreateIterator
allow CreateIterator to access internal data
Definition: bcrsmatrix.hh:1052
void copyWindowStructure(const BCRSMatrix &Mat)
Copy the window structure from another matrix.
Definition: bcrsmatrix.hh:2080
B & entry(size_type row, size_type col)
Returns reference to entry (row,col) of the matrix.
Definition: bcrsmatrix.hh:1254
void mmhv(const X &x, Y &y) const
y -= A^H x
Definition: bcrsmatrix.hh:1770
void usmhv(const field_type &alpha, const X &x, Y &y) const
y += alpha A^H x
Definition: bcrsmatrix.hh:1793
void umtv(const X &x, Y &y) const
y += A^T x
Definition: bcrsmatrix.hh:1680
RealRowIterator< const row_type > const_iterator
The const iterator over the matrix rows.
Definition: bcrsmatrix.hh:665
static constexpr unsigned int blocklevel
increment block level counter
Definition: bcrsmatrix.hh:465
BCRSMatrix(size_type _n, size_type _m, size_type _avg, double _overflowsize, BuildMode bm)
construct matrix with a known average number of entries per row
Definition: bcrsmatrix.hh:742
void allocate(size_type rows, size_type columns, size_type allocationSize, bool allocateRows, bool allocate_data)
Allocate memory for the matrix structure.
Definition: bcrsmatrix.hh:2177
BCRSMatrix & axpy(field_type alpha, const BCRSMatrix &b)
Add the scaled entries of another matrix to this one.
Definition: bcrsmatrix.hh:1550
FieldTraits< ft >::real_type infinity_norm_real() const
simplified infinity norm (uses Manhattan norm for complex values)
Definition: bcrsmatrix.hh:1863
~BCRSMatrix()
destructor
Definition: bcrsmatrix.hh:782
void deallocate(bool deallocateRows=true)
deallocate memory of the matrix.
Definition: bcrsmatrix.hh:2097
Iterator RowIterator
rename the iterators for easier access
Definition: bcrsmatrix.hh:659
row_type & operator[](size_type i)
random access to the rows
Definition: bcrsmatrix.hh:507
BCRSMatrix()
an empty matrix
Definition: bcrsmatrix.hh:707
void endrowsizes()
indicate that size of all rows is defined
Definition: bcrsmatrix.hh:1107
void incrementrowsize(size_type i, size_type s=1)
increment size of row i by s (1 by default)
Definition: bcrsmatrix.hh:1096
void mtv(const X &x, Y &y) const
y = A^T x
Definition: bcrsmatrix.hh:1665
void umhv(const X &x, Y &y) const
y += A^H x
Definition: bcrsmatrix.hh:1747
size_type nonzeroes() const
number of blocks that are stored (the number of blocks that possibly are nonzero)
Definition: bcrsmatrix.hh:1942
ConstIterator ConstRowIterator
rename the const row iterator for easier access
Definition: bcrsmatrix.hh:696
void setrowsize(size_type i, size_type s)
Set number of indices in row i to s.
Definition: bcrsmatrix.hh:1075
BCRSMatrix & operator*=(const field_type &k)
vector space multiplication with scalar
Definition: bcrsmatrix.hh:1442
RealRowIterator< row_type > iterator
The iterator over the (mutable matrix rows.
Definition: bcrsmatrix.hh:629
void usmtv(const field_type &alpha, const X &x, Y &y) const
y += alpha A^T x
Definition: bcrsmatrix.hh:1724
ConstIterator beforeBegin() const
Definition: bcrsmatrix.hh:690
Iterator beforeBegin()
Definition: bcrsmatrix.hh:653
BuildMode
we support two modes
Definition: bcrsmatrix.hh:468
@ implicit
Build entries randomly with an educated guess for the number of entries per row.
Definition: bcrsmatrix.hh:497
@ unknown
Build mode not set!
Definition: bcrsmatrix.hh:501
@ random
Build entries randomly.
Definition: bcrsmatrix.hh:488
@ row_wise
Build in a row-wise manner.
Definition: bcrsmatrix.hh:479
BCRSMatrix(size_type _n, size_type _m, size_type _nnz, BuildMode bm)
matrix with known number of nonzeroes
Definition: bcrsmatrix.hh:714
::Dune::CompressionStatistics< size_type > CompressionStatistics
The type for the statistics object returned by compress()
Definition: bcrsmatrix.hh:462
BCRSMatrix & operator-=(const BCRSMatrix &b)
Subtract the entries of another matrix from this one.
Definition: bcrsmatrix.hh:1525
BCRSMatrix(const BCRSMatrix &Mat)
copy constructor
Definition: bcrsmatrix.hh:763
Iterator end()
Get iterator to one beyond last row.
Definition: bcrsmatrix.hh:639
void setIndices(size_type row, It begin, It end)
Set all column indices for row from the given iterator range.
Definition: bcrsmatrix.hh:1192
void addindex(size_type row, size_type col)
add index (row,col) to the matrix
Definition: bcrsmatrix.hh:1149
row_type::Iterator ColIterator
Iterator for the entries of each row.
Definition: bcrsmatrix.hh:662
FieldTraits< field_type >::real_type frobenius_norm() const
frobenius norm: sqrt(sum over squared values of entries)
Definition: bcrsmatrix.hh:1835
A::size_type size_type
The type for the index access and the size.
Definition: bcrsmatrix.hh:459
BCRSMatrix & operator/=(const field_type &k)
vector space division by scalar
Definition: bcrsmatrix.hh:1470
BCRSMatrix & operator+=(const BCRSMatrix &b)
Add the entries of another matrix to this one.
Definition: bcrsmatrix.hh:1503
CreateIterator createend()
get create iterator pointing to one after the last block
Definition: bcrsmatrix.hh:1061
FieldTraits< field_type >::real_type frobenius_norm2() const
square of frobenius norm, need for block recursion
Definition: bcrsmatrix.hh:1818
Iterator beforeEnd()
Definition: bcrsmatrix.hh:646
row_type::ConstIterator ConstColIterator
Const iterator to the entries of a row.
Definition: bcrsmatrix.hh:699
void usmv(F &&alpha, const X &x, Y &y) const
y += alpha A x
Definition: bcrsmatrix.hh:1642
size_type getrowsize(size_type i) const
get current number of indices in row i
Definition: bcrsmatrix.hh:1086
size_type M() const
number of columns (counted in blocks)
Definition: bcrsmatrix.hh:1936
Imp::CompressedBlockVectorWindow< B, A > row_type
implement row_type with compressed vector
Definition: bcrsmatrix.hh:456
CreateIterator createbegin()
get initial create iterator
Definition: bcrsmatrix.hh:1055
BuildStage
Definition: bcrsmatrix.hh:428
@ rowSizesBuilt
The row sizes of the matrix are known.
Definition: bcrsmatrix.hh:439
@ built
The matrix structure is fully built.
Definition: bcrsmatrix.hh:441
@ notbuilt
Matrix is not built at all, no memory has been allocated, build mode and size can still be set.
Definition: bcrsmatrix.hh:430
@ notAllocated
Matrix is not built at all, no memory has been allocated, build mode and size can still be set.
Definition: bcrsmatrix.hh:432
@ building
Matrix is currently being built, some memory has been allocated, build mode and size are fixed.
Definition: bcrsmatrix.hh:434
BuildMode buildMode() const
The currently selected build mode of the matrix.
Definition: bcrsmatrix.hh:1957
void mmv(const X &x, Y &y) const
y -= A x
Definition: bcrsmatrix.hh:1619
FieldTraits< ft >::real_type infinity_norm() const
infinity norm (row sum norm, how to generalize for blocks?)
Definition: bcrsmatrix.hh:1843
void mv(const X &x, Y &y) const
y = A x
Definition: bcrsmatrix.hh:1570
B block_type
export the type representing the components
Definition: bcrsmatrix.hh:450
void mmtv(const X &x, Y &y) const
y -= A^T x
Definition: bcrsmatrix.hh:1703
void umv(const X &x, Y &y) const
y += A x
Definition: bcrsmatrix.hh:1596
void implicit_allocate(size_type _n, size_type _m)
organizes allocation implicit mode calculates correct array size to be allocated and sets the the win...
Definition: bcrsmatrix.hh:2233
BCRSMatrix(size_type _n, size_type _m, BuildMode bm)
matrix with unknown number of nonzeroes
Definition: bcrsmatrix.hh:723
void endindices()
indicate that all indices are defined, check consistency
Definition: bcrsmatrix.hh:1206
CompressionStatistics compress()
Finishes the buildstage in implicit mode.
Definition: bcrsmatrix.hh:1318
void setDataPointers()
Set data pointers for all rows.
Definition: bcrsmatrix.hh:2061
size_type N() const
number of rows (counted in blocks)
Definition: bcrsmatrix.hh:1930
void setBuildMode(BuildMode bm)
Sets the build mode of the matrix.
Definition: bcrsmatrix.hh:791
void setSize(size_type rows, size_type columns, size_type nnz=0)
Set the size of the matrix.
Definition: bcrsmatrix.hh:819
BCRSMatrix & operator=(const BCRSMatrix &Mat)
assignment
Definition: bcrsmatrix.hh:869
void setColumnPointers(ConstRowIterator row)
Copy row sizes from iterator range starting at row and set column index pointers for all rows.
Definition: bcrsmatrix.hh:2035
ConstIterator end() const
Get const iterator to one beyond last row.
Definition: bcrsmatrix.hh:676
ConstIterator begin() const
Get const iterator to first row.
Definition: bcrsmatrix.hh:670
void setImplicitBuildModeParameters(size_type _avg, double _overflow)
Set parameters needed for creation in implicit build mode.
Definition: bcrsmatrix.hh:847
A allocator_type
export the allocator type
Definition: bcrsmatrix.hh:453
ConstIterator beforeEnd() const
Definition: bcrsmatrix.hh:683
Proxy row object for entry access.
Definition: bcrsmatrix.hh:132
block_type & operator[](size_type j) const
Returns entry in column j.
Definition: bcrsmatrix.hh:137
A wrapper for uniform access to the BCRSMatrix during and after the build stage in implicit build mod...
Definition: bcrsmatrix.hh:112
Matrix::block_type block_type
The block_type of the underlying matrix.
Definition: bcrsmatrix.hh:120
ImplicitMatrixBuilder(Matrix &m)
Creates an ImplicitMatrixBuilder for matrix m.
Definition: bcrsmatrix.hh:165
M_ Matrix
The underlying matrix.
Definition: bcrsmatrix.hh:117
ImplicitMatrixBuilder(Matrix &m, size_type rows, size_type cols, size_type avg_cols_per_row, double overflow_fraction)
Sets up matrix m for implicit construction using the given parameters and creates an ImplicitBmatrixu...
Definition: bcrsmatrix.hh:189
size_type M() const
The number of columns in the matrix.
Definition: bcrsmatrix.hh:212
Matrix::size_type size_type
The size_type of the underlying matrix.
Definition: bcrsmatrix.hh:123
row_object operator[](size_type i) const
Returns a proxy for entries in row i.
Definition: bcrsmatrix.hh:200
size_type N() const
The number of rows in the matrix.
Definition: bcrsmatrix.hh:206
The overflow error used during implicit BCRSMatrix construction was exhausted.
Definition: istlexception.hh:34
Default exception if a function was called while the object is not in a valid state for that function...
Definition: exceptions.hh:279
A generic dynamic dense matrix.
Definition: matrix.hh:558
A::size_type size_type
Type for indices and sizes.
Definition: matrix.hh:574
T block_type
Export the type representing the components.
Definition: matrix.hh:565
Base class for stl conformant forward iterators.
Definition: iteratorfacades.hh:432
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_THROW(E, m)
Definition: exceptions.hh:216
T accumulate(Range &&range, T value, F &&f)
Accumulate values.
Definition: hybridutilities.hh:290
constexpr decltype(auto) elementAt(Container &&c, Index &&i)
Get element at given position from container.
Definition: hybridutilities.hh:134
auto max(ADLTag< 0 >, const V &v1, const V &v2)
implements binary Simd::max()
Definition: defaults.hh:79
DWarnType dwarn(std::cerr)
Stream for warnings indicating problems.
Definition: stdstreams.hh:159
This file implements iterator facade classes for writing stl conformant iterators.
Some handy generic functions for ISTL matrices.
Dune namespace.
Definition: alignedallocator.hh:14
Implements a scalar matrix view wrapper around an existing scalar.
Implements a scalar vector view wrapper around an existing scalar.
Standard Dune debug streams.
Statistics about compression achieved in implicit mode.
Definition: bcrsmatrix.hh:83
size_type overflow_total
total number of elements written to the overflow area during construction.
Definition: bcrsmatrix.hh:89
size_type maximum
maximum number of non-zeroes per row.
Definition: bcrsmatrix.hh:87
double avg
average number of non-zeroes per row.
Definition: bcrsmatrix.hh:85
double mem_ratio
fraction of wasted memory resulting from non-used overflow area.
Definition: bcrsmatrix.hh:94
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)