Dune Core Modules (unstable)

bcrsmatrix.hh
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright © DUNE Project contributors, see file LICENSE.md in module root
2// SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
3// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
4// vi: set et ts=4 sw=2 sts=2:
5
6#ifndef DUNE_ISTL_BCRSMATRIX_HH
7#define DUNE_ISTL_BCRSMATRIX_HH
8
9#include <cmath>
10#include <complex>
11#include <set>
12#include <iostream>
13#include <algorithm>
14#include <numeric>
15#include <vector>
16#include <map>
17#include <memory>
18
19#include "istlexception.hh"
20#include "bvector.hh"
21#include "matrixutils.hh"
22#include "matrixindexset.hh"
29
31
36namespace Dune {
37
77 template<typename M>
78 struct MatrixDimension;
79
81
87 template<typename size_type>
89 {
91 double avg;
93 size_type maximum;
95 size_type overflow_total;
97
100 double mem_ratio;
101 };
102
104
116 template<class M_>
118 {
119
120 public:
121
123 typedef M_ Matrix;
124
127
129 typedef typename Matrix::size_type size_type;
130
132
138 {
139
140 public:
141
144 {
145 return _m.entry(_i,j);
146 }
147
148#ifndef DOXYGEN
149
151 : _m(m)
152 , _i(i)
153 {}
154
155#endif
156
157 private:
158
159 Matrix& _m;
160 size_type _i;
161
162 };
163
165
172 : _m(m)
173 {
174 if (m.buildMode() != Matrix::implicit)
175 DUNE_THROW(BCRSMatrixError,"You can only create an ImplicitBuilder for a matrix in implicit build mode");
176 if (m.buildStage() != Matrix::building)
177 DUNE_THROW(BCRSMatrixError,"You can only create an ImplicitBuilder for a matrix with set size that has not been compressed() yet");
178 }
179
181
195 ImplicitMatrixBuilder(Matrix& m, size_type rows, size_type cols, size_type avg_cols_per_row, double overflow_fraction)
196 : _m(m)
197 {
198 if (m.buildStage() != Matrix::notAllocated)
199 DUNE_THROW(BCRSMatrixError,"You can only set up a matrix for this ImplicitBuilder if it has no memory allocated yet");
200 m.setBuildMode(Matrix::implicit);
201 m.setImplicitBuildModeParameters(avg_cols_per_row,overflow_fraction);
202 m.setSize(rows,cols);
203 }
204
207 {
208 return row_object(_m,i);
209 }
210
212 size_type N() const
213 {
214 return _m.N();
215 }
216
218 size_type M() const
219 {
220 return _m.M();
221 }
222
223 private:
224
225 Matrix& _m;
226
227 };
228
465 template<class B, class A=std::allocator<B> >
467 {
468 friend struct MatrixDimension<BCRSMatrix>;
469 public:
483 built=3
484 };
485
486 //===== type definitions and constants
487
489 using field_type = typename Imp::BlockTraits<B>::field_type;
490
492 typedef B block_type;
493
495 typedef A allocator_type;
496
498 typedef typename A::size_type size_type;
499
501 typedef Imp::CompressedBlockVectorWindow<B,size_type> row_type;
502
504 typedef ::Dune::CompressionStatistics<size_type> CompressionStatistics;
505
540 unknown
541 };
542
543 //===== random access interface to rows of the matrix
544
547 {
548#ifdef DUNE_ISTL_WITH_CHECKING
549 if (build_mode == implicit && ready != built)
550 DUNE_THROW(BCRSMatrixError,"You cannot use operator[] in implicit build mode before calling compress()");
551 if (r==0) DUNE_THROW(BCRSMatrixError,"row not initialized yet");
552 if (i>=n) DUNE_THROW(BCRSMatrixError,"index out of range");
553#endif
554 return r[i];
555 }
556
559 {
560#ifdef DUNE_ISTL_WITH_CHECKING
561 if (build_mode == implicit && ready != built)
562 DUNE_THROW(BCRSMatrixError,"You cannot use operator[] in implicit build mode before calling compress()");
563 if (built!=ready) DUNE_THROW(BCRSMatrixError,"row not initialized yet");
564 if (i>=n) DUNE_THROW(BCRSMatrixError,"index out of range");
565#endif
566 return r[i];
567 }
568
569
570 //===== iterator interface to rows of the matrix
571
573 template<class T>
575 : public RandomAccessIteratorFacade<RealRowIterator<T>, T>
576 {
577
578 public:
580 typedef typename std::remove_const<T>::type ValueType;
581
584 friend class RealRowIterator<const ValueType>;
585 friend class RealRowIterator<ValueType>;
586
589 : p(_p), i(_i)
590 {}
591
594 : p(0), i(0)
595 {}
596
598 : p(it.p), i(it.i)
599 {}
600
601
604 {
605 return i;
606 }
607
608 std::ptrdiff_t distanceTo(const RealRowIterator<ValueType>& other) const
609 {
610 assert(other.p==p);
611 return (other.i-i);
612 }
613
614 std::ptrdiff_t distanceTo(const RealRowIterator<const ValueType>& other) const
615 {
616 assert(other.p==p);
617 return (other.i-i);
618 }
619
621 bool equals (const RealRowIterator<ValueType>& other) const
622 {
623 assert(other.p==p);
624 return i==other.i;
625 }
626
628 bool equals (const RealRowIterator<const ValueType>& other) const
629 {
630 assert(other.p==p);
631 return i==other.i;
632 }
633
634 private:
636 void increment()
637 {
638 ++i;
639 }
640
642 void decrement()
643 {
644 --i;
645 }
646
647 void advance(std::ptrdiff_t diff)
648 {
649 i+=diff;
650 }
651
652 T& elementAt(std::ptrdiff_t diff) const
653 {
654 return p[i+diff];
655 }
656
658 row_type& dereference () const
659 {
660 return p[i];
661 }
662
663 row_type* p;
664 size_type i;
665 };
666
670
673 {
674 return Iterator(r,0);
675 }
676
679 {
680 return Iterator(r,n);
681 }
682
686 {
687 return Iterator(r,n-1);
688 }
689
693 {
694 return Iterator(r,-1);
695 }
696
699
701 typedef typename row_type::Iterator ColIterator;
702
706
707
710 {
711 return ConstIterator(r,0);
712 }
713
716 {
717 return ConstIterator(r,n);
718 }
719
723 {
724 return ConstIterator(r,n-1);
725 }
726
730 {
731 return ConstIterator(r,-1);
732 }
733
736
738 typedef typename row_type::ConstIterator ConstColIterator;
739
740 //===== constructors & resizers
741
742 // we use a negative compressionBufferSize to indicate that the implicit
743 // mode parameters have not been set yet
744
747 : build_mode(unknown), ready(notAllocated), n(0), m(0), nnz_(0),
748 allocationSize_(0), r(0), a(0),
749 avg(0), compressionBufferSize_(-1.0)
750 {}
751
754 : build_mode(bm), ready(notAllocated), n(0), m(0), nnz_(0),
755 allocationSize_(0), r(0), a(0),
756 avg(0), compressionBufferSize_(-1.0)
757 {
758 allocate(_n, _m, _nnz,true,false);
759 }
760
763 : build_mode(bm), ready(notAllocated), n(0), m(0), nnz_(0),
764 allocationSize_(0), r(0), a(0),
765 avg(0), compressionBufferSize_(-1.0)
766 {
767 allocate(_n, _m,0,true,false);
768 }
769
771
782 BCRSMatrix (size_type _n, size_type _m, size_type _avg, double compressionBufferSize, BuildMode bm)
783 : build_mode(bm), ready(notAllocated), n(0), m(0), nnz_(0),
784 allocationSize_(0), r(0), a(0),
785 avg(_avg), compressionBufferSize_(compressionBufferSize)
786 {
787 if (bm != implicit)
788 DUNE_THROW(BCRSMatrixError,"Only call this constructor when using the implicit build mode");
789 // Prevent user from setting a negative compression buffer size:
790 // 1) It doesn't make sense
791 // 2) We use a negative value to indicate that the parameters
792 // have not been set yet
793 if (compressionBufferSize_ < 0.0)
794 DUNE_THROW(BCRSMatrixError,"You cannot set a negative overflow fraction");
795 implicit_allocate(_n,_m);
796 }
797
804 : build_mode(Mat.build_mode), ready(notAllocated), n(0), m(0), nnz_(0),
805 allocationSize_(0), r(0), a(0),
806 avg(Mat.avg), compressionBufferSize_(Mat.compressionBufferSize_)
807 {
808 if (!(Mat.ready == notAllocated || Mat.ready == built))
809 DUNE_THROW(InvalidStateException,"BCRSMatrix can only be copy-constructed when source matrix is completely empty (size not set) or fully built)");
810
811 // deep copy in global array
812 size_type _nnz = Mat.nonzeroes();
813
814 j_ = Mat.j_; // enable column index sharing, release array in case of row-wise allocation
815 allocate(Mat.n, Mat.m, _nnz, true, true);
816
817 // build window structure
819 }
820
823 {
824 deallocate();
825 }
826
832 {
833 if (ready == notAllocated)
834 {
835 build_mode = bm;
836 return;
837 }
838 if (ready == building && (build_mode == unknown || build_mode == random || build_mode == row_wise) && (bm == row_wise || bm == random))
839 build_mode = bm;
840 else
841 DUNE_THROW(InvalidStateException, "Matrix structure cannot be changed at this stage anymore (ready == "<<ready<<").");
842 }
843
859 void setSize(size_type rows, size_type columns, size_type nnz=0)
860 {
861 // deallocate already setup memory
862 deallocate();
863
864 if (build_mode == implicit)
865 {
866 if (nnz>0)
867 DUNE_THROW(Dune::BCRSMatrixError,"number of non-zeroes may not be set in implicit mode, use setImplicitBuildModeParameters() instead");
868
869 // implicit allocates differently
870 implicit_allocate(rows,columns);
871 }
872 else
873 {
874 // allocate matrix memory
875 allocate(rows, columns, nnz, true, false);
876 }
877 }
878
887 void setImplicitBuildModeParameters(size_type _avg, double compressionBufferSize)
888 {
889 // Prevent user from setting a negative compression buffer size:
890 // 1) It doesn't make sense
891 // 2) We use a negative value to indicate that the parameters
892 // have not been set yet
893 if (compressionBufferSize < 0.0)
894 DUNE_THROW(BCRSMatrixError,"You cannot set a negative compressionBufferSize value");
895
896 // make sure the parameters aren't changed after memory has been allocated
897 if (ready != notAllocated)
898 DUNE_THROW(InvalidStateException,"You cannot modify build mode parameters at this stage anymore");
899 avg = _avg;
900 compressionBufferSize_ = compressionBufferSize;
901 }
902
910 {
911 // return immediately when self-assignment
912 if (&Mat==this) return *this;
913
914 if (!((ready == notAllocated || ready == built) && (Mat.ready == notAllocated || Mat.ready == built)))
915 DUNE_THROW(InvalidStateException,"BCRSMatrix can only be copied when both target and source are empty or fully built)");
916
917 // make it simple: ALWAYS throw away memory for a and j_
918 // and deallocate rows only if n != Mat.n
919 deallocate(n!=Mat.n);
920
921 // reallocate the rows if required
922 if (n>0 && n!=Mat.n) {
923 // free rows
924 for(row_type *riter=r+(n-1), *rend=r-1; riter!=rend; --riter)
925 std::allocator_traits<decltype(rowAllocator_)>::destroy(rowAllocator_, riter);
926 rowAllocator_.deallocate(r,n);
927 }
928
929 nnz_ = Mat.nonzeroes();
930
931 // allocate a, share j_
932 j_ = Mat.j_;
933 allocate(Mat.n, Mat.m, nnz_, n!=Mat.n, true);
934
935 // build window structure
937 return *this;
938 }
939
942 {
943
944 if (!(ready == notAllocated || ready == built))
945 DUNE_THROW(InvalidStateException,"Scalar assignment only works on fully built BCRSMatrix)");
946
947 for (size_type i=0; i<n; i++) r[i] = k;
948 return *this;
949 }
950
951 //===== row-wise creation interface
952
955 {
956 public:
959 : Mat(_Mat), i(_i), nnz(0), current_row(nullptr, Mat.j_.get(), 0)
960 {
961 if (Mat.build_mode == unknown && Mat.ready == building)
962 {
963 Mat.build_mode = row_wise;
964 }
965 if (i==0 && Mat.ready != building)
966 DUNE_THROW(BCRSMatrixError,"creation only allowed for uninitialized matrix");
967 if(Mat.build_mode!=row_wise)
968 DUNE_THROW(BCRSMatrixError,"creation only allowed if row wise allocation was requested in the constructor");
969 if(i==0 && _Mat.N()==0)
970 // empty Matrix is always built.
971 Mat.ready = built;
972 }
973
976 {
977 // this should only be called if matrix is in creation
978 if (Mat.ready != building)
979 DUNE_THROW(BCRSMatrixError,"matrix already built up");
980
981 // row i is defined through the pattern
982 // get memory for the row and initialize the j_ array
983 // this depends on the allocation mode
984
985 // compute size of the row
986 size_type s = pattern.size();
987
988 if(s>0) {
989 // update number of nonzeroes including this row
990 nnz += s;
991
992 // alloc memory / set window
993 if (Mat.nnz_ > 0)
994 {
995 // memory is allocated in one long array
996
997 // check if that memory is sufficient
998 if (nnz > Mat.nnz_)
999 DUNE_THROW(BCRSMatrixError,"allocated nnz too small");
1000
1001 // set row i
1002 Mat.r[i].set(s,nullptr,current_row.getindexptr());
1003 current_row.setindexptr(current_row.getindexptr()+s);
1004 }else{
1005 // memory is allocated individually per row
1006 // allocate and set row i
1007 B* b = Mat.allocator_.allocate(s);
1008 // use placement new to call constructor that allocates
1009 // additional memory.
1010 new (b) B[s];
1011 size_type* j = Mat.sizeAllocator_.allocate(s);
1012 Mat.r[i].set(s,b,j);
1013 }
1014 }else
1015 // setup empty row
1016 Mat.r[i].set(0,nullptr,nullptr);
1017
1018 // initialize the j array for row i from pattern
1019 std::visit([&](const auto& pattern){
1020 std::copy(pattern.cbegin(), pattern.cend(), Mat.r[i].getindexptr());
1021 }, pattern.storage());
1022
1023 // now go to next row
1024 i++;
1025 pattern.clear();
1026
1027 // check if this was last row
1028 if (i==Mat.n)
1029 {
1030 Mat.ready = built;
1031 if(Mat.nnz_ > 0)
1032 {
1033 // Set nnz to the exact number of nonzero blocks inserted
1034 // as some methods rely on it
1035 Mat.nnz_ = nnz;
1036 // allocate data array
1037 Mat.allocateData();
1038 Mat.setDataPointers();
1039 }
1040 }
1041 // done
1042 return *this;
1043 }
1044
1046 bool operator!= (const CreateIterator& it) const
1047 {
1048 return (i!=it.i) || (&Mat!=&it.Mat);
1049 }
1050
1052 bool operator== (const CreateIterator& it) const
1053 {
1054 return (i==it.i) && (&Mat==&it.Mat);
1055 }
1056
1059 {
1060 return i;
1061 }
1062
1065 {
1066 pattern.insert(j);
1067 }
1068
1071 {
1072 return pattern.contains(j);
1073 }
1080 {
1081 return pattern.size();
1082 }
1083
1084 private:
1085 BCRSMatrix& Mat; // the matrix we are defining
1086 size_type i; // current row to be defined
1087 size_type nnz; // count total number of nonzeros
1088 using PatternType = typename Impl::RowIndexSet;
1089 PatternType pattern; // used to compile entries in a row
1090 row_type current_row; // row pointing to the current row to setup
1091 };
1092
1094 friend class CreateIterator;
1095
1098 {
1099 return CreateIterator(*this,0);
1100 }
1101
1104 {
1105 return CreateIterator(*this,n);
1106 }
1107
1108
1109 //===== random creation interface
1110
1118 {
1119 if (build_mode!=random)
1120 DUNE_THROW(BCRSMatrixError,"requires random build mode");
1121 if (ready != building)
1122 DUNE_THROW(BCRSMatrixError,"matrix row sizes already built up");
1123
1124 r[i].setsize(s);
1125 }
1126
1129 {
1130#ifdef DUNE_ISTL_WITH_CHECKING
1131 if (r==0) DUNE_THROW(BCRSMatrixError,"row not initialized yet");
1132 if (i>=n) DUNE_THROW(BCRSMatrixError,"index out of range");
1133#endif
1134 return r[i].getsize();
1135 }
1136
1139 {
1140 if (build_mode!=random)
1141 DUNE_THROW(BCRSMatrixError,"requires random build mode");
1142 if (ready != building)
1143 DUNE_THROW(BCRSMatrixError,"matrix row sizes already built up");
1144
1145 r[i].setsize(r[i].getsize()+s);
1146 }
1147
1150 {
1151 if (build_mode!=random)
1152 DUNE_THROW(BCRSMatrixError,"requires random build mode");
1153 if (ready != building)
1154 DUNE_THROW(BCRSMatrixError,"matrix row sizes already built up");
1155
1156 // compute total size, check positivity
1157 size_type total=0;
1158 for (size_type i=0; i<n; i++)
1159 {
1160 total += r[i].getsize();
1161 }
1162
1163 if(nnz_ == 0)
1164 // allocate/check memory
1165 allocate(n,m,total,false,false);
1166 else if(nnz_ < total)
1167 DUNE_THROW(BCRSMatrixError,"Specified number of nonzeros ("<<nnz_<<") not "
1168 <<"sufficient for calculated nonzeros ("<<total<<"! ");
1169
1170 // set the window pointers correctly
1172
1173 // initialize j_ array with m (an invalid column index)
1174 // this indicates an unused entry
1175 for (size_type k=0; k<nnz_; k++)
1176 j_.get()[k] = m;
1177 ready = rowSizesBuilt;
1178 }
1179
1181
1192 {
1193 if (build_mode!=random)
1194 DUNE_THROW(BCRSMatrixError,"requires random build mode");
1195 if (ready==built)
1196 DUNE_THROW(BCRSMatrixError,"matrix already built up");
1197 if (ready==building)
1198 DUNE_THROW(BCRSMatrixError,"matrix row sizes not built up yet");
1199 if (ready==notAllocated)
1200 DUNE_THROW(BCRSMatrixError,"matrix size not set and no memory allocated yet");
1201
1202 if (col >= m)
1203 DUNE_THROW(BCRSMatrixError,"column index exceeds matrix size");
1204
1205 // get row range
1206 size_type* const first = r[row].getindexptr();
1207 size_type* const last = first + r[row].getsize();
1208
1209 // find correct insertion position for new column index
1210 size_type* pos = std::lower_bound(first,last,col);
1211
1212 // check if index is already in row
1213 if (pos!=last && *pos == col) return;
1214
1215 // find end of already inserted column indices
1216 size_type* end = std::lower_bound(pos,last,m);
1217 if (end==last)
1218 DUNE_THROW(BCRSMatrixError,"row is too small");
1219
1220 // insert new column index at correct position
1221 std::copy_backward(pos,end,end+1);
1222 *pos = col;
1223 }
1224
1226
1234 template<typename It>
1236 {
1237 size_type row_size = r[row].size();
1238 size_type* col_begin = r[row].getindexptr();
1239 size_type* col_end;
1240 // consistency check between allocated row size and number of passed column indices
1241 if ((col_end = std::copy(begin,end,r[row].getindexptr())) != col_begin + row_size)
1242 DUNE_THROW(BCRSMatrixError,"Given size of row " << row
1243 << " (" << row_size
1244 << ") does not match number of passed entries (" << (col_end - col_begin) << ")");
1245 }
1246
1247
1249
1257 template<typename It>
1258 void setIndices(size_type row, It begin, It end)
1259 {
1260 size_type row_size = r[row].size();
1261 size_type* col_begin = r[row].getindexptr();
1262 size_type* col_end;
1263 // consistency check between allocated row size and number of passed column indices
1264 if ((col_end = std::copy(begin,end,r[row].getindexptr())) != col_begin + row_size)
1265 DUNE_THROW(BCRSMatrixError,"Given size of row " << row
1266 << " (" << row_size
1267 << ") does not match number of passed entries (" << (col_end - col_begin) << ")");
1268 std::sort(col_begin,col_end);
1269 }
1270
1273 {
1274 if (build_mode!=random)
1275 DUNE_THROW(BCRSMatrixError,"requires random build mode");
1276 if (ready==built)
1277 DUNE_THROW(BCRSMatrixError,"matrix already built up");
1278 if (ready==building)
1279 DUNE_THROW(BCRSMatrixError,"row sizes are not built up yet");
1280 if (ready==notAllocated)
1281 DUNE_THROW(BCRSMatrixError,"matrix size not set and no memory allocated yet");
1282
1283 // check if there are undefined indices
1284 RowIterator endi=end();
1285 for (RowIterator i=begin(); i!=endi; ++i)
1286 {
1287 ColIterator endj = (*i).end();
1288 for (ColIterator j=(*i).begin(); j!=endj; ++j) {
1289 if (j.index() >= m) {
1290 dwarn << "WARNING: size of row "<< i.index()<<" is "<<j.offset()<<". But was specified as being "<< (*i).end().offset()
1291 <<". This means you are wasting valuable space and creating additional cache misses!"<<std::endl;
1292 nnz_ -= ((*i).end().offset() - j.offset());
1293 r[i.index()].setsize(j.offset());
1294 break;
1295 }
1296 }
1297 }
1298
1299 allocateData();
1301
1302 // if not, set matrix to built
1303 ready = built;
1304 }
1305
1306 //===== implicit creation interface
1307
1309
1321 {
1322#ifdef DUNE_ISTL_WITH_CHECKING
1323 if (build_mode!=implicit)
1324 DUNE_THROW(BCRSMatrixError,"requires implicit build mode");
1325 if (ready==built)
1326 DUNE_THROW(BCRSMatrixError,"matrix already built up, use operator[] for entry access now");
1327 if (ready==notAllocated)
1328 DUNE_THROW(BCRSMatrixError,"matrix size not set and no memory allocated yet");
1329 if (ready!=building)
1330 DUNE_THROW(InvalidStateException,"You may only use entry() during the 'building' stage");
1331
1332 if (row >= n)
1333 DUNE_THROW(BCRSMatrixError,"row index exceeds matrix size");
1334 if (col >= m)
1335 DUNE_THROW(BCRSMatrixError,"column index exceeds matrix size");
1336#endif
1337
1338 size_type* begin = r[row].getindexptr();
1339 size_type* end = begin + r[row].getsize();
1340
1341 size_type* pos = std::find(begin, end, col);
1342
1343 //treat the case that there was a match in the array
1344 if (pos != end)
1345 if (*pos == col)
1346 {
1347 std::ptrdiff_t offset = pos - r[row].getindexptr();
1348 B* aptr = r[row].getptr() + offset;
1349
1350 return *aptr;
1351 }
1352
1353 //determine whether overflow has to be taken into account or not
1354 if (r[row].getsize() == avg)
1355 return overflow[std::make_pair(row,col)];
1356 else
1357 {
1358 //modify index array
1359 *end = col;
1360
1361 //do simultaneous operations on data array a
1362 std::ptrdiff_t offset = end - r[row].getindexptr();
1363 B* apos = r[row].getptr() + offset;
1364
1365 //increase rowsize
1366 r[row].setsize(r[row].getsize()+1);
1367
1368 //return reference to the newly created entry
1369 return *apos;
1370 }
1371 }
1372
1374
1385 {
1386 if (build_mode!=implicit)
1387 DUNE_THROW(BCRSMatrixError,"requires implicit build mode");
1388 if (ready==built)
1389 DUNE_THROW(BCRSMatrixError,"matrix already built up, no more need for compression");
1390 if (ready==notAllocated)
1391 DUNE_THROW(BCRSMatrixError,"matrix size not set and no memory allocated yet");
1392 if (ready!=building)
1393 DUNE_THROW(InvalidStateException,"You may only call compress() at the end of the 'building' stage");
1394
1395 //calculate statistics
1397 stats.overflow_total = overflow.size();
1398 stats.maximum = 0;
1399
1400 //get insertion iterators pointing to one before start (for later use of ++it)
1401 size_type* jiit = j_.get();
1402 B* aiit = a;
1403
1404 //get iterator to the smallest overflow element
1405 typename OverflowType::iterator oit = overflow.begin();
1406
1407 //store a copy of index pointers on which to perform sorting
1408 std::vector<size_type*> perm;
1409
1410 //iterate over all rows and copy elements into their position in the compressed array
1411 for (size_type i=0; i<n; i++)
1412 {
1413 //get old pointers into a and j and size without overflow changes
1414 size_type* begin = r[i].getindexptr();
1415 //B* apos = r[i].getptr();
1416 size_type size = r[i].getsize();
1417
1418 perm.resize(size);
1419
1420 typename std::vector<size_type*>::iterator it = perm.begin();
1421 for (size_type* iit = begin; iit < begin + size; ++iit, ++it)
1422 *it = iit;
1423
1424 //sort permutation array
1425 std::sort(perm.begin(),perm.end(),PointerCompare<size_type>());
1426
1427 //change row window pointer to their new positions
1428 r[i].setindexptr(jiit);
1429 r[i].setptr(aiit);
1430
1431 for (it = perm.begin(); it != perm.end(); ++it)
1432 {
1433 //check whether there are elements in the overflow area which take precedence
1434 while ((oit!=overflow.end()) && (oit->first < std::make_pair(i,**it)))
1435 {
1436 //check whether there is enough memory to write to
1437 if (jiit > begin)
1439 "Allocated memory for BCRSMatrix exhausted during compress()!"
1440 "Please increase either the average number of entries per row or the compressionBufferSize value."
1441 );
1442 //copy an element from the overflow area to the insertion position in a and j
1443 *jiit = oit->first.second;
1444 ++jiit;
1445 *aiit = oit->second;
1446 ++aiit;
1447 ++oit;
1448 r[i].setsize(r[i].getsize()+1);
1449 }
1450
1451 //check whether there is enough memory to write to
1452 if (jiit > begin)
1454 "Allocated memory for BCRSMatrix exhausted during compress()!"
1455 "Please increase either the average number of entries per row or the compressionBufferSize value."
1456 );
1457
1458 //copy element from array
1459 *jiit = **it;
1460 ++jiit;
1461 B* apos = *it - j_.get() + a;
1462 *aiit = *apos;
1463 ++aiit;
1464 }
1465
1466 //copy remaining elements from the overflow area
1467 while ((oit!=overflow.end()) && (oit->first.first == i))
1468 {
1469 //check whether there is enough memory to write to
1470 if (jiit > begin)
1472 "Allocated memory for BCRSMatrix exhausted during compress()!"
1473 "Please increase either the average number of entries per row or the compressionBufferSize value."
1474 );
1475
1476 //copy and element from the overflow area to the insertion position in a and j
1477 *jiit = oit->first.second;
1478 ++jiit;
1479 *aiit = oit->second;
1480 ++aiit;
1481 ++oit;
1482 r[i].setsize(r[i].getsize()+1);
1483 }
1484
1485 // update maximum row size
1486 if (r[i].getsize()>stats.maximum)
1487 stats.maximum = r[i].getsize();
1488 }
1489
1490 // overflow area may be cleared
1491 overflow.clear();
1492
1493 //determine average number of entries and memory usage
1494 if ( n == 0)
1495 {
1496 stats.avg = 0;
1497 stats.mem_ratio = 1;
1498 }
1499 else
1500 {
1501 std::ptrdiff_t diff = (r[n-1].getindexptr() + r[n-1].getsize() - j_.get());
1502 nnz_ = diff;
1503 stats.avg = (double) (nnz_) / (double) n;
1504 stats.mem_ratio = (double) (nnz_) / (double) allocationSize_;
1505 }
1506
1507 //matrix is now built
1508 ready = built;
1509
1510 return stats;
1511 }
1512
1513 //===== vector space arithmetic
1514
1517 {
1518#ifdef DUNE_ISTL_WITH_CHECKING
1519 if (ready != built)
1520 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1521#endif
1522
1523 if (nnz_ > 0)
1524 {
1525 // process 1D array
1526 for (size_type i=0; i<nnz_; i++)
1527 a[i] *= k;
1528 }
1529 else
1530 {
1531 RowIterator endi=end();
1532 for (RowIterator i=begin(); i!=endi; ++i)
1533 {
1534 ColIterator endj = (*i).end();
1535 for (ColIterator j=(*i).begin(); j!=endj; ++j)
1536 (*j) *= k;
1537 }
1538 }
1539
1540 return *this;
1541 }
1542
1545 {
1546#ifdef DUNE_ISTL_WITH_CHECKING
1547 if (ready != built)
1548 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1549#endif
1550
1551 if (nnz_ > 0)
1552 {
1553 // process 1D array
1554 for (size_type i=0; i<nnz_; i++)
1555 a[i] /= k;
1556 }
1557 else
1558 {
1559 RowIterator endi=end();
1560 for (RowIterator i=begin(); i!=endi; ++i)
1561 {
1562 ColIterator endj = (*i).end();
1563 for (ColIterator j=(*i).begin(); j!=endj; ++j)
1564 (*j) /= k;
1565 }
1566 }
1567
1568 return *this;
1569 }
1570
1571
1578 {
1579#ifdef DUNE_ISTL_WITH_CHECKING
1580 if (ready != built || b.ready != built)
1581 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1582 if(N()!=b.N() || M() != b.M())
1583 DUNE_THROW(RangeError, "Matrix sizes do not match!");
1584#endif
1585 RowIterator endi=end();
1586 ConstRowIterator j=b.begin();
1587 for (RowIterator i=begin(); i!=endi; ++i, ++j) {
1588 i->operator+=(*j);
1589 }
1590
1591 return *this;
1592 }
1593
1600 {
1601#ifdef DUNE_ISTL_WITH_CHECKING
1602 if (ready != built || b.ready != built)
1603 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1604 if(N()!=b.N() || M() != b.M())
1605 DUNE_THROW(RangeError, "Matrix sizes do not match!");
1606#endif
1607 RowIterator endi=end();
1608 ConstRowIterator j=b.begin();
1609 for (RowIterator i=begin(); i!=endi; ++i, ++j) {
1610 i->operator-=(*j);
1611 }
1612
1613 return *this;
1614 }
1615
1625 {
1626#ifdef DUNE_ISTL_WITH_CHECKING
1627 if (ready != built || b.ready != built)
1628 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1629 if(N()!=b.N() || M() != b.M())
1630 DUNE_THROW(RangeError, "Matrix sizes do not match!");
1631#endif
1632 RowIterator endi=end();
1633 ConstRowIterator j=b.begin();
1634 for(RowIterator i=begin(); i!=endi; ++i, ++j)
1635 i->axpy(alpha, *j);
1636
1637 return *this;
1638 }
1639
1640 //===== linear maps
1641
1643 template<class X, class Y>
1644 void mv (const X& x, Y& y) const
1645 {
1646#ifdef DUNE_ISTL_WITH_CHECKING
1647 if (ready != built)
1648 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1649 if (x.N()!=M()) DUNE_THROW(BCRSMatrixError,
1650 "Size mismatch: M: " << N() << "x" << M() << " x: " << x.N());
1651 if (y.N()!=N()) DUNE_THROW(BCRSMatrixError,
1652 "Size mismatch: M: " << N() << "x" << M() << " y: " << y.N());
1653#endif
1654 ConstRowIterator endi=end();
1655 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1656 {
1657 y[i.index()]=0;
1658 ConstColIterator endj = (*i).end();
1659 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1660 {
1661 auto&& xj = Impl::asVector(x[j.index()]);
1662 auto&& yi = Impl::asVector(y[i.index()]);
1663 Impl::asMatrix(*j).umv(xj, yi);
1664 }
1665 }
1666 }
1667
1669 template<class X, class Y>
1670 void umv (const X& x, Y& y) const
1671 {
1672#ifdef DUNE_ISTL_WITH_CHECKING
1673 if (ready != built)
1674 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1675 if (x.N()!=M()) DUNE_THROW(BCRSMatrixError,"index out of range");
1676 if (y.N()!=N()) DUNE_THROW(BCRSMatrixError,"index out of range");
1677#endif
1678 ConstRowIterator endi=end();
1679 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1680 {
1681 ConstColIterator endj = (*i).end();
1682 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1683 {
1684 auto&& xj = Impl::asVector(x[j.index()]);
1685 auto&& yi = Impl::asVector(y[i.index()]);
1686 Impl::asMatrix(*j).umv(xj,yi);
1687 }
1688 }
1689 }
1690
1692 template<class X, class Y>
1693 void mmv (const X& x, Y& y) const
1694 {
1695#ifdef DUNE_ISTL_WITH_CHECKING
1696 if (ready != built)
1697 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1698 if (x.N()!=M()) DUNE_THROW(BCRSMatrixError,"index out of range");
1699 if (y.N()!=N()) DUNE_THROW(BCRSMatrixError,"index out of range");
1700#endif
1701 ConstRowIterator endi=end();
1702 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1703 {
1704 ConstColIterator endj = (*i).end();
1705 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1706 {
1707 auto&& xj = Impl::asVector(x[j.index()]);
1708 auto&& yi = Impl::asVector(y[i.index()]);
1709 Impl::asMatrix(*j).mmv(xj,yi);
1710 }
1711 }
1712 }
1713
1715 template<class X, class Y, class F>
1716 void usmv (F&& alpha, const X& x, Y& y) const
1717 {
1718#ifdef DUNE_ISTL_WITH_CHECKING
1719 if (ready != built)
1720 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1721 if (x.N()!=M()) DUNE_THROW(BCRSMatrixError,"index out of range");
1722 if (y.N()!=N()) DUNE_THROW(BCRSMatrixError,"index out of range");
1723#endif
1724 ConstRowIterator endi=end();
1725 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1726 {
1727 ConstColIterator endj = (*i).end();
1728 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1729 {
1730 auto&& xj = Impl::asVector(x[j.index()]);
1731 auto&& yi = Impl::asVector(y[i.index()]);
1732 Impl::asMatrix(*j).usmv(alpha,xj,yi);
1733 }
1734 }
1735 }
1736
1738 template<class X, class Y>
1739 void mtv (const X& x, Y& y) const
1740 {
1741#ifdef DUNE_ISTL_WITH_CHECKING
1742 if (ready != built)
1743 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1744 if (x.N()!=N()) DUNE_THROW(BCRSMatrixError,"index out of range");
1745 if (y.N()!=M()) DUNE_THROW(BCRSMatrixError,"index out of range");
1746#endif
1747 for(size_type i=0; i<y.N(); ++i)
1748 y[i]=0;
1749 umtv(x,y);
1750 }
1751
1753 template<class X, class Y>
1754 void umtv (const X& x, Y& y) const
1755 {
1756#ifdef DUNE_ISTL_WITH_CHECKING
1757 if (ready != built)
1758 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1759 if (x.N()!=N()) DUNE_THROW(BCRSMatrixError,"index out of range");
1760 if (y.N()!=M()) DUNE_THROW(BCRSMatrixError,"index out of range");
1761#endif
1762 ConstRowIterator endi=end();
1763 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1764 {
1765 ConstColIterator endj = (*i).end();
1766 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1767 {
1768 auto&& xi = Impl::asVector(x[i.index()]);
1769 auto&& yj = Impl::asVector(y[j.index()]);
1770 Impl::asMatrix(*j).umtv(xi,yj);
1771 }
1772 }
1773 }
1774
1776 template<class X, class Y>
1777 void mmtv (const X& x, Y& y) const
1778 {
1779#ifdef DUNE_ISTL_WITH_CHECKING
1780 if (x.N()!=N()) DUNE_THROW(BCRSMatrixError,"index out of range");
1781 if (y.N()!=M()) DUNE_THROW(BCRSMatrixError,"index out of range");
1782#endif
1783 ConstRowIterator endi=end();
1784 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1785 {
1786 ConstColIterator endj = (*i).end();
1787 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1788 {
1789 auto&& xi = Impl::asVector(x[i.index()]);
1790 auto&& yj = Impl::asVector(y[j.index()]);
1791 Impl::asMatrix(*j).mmtv(xi,yj);
1792 }
1793 }
1794 }
1795
1797 template<class X, class Y>
1798 void usmtv (const field_type& alpha, const X& x, Y& y) const
1799 {
1800#ifdef DUNE_ISTL_WITH_CHECKING
1801 if (ready != built)
1802 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1803 if (x.N()!=N()) DUNE_THROW(BCRSMatrixError,"index out of range");
1804 if (y.N()!=M()) DUNE_THROW(BCRSMatrixError,"index out of range");
1805#endif
1806 ConstRowIterator endi=end();
1807 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1808 {
1809 ConstColIterator endj = (*i).end();
1810 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1811 {
1812 auto&& xi = Impl::asVector(x[i.index()]);
1813 auto&& yj = Impl::asVector(y[j.index()]);
1814 Impl::asMatrix(*j).usmtv(alpha,xi,yj);
1815 }
1816 }
1817 }
1818
1820 template<class X, class Y>
1821 void umhv (const X& x, Y& y) const
1822 {
1823#ifdef DUNE_ISTL_WITH_CHECKING
1824 if (ready != built)
1825 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1826 if (x.N()!=N()) DUNE_THROW(BCRSMatrixError,"index out of range");
1827 if (y.N()!=M()) DUNE_THROW(BCRSMatrixError,"index out of range");
1828#endif
1829 ConstRowIterator endi=end();
1830 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1831 {
1832 ConstColIterator endj = (*i).end();
1833 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1834 {
1835 auto&& xi = Impl::asVector(x[i.index()]);
1836 auto&& yj = Impl::asVector(y[j.index()]);
1837 Impl::asMatrix(*j).umhv(xi,yj);
1838 }
1839 }
1840 }
1841
1843 template<class X, class Y>
1844 void mmhv (const X& x, Y& y) const
1845 {
1846#ifdef DUNE_ISTL_WITH_CHECKING
1847 if (ready != built)
1848 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1849 if (x.N()!=N()) DUNE_THROW(BCRSMatrixError,"index out of range");
1850 if (y.N()!=M()) DUNE_THROW(BCRSMatrixError,"index out of range");
1851#endif
1852 ConstRowIterator endi=end();
1853 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1854 {
1855 ConstColIterator endj = (*i).end();
1856 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1857 {
1858 auto&& xi = Impl::asVector(x[i.index()]);
1859 auto&& yj = Impl::asVector(y[j.index()]);
1860 Impl::asMatrix(*j).mmhv(xi,yj);
1861 }
1862 }
1863 }
1864
1866 template<class X, class Y>
1867 void usmhv (const field_type& alpha, const X& x, Y& y) const
1868 {
1869#ifdef DUNE_ISTL_WITH_CHECKING
1870 if (ready != built)
1871 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1872 if (x.N()!=N()) DUNE_THROW(BCRSMatrixError,"index out of range");
1873 if (y.N()!=M()) DUNE_THROW(BCRSMatrixError,"index out of range");
1874#endif
1875 ConstRowIterator endi=end();
1876 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1877 {
1878 ConstColIterator endj = (*i).end();
1879 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1880 {
1881 auto&& xi = Impl::asVector(x[i.index()]);
1882 auto&& yj = Impl::asVector(y[j.index()]);
1883 Impl::asMatrix(*j).usmhv(alpha,xi,yj);
1884 }
1885 }
1886 }
1887
1888
1889 //===== norms
1890
1892 typename FieldTraits<field_type>::real_type frobenius_norm2 () const
1893 {
1894#ifdef DUNE_ISTL_WITH_CHECKING
1895 if (ready != built)
1896 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1897#endif
1898
1899 typename FieldTraits<field_type>::real_type sum=0;
1900
1901 for (auto&& row : (*this))
1902 for (auto&& entry : row)
1903 sum += Impl::asMatrix(entry).frobenius_norm2();
1904
1905 return sum;
1906 }
1907
1909 typename FieldTraits<field_type>::real_type frobenius_norm () const
1910 {
1911 return sqrt(frobenius_norm2());
1912 }
1913
1915 template <typename ft = field_type,
1916 typename std::enable_if<!HasNaN<ft>::value, int>::type = 0>
1917 typename FieldTraits<ft>::real_type infinity_norm() const {
1918 if (ready != built)
1919 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1920
1921 using real_type = typename FieldTraits<ft>::real_type;
1922 using std::max;
1923
1924 real_type norm = 0;
1925 for (auto const &x : *this) {
1926 real_type sum = 0;
1927 for (auto const &y : x)
1928 sum += Impl::asMatrix(y).infinity_norm();
1929 norm = max(sum, norm);
1930 }
1931 return norm;
1932 }
1933
1935 template <typename ft = field_type,
1936 typename std::enable_if<!HasNaN<ft>::value, int>::type = 0>
1937 typename FieldTraits<ft>::real_type infinity_norm_real() const {
1938 if (ready != built)
1939 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1940
1941 using real_type = typename FieldTraits<ft>::real_type;
1942 using std::max;
1943
1944 real_type norm = 0;
1945 for (auto const &x : *this) {
1946 real_type sum = 0;
1947 for (auto const &y : x)
1948 sum += Impl::asMatrix(y).infinity_norm_real();
1949 norm = max(sum, norm);
1950 }
1951 return norm;
1952 }
1953
1955 template <typename ft = field_type,
1956 typename std::enable_if<HasNaN<ft>::value, int>::type = 0>
1957 typename FieldTraits<ft>::real_type infinity_norm() const {
1958 if (ready != built)
1959 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1960
1961 using real_type = typename FieldTraits<ft>::real_type;
1962 using std::max;
1963
1964 real_type norm = 0;
1965 real_type isNaN = 1;
1966 for (auto const &x : *this) {
1967 real_type sum = 0;
1968 for (auto const &y : x)
1969 sum += Impl::asMatrix(y).infinity_norm();
1970 norm = max(sum, norm);
1971 isNaN += sum;
1972 }
1973
1974 return norm * (isNaN / isNaN);
1975 }
1976
1978 template <typename ft = field_type,
1979 typename std::enable_if<HasNaN<ft>::value, int>::type = 0>
1980 typename FieldTraits<ft>::real_type infinity_norm_real() const {
1981 if (ready != built)
1982 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1983
1984 using real_type = typename FieldTraits<ft>::real_type;
1985 using std::max;
1986
1987 real_type norm = 0;
1988 real_type isNaN = 1;
1989
1990 for (auto const &x : *this) {
1991 real_type sum = 0;
1992 for (auto const &y : x)
1993 sum += Impl::asMatrix(y).infinity_norm_real();
1994 norm = max(sum, norm);
1995 isNaN += sum;
1996 }
1997
1998 return norm * (isNaN / isNaN);
1999 }
2000
2001 //===== sizes
2002
2004 size_type N () const
2005 {
2006 return n;
2007 }
2008
2010 size_type M () const
2011 {
2012 return m;
2013 }
2014
2017 {
2018 // in case of row-wise allocation
2019 if( nnz_ <= 0 )
2020 nnz_ = std::accumulate( begin(), end(), size_type( 0 ), [] ( size_type s, const row_type &row ) { return s+row.getsize(); } );
2021 return nnz_;
2022 }
2023
2026 {
2027 return ready;
2028 }
2029
2032 {
2033 return build_mode;
2034 }
2035
2036 //===== query
2037
2039 bool exists (size_type i, size_type j) const
2040 {
2041#ifdef DUNE_ISTL_WITH_CHECKING
2042 if (i<0 || i>=n) DUNE_THROW(BCRSMatrixError,"row index out of range");
2043 if (j<0 || j>=m) DUNE_THROW(BCRSMatrixError,"column index out of range");
2044#endif
2045 return (r[i].size() && r[i].find(j) != r[i].end());
2046 }
2047
2048
2049 protected:
2050 // state information
2051 BuildMode build_mode; // row wise or whole matrix
2052 BuildStage ready; // indicate the stage the matrix building is in
2053
2054 // The allocator used for memory management
2055 typename std::allocator_traits<A>::template rebind_alloc<B> allocator_;
2056
2057 typename std::allocator_traits<A>::template rebind_alloc<row_type> rowAllocator_;
2058
2059 typename std::allocator_traits<A>::template rebind_alloc<size_type> sizeAllocator_;
2060
2061 // size of the matrix
2062 size_type n; // number of rows
2063 size_type m; // number of columns
2064 mutable size_type nnz_; // number of nonzeroes contained in the matrix
2065 size_type allocationSize_; //allocated size of a and j arrays, except in implicit mode: nnz_==allocationsSize_
2066 // zero means that memory is allocated separately for each row.
2067
2068 // the rows are dynamically allocated
2069 row_type* r; // [n] the individual rows having pointers into a,j arrays
2070
2071 // dynamically allocated memory
2072 B* a; // [allocationSize] non-zero entries of the matrix in row-wise ordering
2073 // If a single array of column indices is used, it can be shared
2074 // between different matrices with the same sparsity pattern
2075 std::shared_ptr<size_type> j_; // [allocationSize] column indices of entries
2076
2077 // additional data is needed in implicit buildmode
2078 size_type avg;
2079 double compressionBufferSize_;
2080
2081 typedef std::map<std::pair<size_type,size_type>, B> OverflowType;
2082 OverflowType overflow;
2083
2084 void setWindowPointers(ConstRowIterator row)
2085 {
2086 row_type current_row(a,j_.get(),0); // Pointers to current row data
2087 for (size_type i=0; i<n; i++, ++row) {
2088 // set row i
2089 size_type s = row->getsize();
2090
2091 if (s>0) {
2092 // setup pointers and size
2093 r[i].set(s,current_row.getptr(), current_row.getindexptr());
2094 // update pointer for next row
2095 current_row.setptr(current_row.getptr()+s);
2096 current_row.setindexptr(current_row.getindexptr()+s);
2097 } else{
2098 // empty row
2099 r[i].set(0,nullptr,nullptr);
2100 }
2101 }
2102 }
2103
2105
2110 {
2111 size_type* jptr = j_.get();
2112 for (size_type i=0; i<n; ++i, ++row) {
2113 // set row i
2114 size_type s = row->getsize();
2115
2116 if (s>0) {
2117 // setup pointers and size
2118 r[i].setsize(s);
2119 r[i].setindexptr(jptr);
2120 } else{
2121 // empty row
2122 r[i].set(0,nullptr,nullptr);
2123 }
2124
2125 // advance position in global array
2126 jptr += s;
2127 }
2128 }
2129
2131
2136 {
2137 B* aptr = a;
2138 for (size_type i=0; i<n; ++i) {
2139 // set row i
2140 if (r[i].getsize() > 0) {
2141 // setup pointers and size
2142 r[i].setptr(aptr);
2143 } else{
2144 // empty row
2145 r[i].set(0,nullptr,nullptr);
2146 }
2147
2148 // advance position in global array
2149 aptr += r[i].getsize();
2150 }
2151 }
2152
2155 {
2156 setWindowPointers(Mat.begin());
2157
2158 // copy data
2159 for (size_type i=0; i<n; i++) r[i] = Mat.r[i];
2160
2161 // finish off
2162 build_mode = row_wise; // dummy
2163 ready = built;
2164 }
2165
2171 void deallocate(bool deallocateRows=true)
2172 {
2173
2174 if (notAllocated)
2175 return;
2176
2177 if (allocationSize_>0)
2178 {
2179 // a,j_ have been allocated as one long vector
2180 j_.reset();
2181 if (a)
2182 {
2183 for(B *aiter=a+(allocationSize_-1), *aend=a-1; aiter!=aend; --aiter)
2184 std::allocator_traits<decltype(allocator_)>::destroy(allocator_, aiter);
2185 allocator_.deallocate(a,allocationSize_);
2186 a = nullptr;
2187 }
2188 }
2189 else if (r)
2190 {
2191 // check if memory for rows have been allocated individually
2192 for (size_type i=0; i<n; i++)
2193 if (r[i].getsize()>0)
2194 {
2195 for (B *col=r[i].getptr()+(r[i].getsize()-1),
2196 *colend = r[i].getptr()-1; col!=colend; --col) {
2197 std::allocator_traits<decltype(allocator_)>::destroy(allocator_, col);
2198 }
2199 sizeAllocator_.deallocate(r[i].getindexptr(),1);
2200 allocator_.deallocate(r[i].getptr(),1);
2201 // clear out row data in case we don't want to deallocate the rows
2202 // otherwise we might run into a double free problem here later
2203 r[i].set(0,nullptr,nullptr);
2204 }
2205 }
2206
2207 // deallocate the rows
2208 if (n>0 && deallocateRows && r) {
2209 for(row_type *riter=r+(n-1), *rend=r-1; riter!=rend; --riter)
2210 std::allocator_traits<decltype(rowAllocator_)>::destroy(rowAllocator_, riter);
2211 rowAllocator_.deallocate(r,n);
2212 r = nullptr;
2213 }
2214
2215 // Mark matrix as not built at all.
2216 ready=notAllocated;
2217
2218 }
2219
2238 void allocate(size_type rows, size_type columns, size_type allocationSize, bool allocateRows, bool allocate_data)
2239 {
2240 // Store size
2241 n = rows;
2242 m = columns;
2243 nnz_ = allocationSize;
2244 allocationSize_ = allocationSize;
2245
2246 // allocate rows
2247 if(allocateRows) {
2248 if (n>0) {
2249 if (r)
2250 DUNE_THROW(InvalidStateException,"Rows have already been allocated, cannot allocate a second time");
2251 r = rowAllocator_.allocate(rows);
2252 // initialize row entries
2253 for(row_type* ri=r; ri!=r+rows; ++ri)
2254 std::allocator_traits<decltype(rowAllocator_)>::construct(rowAllocator_, ri, row_type());
2255 }else{
2256 r = 0;
2257 }
2258 }
2259
2260 // allocate a and j_ array
2261 if (allocate_data)
2262 allocateData();
2263 // allocate column indices only if not yet present (enable sharing)
2264 if (allocationSize_>0) {
2265 // we copy allocator and size to the deleter since _j may outlive this class
2266 if (!j_.get())
2267 j_.reset(sizeAllocator_.allocate(allocationSize_),
2268 [alloc = sizeAllocator_, size = allocationSize_](auto ptr) mutable {
2269 alloc.deallocate(ptr, size);
2270 });
2271 }else{
2272 j_.reset();
2273 }
2274
2275 // Mark the matrix as not built.
2276 ready = building;
2277 }
2278
2279 void allocateData()
2280 {
2281 if (a)
2282 DUNE_THROW(InvalidStateException,"Cannot allocate data array (already allocated)");
2283 if (allocationSize_>0) {
2284 a = allocator_.allocate(allocationSize_);
2285 // use placement new to call constructor that allocates
2286 // additional memory.
2287 new (a) B[allocationSize_];
2288 } else {
2289 a = nullptr;
2290 }
2291 }
2292
2299 {
2300 if (build_mode != implicit)
2301 DUNE_THROW(InvalidStateException,"implicit_allocate() may only be called in implicit build mode");
2302 if (ready != notAllocated)
2303 DUNE_THROW(InvalidStateException,"memory has already been allocated");
2304
2305 // check to make sure the user has actually set the parameters
2306 if (compressionBufferSize_ < 0)
2307 DUNE_THROW(InvalidStateException,"You have to set the implicit build mode parameters before starting to build the matrix");
2308 //calculate size of overflow region, add buffer for row sort!
2309 size_type osize = (size_type) (_n*avg)*compressionBufferSize_ + 4*avg;
2310 allocationSize_ = _n*avg + osize;
2311
2312 allocate(_n, _m, allocationSize_,true,true);
2313
2314 //set row pointers correctly
2315 size_type* jptr = j_.get() + osize;
2316 B* aptr = a + osize;
2317 for (size_type i=0; i<n; i++)
2318 {
2319 r[i].set(0,aptr,jptr);
2320 jptr = jptr + avg;
2321 aptr = aptr + avg;
2322 }
2323
2324 ready = building;
2325 }
2326 };
2327
2328
2329 template<class B, class A>
2330 struct FieldTraits< BCRSMatrix<B, A> >
2331 {
2332 using field_type = typename BCRSMatrix<B, A>::field_type;
2333 using real_type = typename FieldTraits<field_type>::real_type;
2334 };
2335
2338} // end namespace
2339
2340#endif
Helper functions for determining the vector/matrix block level.
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:24
Iterator class for sequential creation of blocks
Definition: bcrsmatrix.hh:955
bool operator==(const CreateIterator &it) const
equality
Definition: bcrsmatrix.hh:1052
CreateIterator & operator++()
prefix increment
Definition: bcrsmatrix.hh:975
size_type index() const
The number of the row that the iterator currently points to.
Definition: bcrsmatrix.hh:1058
bool operator!=(const CreateIterator &it) const
inequality
Definition: bcrsmatrix.hh:1046
CreateIterator(BCRSMatrix &_Mat, size_type _i)
constructor
Definition: bcrsmatrix.hh:958
void insert(size_type j)
put column index in row
Definition: bcrsmatrix.hh:1064
bool contains(size_type j)
return true if column index is in row
Definition: bcrsmatrix.hh:1070
size_type size() const
Get the current row size.
Definition: bcrsmatrix.hh:1079
Iterator access to matrix rows
Definition: bcrsmatrix.hh:576
RealRowIterator()
empty constructor, use with care!
Definition: bcrsmatrix.hh:593
bool equals(const RealRowIterator< ValueType > &other) const
equality
Definition: bcrsmatrix.hh:621
std::remove_const< T >::type ValueType
The unqualified value type.
Definition: bcrsmatrix.hh:580
bool equals(const RealRowIterator< const ValueType > &other) const
equality
Definition: bcrsmatrix.hh:628
RealRowIterator(row_type *_p, size_type _i)
constructor
Definition: bcrsmatrix.hh:588
size_type index() const
return index
Definition: bcrsmatrix.hh:603
A sparse block matrix with compressed row storage.
Definition: bcrsmatrix.hh:467
typename Imp::BlockTraits< B >::field_type field_type
export the type representing the field
Definition: bcrsmatrix.hh:489
bool exists(size_type i, size_type j) const
return true if (i,j) is in pattern
Definition: bcrsmatrix.hh:2039
BuildStage buildStage() const
The current build stage of the matrix.
Definition: bcrsmatrix.hh:2025
Iterator begin()
Get iterator to first row.
Definition: bcrsmatrix.hh:672
friend class CreateIterator
allow CreateIterator to access internal data
Definition: bcrsmatrix.hh:1094
void copyWindowStructure(const BCRSMatrix &Mat)
Copy the window structure from another matrix.
Definition: bcrsmatrix.hh:2154
B & entry(size_type row, size_type col)
Returns reference to entry (row,col) of the matrix.
Definition: bcrsmatrix.hh:1320
void mmhv(const X &x, Y &y) const
y -= A^H x
Definition: bcrsmatrix.hh:1844
void usmhv(const field_type &alpha, const X &x, Y &y) const
y += alpha A^H x
Definition: bcrsmatrix.hh:1867
void umtv(const X &x, Y &y) const
y += A^T x
Definition: bcrsmatrix.hh:1754
RealRowIterator< const row_type > const_iterator
The const iterator over the matrix rows.
Definition: bcrsmatrix.hh:704
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:2238
BCRSMatrix & axpy(field_type alpha, const BCRSMatrix &b)
Add the scaled entries of another matrix to this one.
Definition: bcrsmatrix.hh:1624
FieldTraits< ft >::real_type infinity_norm_real() const
simplified infinity norm (uses Manhattan norm for complex values)
Definition: bcrsmatrix.hh:1937
~BCRSMatrix()
destructor
Definition: bcrsmatrix.hh:822
void deallocate(bool deallocateRows=true)
deallocate memory of the matrix.
Definition: bcrsmatrix.hh:2171
Iterator RowIterator
rename the iterators for easier access
Definition: bcrsmatrix.hh:698
row_type & operator[](size_type i)
random access to the rows
Definition: bcrsmatrix.hh:546
BCRSMatrix()
an empty matrix
Definition: bcrsmatrix.hh:746
void endrowsizes()
indicate that size of all rows is defined
Definition: bcrsmatrix.hh:1149
void incrementrowsize(size_type i, size_type s=1)
increment size of row i by s (1 by default)
Definition: bcrsmatrix.hh:1138
void mtv(const X &x, Y &y) const
y = A^T x
Definition: bcrsmatrix.hh:1739
void umhv(const X &x, Y &y) const
y += A^H x
Definition: bcrsmatrix.hh:1821
size_type nonzeroes() const
number of blocks that are stored (the number of blocks that possibly are nonzero)
Definition: bcrsmatrix.hh:2016
ConstIterator ConstRowIterator
rename the const row iterator for easier access
Definition: bcrsmatrix.hh:735
void setIndicesNoSort(size_type row, It begin, It end)
Set all column indices for row from the given iterator range.
Definition: bcrsmatrix.hh:1235
void setrowsize(size_type i, size_type s)
Set number of indices in row i to s.
Definition: bcrsmatrix.hh:1117
BCRSMatrix & operator*=(const field_type &k)
vector space multiplication with scalar
Definition: bcrsmatrix.hh:1516
RealRowIterator< row_type > iterator
The iterator over the (mutable matrix rows.
Definition: bcrsmatrix.hh:668
void usmtv(const field_type &alpha, const X &x, Y &y) const
y += alpha A^T x
Definition: bcrsmatrix.hh:1798
ConstIterator beforeBegin() const
Definition: bcrsmatrix.hh:729
Iterator beforeBegin()
Definition: bcrsmatrix.hh:692
BuildMode
we support two modes
Definition: bcrsmatrix.hh:507
@ implicit
Build entries randomly with an educated guess for the number of entries per row.
Definition: bcrsmatrix.hh:536
@ unknown
Build mode not set!
Definition: bcrsmatrix.hh:540
@ random
Build entries randomly.
Definition: bcrsmatrix.hh:527
@ row_wise
Build in a row-wise manner.
Definition: bcrsmatrix.hh:518
BCRSMatrix(size_type _n, size_type _m, size_type _nnz, BuildMode bm)
matrix with known number of nonzeroes
Definition: bcrsmatrix.hh:753
Imp::CompressedBlockVectorWindow< B, size_type > row_type
implement row_type with compressed vector
Definition: bcrsmatrix.hh:501
::Dune::CompressionStatistics< size_type > CompressionStatistics
The type for the statistics object returned by compress()
Definition: bcrsmatrix.hh:504
BCRSMatrix & operator-=(const BCRSMatrix &b)
Subtract the entries of another matrix from this one.
Definition: bcrsmatrix.hh:1599
BCRSMatrix(const BCRSMatrix &Mat)
copy constructor
Definition: bcrsmatrix.hh:803
Iterator end()
Get iterator to one beyond last row.
Definition: bcrsmatrix.hh:678
void setIndices(size_type row, It begin, It end)
Set all column indices for row from the given iterator range.
Definition: bcrsmatrix.hh:1258
void addindex(size_type row, size_type col)
add index (row,col) to the matrix
Definition: bcrsmatrix.hh:1191
row_type::Iterator ColIterator
Iterator for the entries of each row.
Definition: bcrsmatrix.hh:701
FieldTraits< field_type >::real_type frobenius_norm() const
frobenius norm: sqrt(sum over squared values of entries)
Definition: bcrsmatrix.hh:1909
A::size_type size_type
The type for the index access and the size.
Definition: bcrsmatrix.hh:498
BCRSMatrix & operator/=(const field_type &k)
vector space division by scalar
Definition: bcrsmatrix.hh:1544
BCRSMatrix & operator+=(const BCRSMatrix &b)
Add the entries of another matrix to this one.
Definition: bcrsmatrix.hh:1577
BCRSMatrix(size_type _n, size_type _m, size_type _avg, double compressionBufferSize, BuildMode bm)
construct matrix with a known average number of entries per row
Definition: bcrsmatrix.hh:782
CreateIterator createend()
get create iterator pointing to one after the last block
Definition: bcrsmatrix.hh:1103
FieldTraits< field_type >::real_type frobenius_norm2() const
square of frobenius norm, need for block recursion
Definition: bcrsmatrix.hh:1892
Iterator beforeEnd()
Definition: bcrsmatrix.hh:685
row_type::ConstIterator ConstColIterator
Const iterator to the entries of a row.
Definition: bcrsmatrix.hh:738
void usmv(F &&alpha, const X &x, Y &y) const
y += alpha A x
Definition: bcrsmatrix.hh:1716
size_type getrowsize(size_type i) const
get current number of indices in row i
Definition: bcrsmatrix.hh:1128
size_type M() const
number of columns (counted in blocks)
Definition: bcrsmatrix.hh:2010
CreateIterator createbegin()
get initial create iterator
Definition: bcrsmatrix.hh:1097
BuildStage
Definition: bcrsmatrix.hh:470
@ rowSizesBuilt
The row sizes of the matrix are known.
Definition: bcrsmatrix.hh:481
@ built
The matrix structure is fully built.
Definition: bcrsmatrix.hh:483
@ notbuilt
Matrix is not built at all, no memory has been allocated, build mode and size can still be set.
Definition: bcrsmatrix.hh:472
@ notAllocated
Matrix is not built at all, no memory has been allocated, build mode and size can still be set.
Definition: bcrsmatrix.hh:474
@ building
Matrix is currently being built, some memory has been allocated, build mode and size are fixed.
Definition: bcrsmatrix.hh:476
BuildMode buildMode() const
The currently selected build mode of the matrix.
Definition: bcrsmatrix.hh:2031
void mmv(const X &x, Y &y) const
y -= A x
Definition: bcrsmatrix.hh:1693
FieldTraits< ft >::real_type infinity_norm() const
infinity norm (row sum norm, how to generalize for blocks?)
Definition: bcrsmatrix.hh:1917
void mv(const X &x, Y &y) const
y = A x
Definition: bcrsmatrix.hh:1644
B block_type
export the type representing the components
Definition: bcrsmatrix.hh:492
void mmtv(const X &x, Y &y) const
y -= A^T x
Definition: bcrsmatrix.hh:1777
void umv(const X &x, Y &y) const
y += A x
Definition: bcrsmatrix.hh:1670
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:2298
void setImplicitBuildModeParameters(size_type _avg, double compressionBufferSize)
Set parameters needed for creation in implicit build mode.
Definition: bcrsmatrix.hh:887
BCRSMatrix(size_type _n, size_type _m, BuildMode bm)
matrix with unknown number of nonzeroes
Definition: bcrsmatrix.hh:762
void endindices()
indicate that all indices are defined, check consistency
Definition: bcrsmatrix.hh:1272
CompressionStatistics compress()
Finishes the buildstage in implicit mode.
Definition: bcrsmatrix.hh:1384
void setDataPointers()
Set data pointers for all rows.
Definition: bcrsmatrix.hh:2135
size_type N() const
number of rows (counted in blocks)
Definition: bcrsmatrix.hh:2004
void setBuildMode(BuildMode bm)
Sets the build mode of the matrix.
Definition: bcrsmatrix.hh:831
void setSize(size_type rows, size_type columns, size_type nnz=0)
Set the size of the matrix.
Definition: bcrsmatrix.hh:859
BCRSMatrix & operator=(const BCRSMatrix &Mat)
assignment
Definition: bcrsmatrix.hh:909
void setColumnPointers(ConstRowIterator row)
Copy row sizes from iterator range starting at row and set column index pointers for all rows.
Definition: bcrsmatrix.hh:2109
ConstIterator end() const
Get const iterator to one beyond last row.
Definition: bcrsmatrix.hh:715
ConstIterator begin() const
Get const iterator to first row.
Definition: bcrsmatrix.hh:709
A allocator_type
export the allocator type
Definition: bcrsmatrix.hh:495
ConstIterator beforeEnd() const
Definition: bcrsmatrix.hh:722
Proxy row object for entry access.
Definition: bcrsmatrix.hh:138
block_type & operator[](size_type j) const
Returns entry in column j.
Definition: bcrsmatrix.hh:143
A wrapper for uniform access to the BCRSMatrix during and after the build stage in implicit build mod...
Definition: bcrsmatrix.hh:118
Matrix::block_type block_type
The block_type of the underlying matrix.
Definition: bcrsmatrix.hh:126
ImplicitMatrixBuilder(Matrix &m)
Creates an ImplicitMatrixBuilder for matrix m.
Definition: bcrsmatrix.hh:171
M_ Matrix
The underlying matrix.
Definition: bcrsmatrix.hh:123
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:195
size_type M() const
The number of columns in the matrix.
Definition: bcrsmatrix.hh:218
Matrix::size_type size_type
The size_type of the underlying matrix.
Definition: bcrsmatrix.hh:129
row_object operator[](size_type i) const
Returns a proxy for entries in row i.
Definition: bcrsmatrix.hh:206
size_type N() const
The number of rows in the matrix.
Definition: bcrsmatrix.hh:212
Thrown when the compression buffer used by the implicit BCRSMatrix construction is exhausted.
Definition: istlexception.hh:37
Default exception if a function was called while the object is not in a valid state for that function...
Definition: exceptions.hh:375
A generic dynamic dense matrix.
Definition: matrix.hh:561
A::size_type size_type
Type for indices and sizes.
Definition: matrix.hh:577
T block_type
Export the type representing the components.
Definition: matrix.hh:568
Base class for stl conformant forward iterators.
Definition: iteratorfacades.hh:435
Default exception class for range errors.
Definition: exceptions.hh:348
Type traits to determine the type of reals (when working with complex numbers)
#define DUNE_THROW(E,...)
Definition: exceptions.hh:314
constexpr auto max
Function object that returns the greater of the given values.
Definition: hybridutilities.hh:485
constexpr decltype(auto) elementAt(Container &&c, Index &&i)
Get element at given position from container.
Definition: hybridutilities.hh:127
constexpr T accumulate(Range &&range, T value, F &&f)
Accumulate values.
Definition: hybridutilities.hh:280
DWarnType dwarn(std::cerr)
Stream for warnings indicating problems.
Definition: stdstreams.hh:162
This file implements iterator facade classes for writing stl conformant iterators.
Some handy generic functions for ISTL matrices.
bool isNaN(const FieldVector< K, SIZE > &b, PriorityTag< 2 >, ADLTag)
Returns whether any entry is NaN.
Definition: fvector.hh:458
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
constexpr auto get(std::integer_sequence< T, II... >, std::integral_constant< std::size_t, pos >={})
Return the entry at position pos of the given sequence.
Definition: integersequence.hh:22
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:89
size_type overflow_total
total number of elements written to the overflow area during construction.
Definition: bcrsmatrix.hh:95
size_type maximum
maximum number of non-zeroes per row.
Definition: bcrsmatrix.hh:93
double avg
average number of non-zeroes per row.
Definition: bcrsmatrix.hh:91
double mem_ratio
fraction of wasted memory resulting from non-used overflow area.
Definition: bcrsmatrix.hh:100
Traits for type conversions and type information.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden & Uni Heidelberg  |  generated with Hugo v0.111.3 (Jul 8, 22:38, 2025)