DUNE-FEM (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#include <dune/common/std/no_unique_address.hh>
30
32
37namespace Dune {
38
78 template<typename M>
79 struct MatrixDimension;
80
82
88 template<typename size_type>
90 {
92 double avg;
94 size_type maximum;
96 size_type overflow_total;
98
101 double mem_ratio;
102 };
103
105
117 template<class M_>
119 {
120
121 public:
122
124 typedef M_ Matrix;
125
128
130 typedef typename Matrix::size_type size_type;
131
133
139 {
140
141 public:
142
145 {
146 return _m.entry(_i,j);
147 }
148
149#ifndef DOXYGEN
150
152 : _m(m)
153 , _i(i)
154 {}
155
156#endif
157
158 private:
159
160 Matrix& _m;
161 size_type _i;
162
163 };
164
166
173 : _m(m)
174 {
175 if (m.buildMode() != Matrix::implicit)
176 DUNE_THROW(BCRSMatrixError,"You can only create an ImplicitBuilder for a matrix in implicit build mode");
177 if (m.buildStage() != Matrix::building)
178 DUNE_THROW(BCRSMatrixError,"You can only create an ImplicitBuilder for a matrix with set size that has not been compressed() yet");
179 }
180
182
196 ImplicitMatrixBuilder(Matrix& m, size_type rows, size_type cols, size_type avg_cols_per_row, double overflow_fraction)
197 : _m(m)
198 {
199 if (m.buildStage() != Matrix::notAllocated)
200 DUNE_THROW(BCRSMatrixError,"You can only set up a matrix for this ImplicitBuilder if it has no memory allocated yet");
201 m.setBuildMode(Matrix::implicit);
202 m.setImplicitBuildModeParameters(avg_cols_per_row,overflow_fraction);
203 m.setSize(rows,cols);
204 }
205
208 {
209 return row_object(_m,i);
210 }
211
213 size_type N() const
214 {
215 return _m.N();
216 }
217
219 size_type M() const
220 {
221 return _m.M();
222 }
223
224 private:
225
226 Matrix& _m;
227
228 };
229
466 template<class B, class A=std::allocator<B> >
468 {
469 friend struct MatrixDimension<BCRSMatrix>;
470 public:
484 built=3
485 };
486
487 //===== type definitions and constants
488
490 using field_type = typename Imp::BlockTraits<B>::field_type;
491
493 typedef B block_type;
494
496 typedef typename std::allocator_traits<A>::size_type size_type;
497
499 typedef Imp::CompressedBlockVectorWindow<B,size_type> row_type;
500
502 using allocator_type = A;
503
505 typedef ::Dune::CompressionStatistics<size_type> CompressionStatistics;
506
541 unknown
542 };
543
544 //===== random access interface to rows of the matrix
545
548 {
549#ifdef DUNE_ISTL_WITH_CHECKING
550 if (build_mode == implicit && ready != built)
551 DUNE_THROW(BCRSMatrixError,"You cannot use operator[] in implicit build mode before calling compress()");
552 if (r==0) DUNE_THROW(BCRSMatrixError,"row not initialized yet");
553 if (i>=n) DUNE_THROW(BCRSMatrixError,"index out of range");
554#endif
555 return r[i];
556 }
557
560 {
561#ifdef DUNE_ISTL_WITH_CHECKING
562 if (build_mode == implicit && ready != built)
563 DUNE_THROW(BCRSMatrixError,"You cannot use operator[] in implicit build mode before calling compress()");
564 if (built!=ready) DUNE_THROW(BCRSMatrixError,"row not initialized yet");
565 if (i>=n) DUNE_THROW(BCRSMatrixError,"index out of range");
566#endif
567 return r[i];
568 }
569
570
571 //===== iterator interface to rows of the matrix
572
574 template<class T>
576 : public RandomAccessIteratorFacade<RealRowIterator<T>, T>
577 {
578
579 public:
581 typedef typename std::remove_const<T>::type ValueType;
582
585 friend class RealRowIterator<const ValueType>;
586 friend class RealRowIterator<ValueType>;
587
590 : p(_p), i(_i)
591 {}
592
595 : p(0), i(0)
596 {}
597
599 : p(it.p), i(it.i)
600 {}
601
602
605 {
606 return i;
607 }
608
609 std::ptrdiff_t distanceTo(const RealRowIterator<ValueType>& other) const
610 {
611 assert(other.p==p);
612 return (other.i-i);
613 }
614
615 std::ptrdiff_t distanceTo(const RealRowIterator<const ValueType>& other) const
616 {
617 assert(other.p==p);
618 return (other.i-i);
619 }
620
622 bool equals (const RealRowIterator<ValueType>& other) const
623 {
624 assert(other.p==p);
625 return i==other.i;
626 }
627
629 bool equals (const RealRowIterator<const ValueType>& other) const
630 {
631 assert(other.p==p);
632 return i==other.i;
633 }
634
635 private:
637 void increment()
638 {
639 ++i;
640 }
641
643 void decrement()
644 {
645 --i;
646 }
647
648 void advance(std::ptrdiff_t diff)
649 {
650 i+=diff;
651 }
652
653 T& elementAt(std::ptrdiff_t diff) const
654 {
655 return p[i+diff];
656 }
657
659 row_type& dereference () const
660 {
661 return p[i];
662 }
663
664 row_type* p;
665 size_type i;
666 };
667
671
674 {
675 return Iterator(r,0);
676 }
677
680 {
681 return Iterator(r,n);
682 }
683
687 {
688 return Iterator(r,n-1);
689 }
690
694 {
695 return Iterator(r,-1);
696 }
697
700
702 typedef typename row_type::Iterator ColIterator;
703
707
708
711 {
712 return ConstIterator(r,0);
713 }
714
717 {
718 return ConstIterator(r,n);
719 }
720
724 {
725 return ConstIterator(r,n-1);
726 }
727
731 {
732 return ConstIterator(r,-1);
733 }
734
737
739 typedef typename row_type::ConstIterator ConstColIterator;
740
741 //===== constructors & resizers
742
743 // we use a negative compressionBufferSize to indicate that the implicit
744 // mode parameters have not been set yet
745
748 : build_mode(unknown), ready(notAllocated), n(0), m(0), nnz_(0),
749 allocationSize_(0), r(0), a(0),
750 avg(0), compressionBufferSize_(-1.0)
751 {}
752
755 : build_mode(bm), ready(notAllocated), n(0), m(0), nnz_(0),
756 allocationSize_(0), r(0), a(0),
757 avg(0), compressionBufferSize_(-1.0)
758 {
759 allocate(_n, _m, _nnz,true,false);
760 }
761
764 : build_mode(bm), ready(notAllocated), n(0), m(0), nnz_(0),
765 allocationSize_(0), r(0), a(0),
766 avg(0), compressionBufferSize_(-1.0)
767 {
768 allocate(_n, _m,0,true,false);
769 }
770
772
783 BCRSMatrix (size_type _n, size_type _m, size_type _avg, double compressionBufferSize, BuildMode bm)
784 : build_mode(bm), ready(notAllocated), n(0), m(0), nnz_(0),
785 allocationSize_(0), r(0), a(0),
786 avg(_avg), compressionBufferSize_(compressionBufferSize)
787 {
788 if (bm != implicit)
789 DUNE_THROW(BCRSMatrixError,"Only call this constructor when using the implicit build mode");
790 // Prevent user from setting a negative compression buffer size:
791 // 1) It doesn't make sense
792 // 2) We use a negative value to indicate that the parameters
793 // have not been set yet
794 if (compressionBufferSize_ < 0.0)
795 DUNE_THROW(BCRSMatrixError,"You cannot set a negative overflow fraction");
796 implicit_allocate(_n,_m);
797 }
798
805 : build_mode(Mat.build_mode), ready(notAllocated), n(0), m(0), nnz_(0),
806 allocationSize_(0), r(0), a(0),
807 avg(Mat.avg), compressionBufferSize_(Mat.compressionBufferSize_)
808 {
809 if (!(Mat.ready == notAllocated || Mat.ready == built))
810 DUNE_THROW(InvalidStateException,"BCRSMatrix can only be copy-constructed when source matrix is completely empty (size not set) or fully built)");
811
812 // deep copy in global array
813 size_type _nnz = Mat.nonzeroes();
814
815 j_ = Mat.j_; // enable column index sharing, release array in case of row-wise allocation
816 allocate(Mat.n, Mat.m, _nnz, true, true);
817
818 // build window structure
820 }
821
824 {
825 deallocate();
826 }
827
833 {
834 if (ready == notAllocated)
835 {
836 build_mode = bm;
837 return;
838 }
839 if (ready == building && (build_mode == unknown || build_mode == random || build_mode == row_wise) && (bm == row_wise || bm == random))
840 build_mode = bm;
841 else
842 DUNE_THROW(InvalidStateException, "Matrix structure cannot be changed at this stage anymore (ready == "<<ready<<").");
843 }
844
860 void setSize(size_type rows, size_type columns, size_type nnz=0)
861 {
862 // deallocate already setup memory
863 deallocate();
864
865 if (build_mode == implicit)
866 {
867 if (nnz>0)
868 DUNE_THROW(Dune::BCRSMatrixError,"number of non-zeroes may not be set in implicit mode, use setImplicitBuildModeParameters() instead");
869
870 // implicit allocates differently
871 implicit_allocate(rows,columns);
872 }
873 else
874 {
875 // allocate matrix memory
876 allocate(rows, columns, nnz, true, false);
877 }
878 }
879
888 void setImplicitBuildModeParameters(size_type _avg, double compressionBufferSize)
889 {
890 // Prevent user from setting a negative compression buffer size:
891 // 1) It doesn't make sense
892 // 2) We use a negative value to indicate that the parameters
893 // have not been set yet
894 if (compressionBufferSize < 0.0)
895 DUNE_THROW(BCRSMatrixError,"You cannot set a negative compressionBufferSize value");
896
897 // make sure the parameters aren't changed after memory has been allocated
898 if (ready != notAllocated)
899 DUNE_THROW(InvalidStateException,"You cannot modify build mode parameters at this stage anymore");
900 avg = _avg;
901 compressionBufferSize_ = compressionBufferSize;
902 }
903
911 {
912 // return immediately when self-assignment
913 if (&Mat==this) return *this;
914
915 if (!((ready == notAllocated || ready == built) && (Mat.ready == notAllocated || Mat.ready == built)))
916 DUNE_THROW(InvalidStateException,"BCRSMatrix can only be copied when both target and source are empty or fully built)");
917
918 // make it simple: ALWAYS throw away memory for a and j_
919 // and deallocate rows only if n != Mat.n
920 deallocate(n!=Mat.n);
921
922 // reallocate the rows if required
923 if (n>0 && n!=Mat.n) {
924 // free rows
925 row_allocator_t rowAllocator_(allocator_);
926 for(row_type *riter=r+(n-1), *rend=r-1; riter!=rend; --riter)
927 std::allocator_traits<row_allocator_t>::destroy(rowAllocator_, riter);
928 rowAllocator_.deallocate(r,n);
929 }
930
931 nnz_ = Mat.nonzeroes();
932
933 // allocate a, share j_
934 j_ = Mat.j_;
935 allocate(Mat.n, Mat.m, nnz_, n!=Mat.n, true);
936
937 // build window structure
939 return *this;
940 }
941
944 {
945
946 if (!(ready == notAllocated || ready == built))
947 DUNE_THROW(InvalidStateException,"Scalar assignment only works on fully built BCRSMatrix)");
948
949 for (size_type i=0; i<n; i++) r[i] = k;
950 return *this;
951 }
952
953 //===== row-wise creation interface
954
957 {
958 public:
961 : Mat(_Mat), i(_i), nnz(0), current_row(Mat.j_.get())
962 {
963 if (Mat.build_mode == unknown && Mat.ready == building)
964 {
965 Mat.build_mode = row_wise;
966 }
967 if (i==0 && Mat.ready != building)
968 DUNE_THROW(BCRSMatrixError,"creation only allowed for uninitialized matrix");
969 if(Mat.build_mode!=row_wise)
970 DUNE_THROW(BCRSMatrixError,"creation only allowed if row wise allocation was requested in the constructor");
971 if(i==0 && _Mat.N()==0)
972 // empty Matrix is always built.
973 Mat.ready = built;
974 }
975
978 {
979 // this should only be called if matrix is in creation
980 if (Mat.ready != building)
981 DUNE_THROW(BCRSMatrixError,"matrix already built up");
982
983 // row i is defined through the pattern
984 // get memory for the row and initialize the j_ array
985 // this depends on the allocation mode
986
987 // compute size of the row
988 size_type s = pattern.size();
989
990 if(s>0) {
991 // update number of nonzeroes including this row
992 nnz += s;
993
994 // alloc memory / set window
995 if (Mat.nnz_ > 0)
996 {
997 // memory is allocated in one long array
998
999 // check if that memory is sufficient
1000 if (nnz > Mat.nnz_)
1001 DUNE_THROW(BCRSMatrixError,"allocated nnz too small");
1002
1003 // set row i
1004 Mat.r[i].set(s,nullptr, current_row);
1005 current_row += s;
1006 }else{
1007 // memory is allocated individually per row
1008 // allocate and set row i
1009 block_allocator_t blockAllocator_(Mat.allocator_);
1010 B* b = blockAllocator_.allocate(s);
1011 // use placement new to call constructor that allocates
1012 // additional memory.
1013 new (b) B[s];
1014 size_allocator_t sizeAllocator_(Mat.allocator_);
1015 size_type* j = sizeAllocator_.allocate(s);
1016 Mat.r[i].set(s,b,j);
1017 }
1018 }else
1019 // setup empty row
1020 Mat.r[i].set(0,nullptr,nullptr);
1021
1022 // initialize the j array for row i from pattern
1023 std::visit([&](const auto& row_indices){
1024 std::copy(row_indices.cbegin(), row_indices.cend(), Mat.r[i].getindexptr());
1025 }, pattern.storage());
1026
1027 // now go to next row
1028 i++;
1029 pattern.clear();
1030
1031 // check if this was last row
1032 if (i==Mat.n)
1033 {
1034 Mat.ready = built;
1035 if(Mat.nnz_ > 0)
1036 {
1037 // Set nnz to the exact number of nonzero blocks inserted
1038 // as some methods rely on it
1039 Mat.nnz_ = nnz;
1040 // allocate data array
1041 Mat.allocateData();
1042 Mat.setDataPointers();
1043 }
1044 }
1045 // done
1046 return *this;
1047 }
1048
1050 bool operator!= (const CreateIterator& it) const
1051 {
1052 return (i!=it.i) || (&Mat!=&it.Mat);
1053 }
1054
1056 bool operator== (const CreateIterator& it) const
1057 {
1058 return (i==it.i) && (&Mat==&it.Mat);
1059 }
1060
1063 {
1064 return i;
1065 }
1066
1069 {
1070 pattern.insert(j);
1071 }
1072
1075 {
1076 return pattern.contains(j);
1077 }
1084 {
1085 return pattern.size();
1086 }
1087
1088 private:
1089 BCRSMatrix& Mat; // the matrix we are defining
1090 size_type i; // current row to be defined
1091 size_type nnz; // count total number of nonzeros
1092 using PatternType = typename Impl::RowIndexSet;
1093 PatternType pattern; // used to compile entries in a row
1094 size_type * current_row; // row pointing to the current row to setup
1095 };
1096
1098 friend class CreateIterator;
1099
1102 {
1103 return CreateIterator(*this,0);
1104 }
1105
1108 {
1109 return CreateIterator(*this,n);
1110 }
1111
1112
1113 //===== random creation interface
1114
1122 {
1123 if (build_mode!=random)
1124 DUNE_THROW(BCRSMatrixError,"requires random build mode");
1125 if (ready != building)
1126 DUNE_THROW(BCRSMatrixError,"matrix row sizes already built up");
1127
1128 r[i].setsize(s);
1129 }
1130
1133 {
1134#ifdef DUNE_ISTL_WITH_CHECKING
1135 if (r==0) DUNE_THROW(BCRSMatrixError,"row not initialized yet");
1136 if (i>=n) DUNE_THROW(BCRSMatrixError,"index out of range");
1137#endif
1138 return r[i].getsize();
1139 }
1140
1143 {
1144 if (build_mode!=random)
1145 DUNE_THROW(BCRSMatrixError,"requires random build mode");
1146 if (ready != building)
1147 DUNE_THROW(BCRSMatrixError,"matrix row sizes already built up");
1148
1149 r[i].setsize(r[i].getsize()+s);
1150 }
1151
1154 {
1155 if (build_mode!=random)
1156 DUNE_THROW(BCRSMatrixError,"requires random build mode");
1157 if (ready != building)
1158 DUNE_THROW(BCRSMatrixError,"matrix row sizes already built up");
1159
1160 // compute total size, check positivity
1161 size_type total=0;
1162 for (size_type i=0; i<n; i++)
1163 {
1164 total += r[i].getsize();
1165 }
1166
1167 if(nnz_ == 0)
1168 // allocate/check memory
1169 allocate(n,m,total,false,false);
1170 else if(nnz_ < total)
1171 DUNE_THROW(BCRSMatrixError,"Specified number of nonzeros ("<<nnz_<<") not "
1172 <<"sufficient for calculated nonzeros ("<<total<<"! ");
1173
1174 // set the window pointers correctly
1176
1177 // initialize j_ array with m (an invalid column index)
1178 // this indicates an unused entry
1179 for (size_type k=0; k<nnz_; k++)
1180 j_.get()[k] = m;
1181 ready = rowSizesBuilt;
1182 }
1183
1185
1196 {
1197 if (build_mode!=random)
1198 DUNE_THROW(BCRSMatrixError,"requires random build mode");
1199 if (ready==built)
1200 DUNE_THROW(BCRSMatrixError,"matrix already built up");
1201 if (ready==building)
1202 DUNE_THROW(BCRSMatrixError,"matrix row sizes not built up yet");
1203 if (ready==notAllocated)
1204 DUNE_THROW(BCRSMatrixError,"matrix size not set and no memory allocated yet");
1205
1206 if (col >= m)
1207 DUNE_THROW(BCRSMatrixError,"column index exceeds matrix size");
1208
1209 // get row range
1210 size_type* const first = r[row].getindexptr();
1211 size_type* const last = first + r[row].getsize();
1212
1213 // find correct insertion position for new column index
1214 size_type* pos = std::lower_bound(first,last,col);
1215
1216 // check if index is already in row
1217 if (pos!=last && *pos == col) return;
1218
1219 // find end of already inserted column indices
1220 size_type* end = std::lower_bound(pos,last,m);
1221 if (end==last)
1222 DUNE_THROW(BCRSMatrixError,"row is too small");
1223
1224 // insert new column index at correct position
1225 std::copy_backward(pos,end,end+1);
1226 *pos = col;
1227 }
1228
1230
1238 template<typename It>
1240 {
1241 size_type row_size = r[row].size();
1242 size_type* col_begin = r[row].getindexptr();
1243 size_type* col_end;
1244 // consistency check between allocated row size and number of passed column indices
1245 if ((col_end = std::copy(begin,end,r[row].getindexptr())) != col_begin + row_size)
1246 DUNE_THROW(BCRSMatrixError,"Given size of row " << row
1247 << " (" << row_size
1248 << ") does not match number of passed entries (" << (col_end - col_begin) << ")");
1249 }
1250
1251
1253
1261 template<typename It>
1262 void setIndices(size_type row, It begin, It end)
1263 {
1264 size_type row_size = r[row].size();
1265 size_type* col_begin = r[row].getindexptr();
1266 size_type* col_end;
1267 // consistency check between allocated row size and number of passed column indices
1268 if ((col_end = std::copy(begin,end,r[row].getindexptr())) != col_begin + row_size)
1269 DUNE_THROW(BCRSMatrixError,"Given size of row " << row
1270 << " (" << row_size
1271 << ") does not match number of passed entries (" << (col_end - col_begin) << ")");
1272 std::sort(col_begin,col_end);
1273 }
1274
1277 {
1278 if (build_mode!=random)
1279 DUNE_THROW(BCRSMatrixError,"requires random build mode");
1280 if (ready==built)
1281 DUNE_THROW(BCRSMatrixError,"matrix already built up");
1282 if (ready==building)
1283 DUNE_THROW(BCRSMatrixError,"row sizes are not built up yet");
1284 if (ready==notAllocated)
1285 DUNE_THROW(BCRSMatrixError,"matrix size not set and no memory allocated yet");
1286
1287 // check if there are undefined indices
1288 RowIterator endi=end();
1289 for (RowIterator i=begin(); i!=endi; ++i)
1290 {
1291 ColIterator endj = (*i).end();
1292 for (ColIterator j=(*i).begin(); j!=endj; ++j) {
1293 if (j.index() >= m) {
1294 dwarn << "WARNING: size of row "<< i.index()<<" is "<<j.offset()<<". But was specified as being "<< (*i).end().offset()
1295 <<". This means you are wasting valuable space and creating additional cache misses!"<<std::endl;
1296 nnz_ -= ((*i).end().offset() - j.offset());
1297 r[i.index()].setsize(j.offset());
1298 break;
1299 }
1300 }
1301 }
1302
1303 allocateData();
1305
1306 // if not, set matrix to built
1307 ready = built;
1308 }
1309
1310 //===== implicit creation interface
1311
1313
1325 {
1326#ifdef DUNE_ISTL_WITH_CHECKING
1327 if (build_mode!=implicit)
1328 DUNE_THROW(BCRSMatrixError,"requires implicit build mode");
1329 if (ready==built)
1330 DUNE_THROW(BCRSMatrixError,"matrix already built up, use operator[] for entry access now");
1331 if (ready==notAllocated)
1332 DUNE_THROW(BCRSMatrixError,"matrix size not set and no memory allocated yet");
1333 if (ready!=building)
1334 DUNE_THROW(InvalidStateException,"You may only use entry() during the 'building' stage");
1335
1336 if (row >= n)
1337 DUNE_THROW(BCRSMatrixError,"row index exceeds matrix size");
1338 if (col >= m)
1339 DUNE_THROW(BCRSMatrixError,"column index exceeds matrix size");
1340#endif
1341
1342 size_type* begin = r[row].getindexptr();
1343 size_type* end = begin + r[row].getsize();
1344
1345 size_type* pos = std::find(begin, end, col);
1346
1347 //treat the case that there was a match in the array
1348 if (pos != end)
1349 if (*pos == col)
1350 {
1351 std::ptrdiff_t offset = pos - r[row].getindexptr();
1352 B* aptr = r[row].getptr() + offset;
1353
1354 return *aptr;
1355 }
1356
1357 //determine whether overflow has to be taken into account or not
1358 if (r[row].getsize() == avg)
1359 return overflow[std::make_pair(row,col)];
1360 else
1361 {
1362 //modify index array
1363 *end = col;
1364
1365 //do simultaneous operations on data array a
1366 std::ptrdiff_t offset = end - r[row].getindexptr();
1367 B* apos = r[row].getptr() + offset;
1368
1369 //increase rowsize
1370 r[row].setsize(r[row].getsize()+1);
1371
1372 //return reference to the newly created entry
1373 return *apos;
1374 }
1375 }
1376
1378
1389 {
1390 if (build_mode!=implicit)
1391 DUNE_THROW(BCRSMatrixError,"requires implicit build mode");
1392 if (ready==built)
1393 DUNE_THROW(BCRSMatrixError,"matrix already built up, no more need for compression");
1394 if (ready==notAllocated)
1395 DUNE_THROW(BCRSMatrixError,"matrix size not set and no memory allocated yet");
1396 if (ready!=building)
1397 DUNE_THROW(InvalidStateException,"You may only call compress() at the end of the 'building' stage");
1398
1399 //calculate statistics
1401 stats.overflow_total = overflow.size();
1402 stats.maximum = 0;
1403
1404 //get insertion iterators pointing to one before start (for later use of ++it)
1405 size_type* jiit = j_.get();
1406 B* aiit = a;
1407
1408 //get iterator to the smallest overflow element
1409 typename OverflowType::iterator oit = overflow.begin();
1410
1411 //store a copy of index pointers on which to perform sorting
1412 std::vector<size_type*> perm;
1413
1414 //iterate over all rows and copy elements into their position in the compressed array
1415 for (size_type i=0; i<n; i++)
1416 {
1417 //get old pointers into a and j and size without overflow changes
1418 size_type* begin = r[i].getindexptr();
1419 //B* apos = r[i].getptr();
1420 size_type size = r[i].getsize();
1421
1422 perm.resize(size);
1423
1424 typename std::vector<size_type*>::iterator it = perm.begin();
1425 for (size_type* iit = begin; iit < begin + size; ++iit, ++it)
1426 *it = iit;
1427
1428 //sort permutation array
1429 std::sort(perm.begin(),perm.end(), [](const size_type* l, const size_type* r){ return *l < *r; });
1430
1431 //change row window pointer to their new positions
1432 r[i].setindexptr(jiit);
1433 r[i].setptr(aiit);
1434
1435 for (it = perm.begin(); it != perm.end(); ++it)
1436 {
1437 //check whether there are elements in the overflow area which take precedence
1438 while ((oit!=overflow.end()) && (oit->first < std::make_pair(i,**it)))
1439 {
1440 //check whether there is enough memory to write to
1441 if (jiit > begin)
1443 "Allocated memory for BCRSMatrix exhausted during compress()!"
1444 "Please increase either the average number of entries per row or the compressionBufferSize value."
1445 );
1446 //copy an element from the overflow area to the insertion position in a and j
1447 *jiit = oit->first.second;
1448 ++jiit;
1449 *aiit = oit->second;
1450 ++aiit;
1451 ++oit;
1452 r[i].setsize(r[i].getsize()+1);
1453 }
1454
1455 //check whether there is enough memory to write to
1456 if (jiit > begin)
1458 "Allocated memory for BCRSMatrix exhausted during compress()!"
1459 "Please increase either the average number of entries per row or the compressionBufferSize value."
1460 );
1461
1462 //copy element from array
1463 *jiit = **it;
1464 ++jiit;
1465 B* apos = *it - j_.get() + a;
1466 *aiit = *apos;
1467 ++aiit;
1468 }
1469
1470 //copy remaining elements from the overflow area
1471 while ((oit!=overflow.end()) && (oit->first.first == i))
1472 {
1473 //check whether there is enough memory to write to
1474 if (jiit > begin)
1476 "Allocated memory for BCRSMatrix exhausted during compress()!"
1477 "Please increase either the average number of entries per row or the compressionBufferSize value."
1478 );
1479
1480 //copy and element from the overflow area to the insertion position in a and j
1481 *jiit = oit->first.second;
1482 ++jiit;
1483 *aiit = oit->second;
1484 ++aiit;
1485 ++oit;
1486 r[i].setsize(r[i].getsize()+1);
1487 }
1488
1489 // update maximum row size
1490 if (r[i].getsize()>stats.maximum)
1491 stats.maximum = r[i].getsize();
1492 }
1493
1494 // overflow area may be cleared
1495 overflow.clear();
1496
1497 //determine average number of entries and memory usage
1498 if ( n == 0)
1499 {
1500 stats.avg = 0;
1501 stats.mem_ratio = 1;
1502 }
1503 else
1504 {
1505 std::ptrdiff_t diff = (r[n-1].getindexptr() + r[n-1].getsize() - j_.get());
1506 nnz_ = diff;
1507 stats.avg = (double) (nnz_) / (double) n;
1508 stats.mem_ratio = (double) (nnz_) / (double) allocationSize_;
1509 }
1510
1511 //matrix is now built
1512 ready = built;
1513
1514 return stats;
1515 }
1516
1517 //===== vector space arithmetic
1518
1521 {
1522#ifdef DUNE_ISTL_WITH_CHECKING
1523 if (ready != built)
1524 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1525#endif
1526
1527 if (nnz_ > 0)
1528 {
1529 // process 1D array
1530 for (size_type i=0; i<nnz_; i++)
1531 a[i] *= k;
1532 }
1533 else
1534 {
1535 RowIterator endi=end();
1536 for (RowIterator i=begin(); i!=endi; ++i)
1537 {
1538 ColIterator endj = (*i).end();
1539 for (ColIterator j=(*i).begin(); j!=endj; ++j)
1540 (*j) *= k;
1541 }
1542 }
1543
1544 return *this;
1545 }
1546
1549 {
1550#ifdef DUNE_ISTL_WITH_CHECKING
1551 if (ready != built)
1552 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1553#endif
1554
1555 if (nnz_ > 0)
1556 {
1557 // process 1D array
1558 for (size_type i=0; i<nnz_; i++)
1559 a[i] /= k;
1560 }
1561 else
1562 {
1563 RowIterator endi=end();
1564 for (RowIterator i=begin(); i!=endi; ++i)
1565 {
1566 ColIterator endj = (*i).end();
1567 for (ColIterator j=(*i).begin(); j!=endj; ++j)
1568 (*j) /= k;
1569 }
1570 }
1571
1572 return *this;
1573 }
1574
1575
1582 {
1583#ifdef DUNE_ISTL_WITH_CHECKING
1584 if (ready != built || b.ready != built)
1585 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1586 if(N()!=b.N() || M() != b.M())
1587 DUNE_THROW(RangeError, "Matrix sizes do not match!");
1588#endif
1589 RowIterator endi=end();
1590 ConstRowIterator j=b.begin();
1591 for (RowIterator i=begin(); i!=endi; ++i, ++j) {
1592 i->operator+=(*j);
1593 }
1594
1595 return *this;
1596 }
1597
1604 {
1605#ifdef DUNE_ISTL_WITH_CHECKING
1606 if (ready != built || b.ready != built)
1607 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1608 if(N()!=b.N() || M() != b.M())
1609 DUNE_THROW(RangeError, "Matrix sizes do not match!");
1610#endif
1611 RowIterator endi=end();
1612 ConstRowIterator j=b.begin();
1613 for (RowIterator i=begin(); i!=endi; ++i, ++j) {
1614 i->operator-=(*j);
1615 }
1616
1617 return *this;
1618 }
1619
1629 {
1630#ifdef DUNE_ISTL_WITH_CHECKING
1631 if (ready != built || b.ready != built)
1632 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1633 if(N()!=b.N() || M() != b.M())
1634 DUNE_THROW(RangeError, "Matrix sizes do not match!");
1635#endif
1636 RowIterator endi=end();
1637 ConstRowIterator j=b.begin();
1638 for(RowIterator i=begin(); i!=endi; ++i, ++j)
1639 i->axpy(alpha, *j);
1640
1641 return *this;
1642 }
1643
1644 //===== linear maps
1645
1647 template<class X, class Y>
1648 void mv (const X& x, Y& y) const
1649 {
1650#ifdef DUNE_ISTL_WITH_CHECKING
1651 if (ready != built)
1652 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1653 if (x.N()!=M()) DUNE_THROW(BCRSMatrixError,
1654 "Size mismatch: M: " << N() << "x" << M() << " x: " << x.N());
1655 if (y.N()!=N()) DUNE_THROW(BCRSMatrixError,
1656 "Size mismatch: M: " << N() << "x" << M() << " y: " << y.N());
1657#endif
1658 ConstRowIterator endi=end();
1659 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1660 {
1661 y[i.index()]=0;
1662 ConstColIterator endj = (*i).end();
1663 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1664 {
1665 auto&& xj = Impl::asVector(x[j.index()]);
1666 auto&& yi = Impl::asVector(y[i.index()]);
1667 Impl::asMatrix(*j).umv(xj, yi);
1668 }
1669 }
1670 }
1671
1673 template<class X, class Y>
1674 void umv (const X& x, Y& y) const
1675 {
1676#ifdef DUNE_ISTL_WITH_CHECKING
1677 if (ready != built)
1678 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1679 if (x.N()!=M()) DUNE_THROW(BCRSMatrixError,"index out of range");
1680 if (y.N()!=N()) DUNE_THROW(BCRSMatrixError,"index out of range");
1681#endif
1682 ConstRowIterator endi=end();
1683 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1684 {
1685 ConstColIterator endj = (*i).end();
1686 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1687 {
1688 auto&& xj = Impl::asVector(x[j.index()]);
1689 auto&& yi = Impl::asVector(y[i.index()]);
1690 Impl::asMatrix(*j).umv(xj,yi);
1691 }
1692 }
1693 }
1694
1696 template<class X, class Y>
1697 void mmv (const X& x, Y& y) const
1698 {
1699#ifdef DUNE_ISTL_WITH_CHECKING
1700 if (ready != built)
1701 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1702 if (x.N()!=M()) DUNE_THROW(BCRSMatrixError,"index out of range");
1703 if (y.N()!=N()) DUNE_THROW(BCRSMatrixError,"index out of range");
1704#endif
1705 ConstRowIterator endi=end();
1706 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1707 {
1708 ConstColIterator endj = (*i).end();
1709 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1710 {
1711 auto&& xj = Impl::asVector(x[j.index()]);
1712 auto&& yi = Impl::asVector(y[i.index()]);
1713 Impl::asMatrix(*j).mmv(xj,yi);
1714 }
1715 }
1716 }
1717
1719 template<class X, class Y, class F>
1720 void usmv (F&& alpha, const X& x, Y& y) const
1721 {
1722#ifdef DUNE_ISTL_WITH_CHECKING
1723 if (ready != built)
1724 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1725 if (x.N()!=M()) DUNE_THROW(BCRSMatrixError,"index out of range");
1726 if (y.N()!=N()) DUNE_THROW(BCRSMatrixError,"index out of range");
1727#endif
1728 ConstRowIterator endi=end();
1729 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1730 {
1731 ConstColIterator endj = (*i).end();
1732 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1733 {
1734 auto&& xj = Impl::asVector(x[j.index()]);
1735 auto&& yi = Impl::asVector(y[i.index()]);
1736 Impl::asMatrix(*j).usmv(alpha,xj,yi);
1737 }
1738 }
1739 }
1740
1742 template<class X, class Y>
1743 void mtv (const X& x, Y& y) const
1744 {
1745#ifdef DUNE_ISTL_WITH_CHECKING
1746 if (ready != built)
1747 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1748 if (x.N()!=N()) DUNE_THROW(BCRSMatrixError,"index out of range");
1749 if (y.N()!=M()) DUNE_THROW(BCRSMatrixError,"index out of range");
1750#endif
1751 for(size_type i=0; i<y.N(); ++i)
1752 y[i]=0;
1753 umtv(x,y);
1754 }
1755
1757 template<class X, class Y>
1758 void umtv (const X& x, Y& y) const
1759 {
1760#ifdef DUNE_ISTL_WITH_CHECKING
1761 if (ready != built)
1762 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1763 if (x.N()!=N()) DUNE_THROW(BCRSMatrixError,"index out of range");
1764 if (y.N()!=M()) DUNE_THROW(BCRSMatrixError,"index out of range");
1765#endif
1766 ConstRowIterator endi=end();
1767 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1768 {
1769 ConstColIterator endj = (*i).end();
1770 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1771 {
1772 auto&& xi = Impl::asVector(x[i.index()]);
1773 auto&& yj = Impl::asVector(y[j.index()]);
1774 Impl::asMatrix(*j).umtv(xi,yj);
1775 }
1776 }
1777 }
1778
1780 template<class X, class Y>
1781 void mmtv (const X& x, Y& y) const
1782 {
1783#ifdef DUNE_ISTL_WITH_CHECKING
1784 if (x.N()!=N()) DUNE_THROW(BCRSMatrixError,"index out of range");
1785 if (y.N()!=M()) DUNE_THROW(BCRSMatrixError,"index out of range");
1786#endif
1787 ConstRowIterator endi=end();
1788 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1789 {
1790 ConstColIterator endj = (*i).end();
1791 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1792 {
1793 auto&& xi = Impl::asVector(x[i.index()]);
1794 auto&& yj = Impl::asVector(y[j.index()]);
1795 Impl::asMatrix(*j).mmtv(xi,yj);
1796 }
1797 }
1798 }
1799
1801 template<class X, class Y>
1802 void usmtv (const field_type& alpha, const X& x, Y& y) const
1803 {
1804#ifdef DUNE_ISTL_WITH_CHECKING
1805 if (ready != built)
1806 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1807 if (x.N()!=N()) DUNE_THROW(BCRSMatrixError,"index out of range");
1808 if (y.N()!=M()) DUNE_THROW(BCRSMatrixError,"index out of range");
1809#endif
1810 ConstRowIterator endi=end();
1811 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1812 {
1813 ConstColIterator endj = (*i).end();
1814 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1815 {
1816 auto&& xi = Impl::asVector(x[i.index()]);
1817 auto&& yj = Impl::asVector(y[j.index()]);
1818 Impl::asMatrix(*j).usmtv(alpha,xi,yj);
1819 }
1820 }
1821 }
1822
1824 template<class X, class Y>
1825 void umhv (const X& x, Y& y) const
1826 {
1827#ifdef DUNE_ISTL_WITH_CHECKING
1828 if (ready != built)
1829 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1830 if (x.N()!=N()) DUNE_THROW(BCRSMatrixError,"index out of range");
1831 if (y.N()!=M()) DUNE_THROW(BCRSMatrixError,"index out of range");
1832#endif
1833 ConstRowIterator endi=end();
1834 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1835 {
1836 ConstColIterator endj = (*i).end();
1837 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1838 {
1839 auto&& xi = Impl::asVector(x[i.index()]);
1840 auto&& yj = Impl::asVector(y[j.index()]);
1841 Impl::asMatrix(*j).umhv(xi,yj);
1842 }
1843 }
1844 }
1845
1847 template<class X, class Y>
1848 void mmhv (const X& x, Y& y) const
1849 {
1850#ifdef DUNE_ISTL_WITH_CHECKING
1851 if (ready != built)
1852 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1853 if (x.N()!=N()) DUNE_THROW(BCRSMatrixError,"index out of range");
1854 if (y.N()!=M()) DUNE_THROW(BCRSMatrixError,"index out of range");
1855#endif
1856 ConstRowIterator endi=end();
1857 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1858 {
1859 ConstColIterator endj = (*i).end();
1860 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1861 {
1862 auto&& xi = Impl::asVector(x[i.index()]);
1863 auto&& yj = Impl::asVector(y[j.index()]);
1864 Impl::asMatrix(*j).mmhv(xi,yj);
1865 }
1866 }
1867 }
1868
1870 template<class X, class Y>
1871 void usmhv (const field_type& alpha, const X& x, Y& y) const
1872 {
1873#ifdef DUNE_ISTL_WITH_CHECKING
1874 if (ready != built)
1875 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1876 if (x.N()!=N()) DUNE_THROW(BCRSMatrixError,"index out of range");
1877 if (y.N()!=M()) DUNE_THROW(BCRSMatrixError,"index out of range");
1878#endif
1879 ConstRowIterator endi=end();
1880 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1881 {
1882 ConstColIterator endj = (*i).end();
1883 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1884 {
1885 auto&& xi = Impl::asVector(x[i.index()]);
1886 auto&& yj = Impl::asVector(y[j.index()]);
1887 Impl::asMatrix(*j).usmhv(alpha,xi,yj);
1888 }
1889 }
1890 }
1891
1892
1893 //===== norms
1894
1896 typename FieldTraits<field_type>::real_type frobenius_norm2 () const
1897 {
1898#ifdef DUNE_ISTL_WITH_CHECKING
1899 if (ready != built)
1900 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1901#endif
1902
1903 typename FieldTraits<field_type>::real_type sum=0;
1904
1905 for (auto&& row : (*this))
1906 for (auto&& entry : row)
1907 sum += Impl::asMatrix(entry).frobenius_norm2();
1908
1909 return sum;
1910 }
1911
1913 typename FieldTraits<field_type>::real_type frobenius_norm () const
1914 {
1915 return sqrt(frobenius_norm2());
1916 }
1917
1919 template <typename ft = field_type,
1920 typename std::enable_if<!HasNaN<ft>::value, int>::type = 0>
1921 typename FieldTraits<ft>::real_type infinity_norm() const {
1922 if (ready != built)
1923 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1924
1925 using real_type = typename FieldTraits<ft>::real_type;
1926 using std::max;
1927
1928 real_type norm = 0;
1929 for (auto const &x : *this) {
1930 real_type sum = 0;
1931 for (auto const &y : x)
1932 sum += Impl::asMatrix(y).infinity_norm();
1933 norm = max(sum, norm);
1934 }
1935 return norm;
1936 }
1937
1939 template <typename ft = field_type,
1940 typename std::enable_if<!HasNaN<ft>::value, int>::type = 0>
1941 typename FieldTraits<ft>::real_type infinity_norm_real() const {
1942 if (ready != built)
1943 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1944
1945 using real_type = typename FieldTraits<ft>::real_type;
1946 using std::max;
1947
1948 real_type norm = 0;
1949 for (auto const &x : *this) {
1950 real_type sum = 0;
1951 for (auto const &y : x)
1952 sum += Impl::asMatrix(y).infinity_norm_real();
1953 norm = max(sum, norm);
1954 }
1955 return norm;
1956 }
1957
1959 template <typename ft = field_type,
1960 typename std::enable_if<HasNaN<ft>::value, int>::type = 0>
1961 typename FieldTraits<ft>::real_type infinity_norm() const {
1962 if (ready != built)
1963 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1964
1965 using real_type = typename FieldTraits<ft>::real_type;
1966 using std::max;
1967
1968 real_type norm = 0;
1969 real_type is_nan = 1;
1970 for (auto const &x : *this) {
1971 real_type sum = 0;
1972 for (auto const &y : x)
1973 sum += Impl::asMatrix(y).infinity_norm();
1974 norm = max(sum, norm);
1975 is_nan += sum;
1976 }
1977
1978 return norm * (is_nan / is_nan);
1979 }
1980
1982 template <typename ft = field_type,
1983 typename std::enable_if<HasNaN<ft>::value, int>::type = 0>
1984 typename FieldTraits<ft>::real_type infinity_norm_real() const {
1985 if (ready != built)
1986 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1987
1988 using real_type = typename FieldTraits<ft>::real_type;
1989 using std::max;
1990
1991 real_type norm = 0;
1992 real_type is_nan = 1;
1993
1994 for (auto const &x : *this) {
1995 real_type sum = 0;
1996 for (auto const &y : x)
1997 sum += Impl::asMatrix(y).infinity_norm_real();
1998 norm = max(sum, norm);
1999 is_nan += sum;
2000 }
2001
2002 return norm * (is_nan / is_nan);
2003 }
2004
2005 //===== sizes
2006
2008 size_type N () const
2009 {
2010 return n;
2011 }
2012
2014 size_type M () const
2015 {
2016 return m;
2017 }
2018
2021 {
2022 // in case of row-wise allocation
2023 if( nnz_ <= 0 )
2024 nnz_ = std::accumulate( begin(), end(), size_type( 0 ), [] ( size_type s, const row_type &row ) { return s+row.getsize(); } );
2025 return nnz_;
2026 }
2027
2030 {
2031 return ready;
2032 }
2033
2036 {
2037 return build_mode;
2038 }
2039
2040 //===== query
2041
2043 bool exists (size_type i, size_type j) const
2044 {
2045#ifdef DUNE_ISTL_WITH_CHECKING
2046 if (i<0 || i>=n) DUNE_THROW(BCRSMatrixError,"row index out of range");
2047 if (j<0 || j>=m) DUNE_THROW(BCRSMatrixError,"column index out of range");
2048#endif
2049 return (r[i].size() && r[i].find(j) != r[i].end());
2050 }
2051
2052
2053 protected:
2054 // state information
2055 BuildMode build_mode; // row wise or whole matrix
2056 BuildStage ready; // indicate the stage the matrix building is in
2057
2058 // size of the matrix
2059 size_type n; // number of rows
2060 size_type m; // number of columns
2061 mutable size_type nnz_; // number of nonzeroes contained in the matrix
2062 size_type allocationSize_; //allocated size of a and j arrays, except in implicit mode: nnz_==allocationsSize_
2063 // zero means that memory is allocated separately for each row.
2064
2065 // the rows are dynamically allocated
2066 row_type* r; // [n] the individual rows having pointers into a,j arrays
2067
2068 // dynamically allocated memory
2069 B* a; // [allocationSize] non-zero entries of the matrix in row-wise ordering
2070 // If a single array of column indices is used, it can be shared
2071 // between different matrices with the same sparsity pattern
2072 std::shared_ptr<size_type> j_; // [allocationSize] column indices of entries
2073
2074 // additional data is needed in implicit buildmode
2075 size_type avg;
2076 double compressionBufferSize_;
2077
2078 typedef std::map<std::pair<size_type,size_type>, B> OverflowType;
2079 OverflowType overflow;
2080
2081 // The allocator used for memory management
2082 DUNE_NO_UNIQUE_ADDRESS allocator_type allocator_;
2083
2084 // Rebound allocators
2085 // NOTE: If the original allocator A is stateful, the rebound allocators will share the same state as A, as they are rebound from A.
2086 using block_allocator_t = typename std::allocator_traits<A>::template rebind_alloc<block_type>;
2087 using row_allocator_t = typename std::allocator_traits<A>::template rebind_alloc<row_type>;
2088 using size_allocator_t = typename std::allocator_traits<A>::template rebind_alloc<size_type>;
2089
2090 void setWindowPointers(ConstRowIterator row)
2091 {
2092 row_type current_row(a,j_.get(),0); // Pointers to current row data
2093 for (size_type i=0; i<n; i++, ++row) {
2094 // set row i
2095 size_type s = row->getsize();
2096
2097 if (s>0) {
2098 // setup pointers and size
2099 r[i].set(s,current_row.getptr(), current_row.getindexptr());
2100 // update pointer for next row
2101 current_row.setptr(current_row.getptr()+s);
2102 current_row.setindexptr(current_row.getindexptr()+s);
2103 } else{
2104 // empty row
2105 r[i].set(0,nullptr,nullptr);
2106 }
2107 }
2108 }
2109
2111
2116 {
2117 size_type* jptr = j_.get();
2118 for (size_type i=0; i<n; ++i, ++row) {
2119 // set row i
2120 size_type s = row->getsize();
2121
2122 if (s>0) {
2123 // setup pointers and size
2124 r[i].setsize(s);
2125 r[i].setindexptr(jptr);
2126 } else{
2127 // empty row
2128 r[i].set(0,nullptr,nullptr);
2129 }
2130
2131 // advance position in global array
2132 jptr += s;
2133 }
2134 }
2135
2137
2142 {
2143 B* aptr = a;
2144 for (size_type i=0; i<n; ++i) {
2145 // set row i
2146 if (r[i].getsize() > 0) {
2147 // setup pointers and size
2148 r[i].setptr(aptr);
2149 } else{
2150 // empty row
2151 r[i].set(0,nullptr,nullptr);
2152 }
2153
2154 // advance position in global array
2155 aptr += r[i].getsize();
2156 }
2157 }
2158
2161 {
2162 setWindowPointers(Mat.begin());
2163
2164 // copy data
2165 for (size_type i=0; i<n; i++) r[i] = Mat.r[i];
2166
2167 // finish off
2168 build_mode = row_wise; // dummy
2169 ready = built;
2170 }
2171
2177 void deallocate(bool deallocateRows=true)
2178 {
2179
2180 if (notAllocated)
2181 return;
2182
2183 if (allocationSize_>0)
2184 {
2185 // a,j_ have been allocated as one long vector
2186 j_.reset();
2187 if (a)
2188 {
2189 block_allocator_t blockAllocator_(allocator_);
2190 for(B *aiter=a+(allocationSize_-1), *aend=a-1; aiter!=aend; --aiter)
2191 std::allocator_traits<block_allocator_t>::destroy(blockAllocator_, aiter);
2192 blockAllocator_.deallocate(a,allocationSize_);
2193 a = nullptr;
2194 }
2195 }
2196 else if (r)
2197 {
2198 // check if memory for rows have been allocated individually
2199 for (size_type i=0; i<n; i++)
2200 if (r[i].getsize()>0)
2201 {
2202 block_allocator_t blockAllocator_(allocator_);
2203 size_allocator_t sizeAllocator_(allocator_);
2204 for (B *col=r[i].getptr()+(r[i].getsize()-1),
2205 *colend = r[i].getptr()-1; col!=colend; --col) {
2206 std::allocator_traits<block_allocator_t>::destroy(blockAllocator_, col);
2207 }
2208 sizeAllocator_.deallocate(r[i].getindexptr(),1);
2209 blockAllocator_.deallocate(r[i].getptr(),1);
2210 // clear out row data in case we don't want to deallocate the rows
2211 // otherwise we might run into a double free problem here later
2212 r[i].set(0,nullptr,nullptr);
2213 }
2214 }
2215
2216 // deallocate the rows
2217 if (n>0 && deallocateRows && r) {
2218 row_allocator_t rowAllocator_(allocator_);
2219 for(row_type *riter=r+(n-1), *rend=r-1; riter!=rend; --riter)
2220 std::allocator_traits<row_allocator_t>::destroy(rowAllocator_, riter);
2221 rowAllocator_.deallocate(r,n);
2222 r = nullptr;
2223 }
2224
2225 // Mark matrix as not built at all.
2226 ready=notAllocated;
2227
2228 }
2229
2248 void allocate(size_type rows, size_type columns, size_type allocationSize, bool allocateRows, bool allocate_data)
2249 {
2250 // Store size
2251 n = rows;
2252 m = columns;
2253 nnz_ = allocationSize;
2254 allocationSize_ = allocationSize;
2255
2256 // allocate rows
2257 if(allocateRows) {
2258 if (n>0) {
2259 row_allocator_t rowAllocator_(allocator_);
2260 if (r)
2261 DUNE_THROW(InvalidStateException,"Rows have already been allocated, cannot allocate a second time");
2262 r = rowAllocator_.allocate(rows);
2263 // initialize row entries
2264 for(row_type* ri=r; ri!=r+rows; ++ri)
2265 std::allocator_traits<row_allocator_t>::construct(rowAllocator_, ri);
2266 }else{
2267 r = 0;
2268 }
2269 }
2270
2271 // allocate a and j_ array
2272 if (allocate_data)
2273 allocateData();
2274 // allocate column indices only if not yet present (enable sharing)
2275 if (allocationSize_>0) {
2276 // we copy allocator and size to the deleter since _j may outlive this class
2277 if (!j_.get()) {
2278 size_allocator_t sizeAllocator_(allocator_);
2279 j_.reset(sizeAllocator_.allocate(allocationSize_),
2280 [alloc = std::move(sizeAllocator_), size = allocationSize_](auto ptr) mutable {
2281 alloc.deallocate(ptr, size);
2282 });
2283 }
2284 }else{
2285 j_.reset();
2286 }
2287
2288 // Mark the matrix as not built.
2289 ready = building;
2290 }
2291
2292 void allocateData()
2293 {
2294 if (a)
2295 DUNE_THROW(InvalidStateException,"Cannot allocate data array (already allocated)");
2296 if (allocationSize_>0) {
2297 block_allocator_t blockAllocator_(allocator_);
2298 a = blockAllocator_.allocate(allocationSize_);
2299 // use placement new to call constructor that allocates
2300 // additional memory.
2301 new (a) B[allocationSize_];
2302 } else {
2303 a = nullptr;
2304 }
2305 }
2306
2313 {
2314 if (build_mode != implicit)
2315 DUNE_THROW(InvalidStateException,"implicit_allocate() may only be called in implicit build mode");
2316 if (ready != notAllocated)
2317 DUNE_THROW(InvalidStateException,"memory has already been allocated");
2318
2319 // check to make sure the user has actually set the parameters
2320 if (compressionBufferSize_ < 0)
2321 DUNE_THROW(InvalidStateException,"You have to set the implicit build mode parameters before starting to build the matrix");
2322 //calculate size of overflow region, add buffer for row sort!
2323 size_type osize = (size_type) (_n*avg)*compressionBufferSize_ + 4*avg;
2324 allocationSize_ = _n*avg + osize;
2325
2326 allocate(_n, _m, allocationSize_,true,true);
2327
2328 //set row pointers correctly
2329 size_type* jptr = j_.get() + osize;
2330 B* aptr = a + osize;
2331 for (size_type i=0; i<n; i++)
2332 {
2333 r[i].set(0,aptr,jptr);
2334 jptr = jptr + avg;
2335 aptr = aptr + avg;
2336 }
2337
2338 ready = building;
2339 }
2340 };
2341
2342
2343 template<class B, class A>
2344 struct FieldTraits< BCRSMatrix<B, A> >
2345 {
2346 using field_type = typename BCRSMatrix<B, A>::field_type;
2347 using real_type = typename FieldTraits<field_type>::real_type;
2348 };
2349
2352} // end namespace
2353
2354#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:957
bool operator==(const CreateIterator &it) const
equality
Definition: bcrsmatrix.hh:1056
CreateIterator & operator++()
prefix increment
Definition: bcrsmatrix.hh:977
size_type index() const
The number of the row that the iterator currently points to.
Definition: bcrsmatrix.hh:1062
bool operator!=(const CreateIterator &it) const
inequality
Definition: bcrsmatrix.hh:1050
CreateIterator(BCRSMatrix &_Mat, size_type _i)
constructor
Definition: bcrsmatrix.hh:960
void insert(size_type j)
put column index in row
Definition: bcrsmatrix.hh:1068
bool contains(size_type j)
return true if column index is in row
Definition: bcrsmatrix.hh:1074
size_type size() const
Get the current row size.
Definition: bcrsmatrix.hh:1083
Iterator access to matrix rows
Definition: bcrsmatrix.hh:577
RealRowIterator()
empty constructor, use with care!
Definition: bcrsmatrix.hh:594
bool equals(const RealRowIterator< ValueType > &other) const
equality
Definition: bcrsmatrix.hh:622
std::remove_const< T >::type ValueType
The unqualified value type.
Definition: bcrsmatrix.hh:581
bool equals(const RealRowIterator< const ValueType > &other) const
equality
Definition: bcrsmatrix.hh:629
RealRowIterator(row_type *_p, size_type _i)
constructor
Definition: bcrsmatrix.hh:589
size_type index() const
return index
Definition: bcrsmatrix.hh:604
A sparse block matrix with compressed row storage.
Definition: bcrsmatrix.hh:468
typename Imp::BlockTraits< B >::field_type field_type
export the type representing the field
Definition: bcrsmatrix.hh:490
bool exists(size_type i, size_type j) const
return true if (i,j) is in pattern
Definition: bcrsmatrix.hh:2043
BuildStage buildStage() const
The current build stage of the matrix.
Definition: bcrsmatrix.hh:2029
Iterator begin()
Get iterator to first row.
Definition: bcrsmatrix.hh:673
friend class CreateIterator
allow CreateIterator to access internal data
Definition: bcrsmatrix.hh:1098
void copyWindowStructure(const BCRSMatrix &Mat)
Copy the window structure from another matrix.
Definition: bcrsmatrix.hh:2160
B & entry(size_type row, size_type col)
Returns reference to entry (row,col) of the matrix.
Definition: bcrsmatrix.hh:1324
void mmhv(const X &x, Y &y) const
y -= A^H x
Definition: bcrsmatrix.hh:1848
void usmhv(const field_type &alpha, const X &x, Y &y) const
y += alpha A^H x
Definition: bcrsmatrix.hh:1871
void umtv(const X &x, Y &y) const
y += A^T x
Definition: bcrsmatrix.hh:1758
RealRowIterator< const row_type > const_iterator
The const iterator over the matrix rows.
Definition: bcrsmatrix.hh:705
std::allocator_traits< A >::size_type size_type
The type for the index access and the size.
Definition: bcrsmatrix.hh:496
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:2248
BCRSMatrix & axpy(field_type alpha, const BCRSMatrix &b)
Add the scaled entries of another matrix to this one.
Definition: bcrsmatrix.hh:1628
FieldTraits< ft >::real_type infinity_norm_real() const
simplified infinity norm (uses Manhattan norm for complex values)
Definition: bcrsmatrix.hh:1941
~BCRSMatrix()
destructor
Definition: bcrsmatrix.hh:823
void deallocate(bool deallocateRows=true)
deallocate memory of the matrix.
Definition: bcrsmatrix.hh:2177
Iterator RowIterator
rename the iterators for easier access
Definition: bcrsmatrix.hh:699
row_type & operator[](size_type i)
random access to the rows
Definition: bcrsmatrix.hh:547
BCRSMatrix()
an empty matrix
Definition: bcrsmatrix.hh:747
void endrowsizes()
indicate that size of all rows is defined
Definition: bcrsmatrix.hh:1153
void incrementrowsize(size_type i, size_type s=1)
increment size of row i by s (1 by default)
Definition: bcrsmatrix.hh:1142
void mtv(const X &x, Y &y) const
y = A^T x
Definition: bcrsmatrix.hh:1743
void umhv(const X &x, Y &y) const
y += A^H x
Definition: bcrsmatrix.hh:1825
size_type nonzeroes() const
number of blocks that are stored (the number of blocks that possibly are nonzero)
Definition: bcrsmatrix.hh:2020
ConstIterator ConstRowIterator
rename the const row iterator for easier access
Definition: bcrsmatrix.hh:736
void setIndicesNoSort(size_type row, It begin, It end)
Set all column indices for row from the given iterator range.
Definition: bcrsmatrix.hh:1239
void setrowsize(size_type i, size_type s)
Set number of indices in row i to s.
Definition: bcrsmatrix.hh:1121
BCRSMatrix & operator*=(const field_type &k)
vector space multiplication with scalar
Definition: bcrsmatrix.hh:1520
RealRowIterator< row_type > iterator
The iterator over the (mutable matrix rows.
Definition: bcrsmatrix.hh:669
void usmtv(const field_type &alpha, const X &x, Y &y) const
y += alpha A^T x
Definition: bcrsmatrix.hh:1802
ConstIterator beforeBegin() const
Definition: bcrsmatrix.hh:730
Iterator beforeBegin()
Definition: bcrsmatrix.hh:693
BuildMode
we support two modes
Definition: bcrsmatrix.hh:508
@ implicit
Build entries randomly with an educated guess for the number of entries per row.
Definition: bcrsmatrix.hh:537
@ unknown
Build mode not set!
Definition: bcrsmatrix.hh:541
@ random
Build entries randomly.
Definition: bcrsmatrix.hh:528
@ row_wise
Build in a row-wise manner.
Definition: bcrsmatrix.hh:519
A allocator_type
allocator of blocks
Definition: bcrsmatrix.hh:502
BCRSMatrix(size_type _n, size_type _m, size_type _nnz, BuildMode bm)
matrix with known number of nonzeroes
Definition: bcrsmatrix.hh:754
Imp::CompressedBlockVectorWindow< B, size_type > row_type
implement row_type with compressed vector
Definition: bcrsmatrix.hh:499
::Dune::CompressionStatistics< size_type > CompressionStatistics
The type for the statistics object returned by compress()
Definition: bcrsmatrix.hh:505
BCRSMatrix & operator-=(const BCRSMatrix &b)
Subtract the entries of another matrix from this one.
Definition: bcrsmatrix.hh:1603
BCRSMatrix(const BCRSMatrix &Mat)
copy constructor
Definition: bcrsmatrix.hh:804
Iterator end()
Get iterator to one beyond last row.
Definition: bcrsmatrix.hh:679
void setIndices(size_type row, It begin, It end)
Set all column indices for row from the given iterator range.
Definition: bcrsmatrix.hh:1262
void addindex(size_type row, size_type col)
add index (row,col) to the matrix
Definition: bcrsmatrix.hh:1195
row_type::Iterator ColIterator
Iterator for the entries of each row.
Definition: bcrsmatrix.hh:702
FieldTraits< field_type >::real_type frobenius_norm() const
frobenius norm: sqrt(sum over squared values of entries)
Definition: bcrsmatrix.hh:1913
BCRSMatrix & operator/=(const field_type &k)
vector space division by scalar
Definition: bcrsmatrix.hh:1548
BCRSMatrix & operator+=(const BCRSMatrix &b)
Add the entries of another matrix to this one.
Definition: bcrsmatrix.hh:1581
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:783
CreateIterator createend()
get create iterator pointing to one after the last block
Definition: bcrsmatrix.hh:1107
FieldTraits< field_type >::real_type frobenius_norm2() const
square of frobenius norm, need for block recursion
Definition: bcrsmatrix.hh:1896
Iterator beforeEnd()
Definition: bcrsmatrix.hh:686
row_type::ConstIterator ConstColIterator
Const iterator to the entries of a row.
Definition: bcrsmatrix.hh:739
void usmv(F &&alpha, const X &x, Y &y) const
y += alpha A x
Definition: bcrsmatrix.hh:1720
size_type getrowsize(size_type i) const
get current number of indices in row i
Definition: bcrsmatrix.hh:1132
size_type M() const
number of columns (counted in blocks)
Definition: bcrsmatrix.hh:2014
CreateIterator createbegin()
get initial create iterator
Definition: bcrsmatrix.hh:1101
BuildStage
Definition: bcrsmatrix.hh:471
@ rowSizesBuilt
The row sizes of the matrix are known.
Definition: bcrsmatrix.hh:482
@ built
The matrix structure is fully built.
Definition: bcrsmatrix.hh:484
@ notbuilt
Matrix is not built at all, no memory has been allocated, build mode and size can still be set.
Definition: bcrsmatrix.hh:473
@ notAllocated
Matrix is not built at all, no memory has been allocated, build mode and size can still be set.
Definition: bcrsmatrix.hh:475
@ building
Matrix is currently being built, some memory has been allocated, build mode and size are fixed.
Definition: bcrsmatrix.hh:477
BuildMode buildMode() const
The currently selected build mode of the matrix.
Definition: bcrsmatrix.hh:2035
void mmv(const X &x, Y &y) const
y -= A x
Definition: bcrsmatrix.hh:1697
FieldTraits< ft >::real_type infinity_norm() const
infinity norm (row sum norm, how to generalize for blocks?)
Definition: bcrsmatrix.hh:1921
void mv(const X &x, Y &y) const
y = A x
Definition: bcrsmatrix.hh:1648
B block_type
export the type representing the components
Definition: bcrsmatrix.hh:493
void mmtv(const X &x, Y &y) const
y -= A^T x
Definition: bcrsmatrix.hh:1781
void umv(const X &x, Y &y) const
y += A x
Definition: bcrsmatrix.hh:1674
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:2312
void setImplicitBuildModeParameters(size_type _avg, double compressionBufferSize)
Set parameters needed for creation in implicit build mode.
Definition: bcrsmatrix.hh:888
BCRSMatrix(size_type _n, size_type _m, BuildMode bm)
matrix with unknown number of nonzeroes
Definition: bcrsmatrix.hh:763
void endindices()
indicate that all indices are defined, check consistency
Definition: bcrsmatrix.hh:1276
CompressionStatistics compress()
Finishes the buildstage in implicit mode.
Definition: bcrsmatrix.hh:1388
void setDataPointers()
Set data pointers for all rows.
Definition: bcrsmatrix.hh:2141
size_type N() const
number of rows (counted in blocks)
Definition: bcrsmatrix.hh:2008
void setBuildMode(BuildMode bm)
Sets the build mode of the matrix.
Definition: bcrsmatrix.hh:832
void setSize(size_type rows, size_type columns, size_type nnz=0)
Set the size of the matrix.
Definition: bcrsmatrix.hh:860
BCRSMatrix & operator=(const BCRSMatrix &Mat)
assignment
Definition: bcrsmatrix.hh:910
void setColumnPointers(ConstRowIterator row)
Copy row sizes from iterator range starting at row and set column index pointers for all rows.
Definition: bcrsmatrix.hh:2115
ConstIterator end() const
Get const iterator to one beyond last row.
Definition: bcrsmatrix.hh:716
ConstIterator begin() const
Get const iterator to first row.
Definition: bcrsmatrix.hh:710
ConstIterator beforeEnd() const
Definition: bcrsmatrix.hh:723
Proxy row object for entry access.
Definition: bcrsmatrix.hh:139
block_type & operator[](size_type j) const
Returns entry in column j.
Definition: bcrsmatrix.hh:144
A wrapper for uniform access to the BCRSMatrix during and after the build stage in implicit build mod...
Definition: bcrsmatrix.hh:119
Matrix::block_type block_type
The block_type of the underlying matrix.
Definition: bcrsmatrix.hh:127
ImplicitMatrixBuilder(Matrix &m)
Creates an ImplicitMatrixBuilder for matrix m.
Definition: bcrsmatrix.hh:172
M_ Matrix
The underlying matrix.
Definition: bcrsmatrix.hh:124
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:196
size_type M() const
The number of columns in the matrix.
Definition: bcrsmatrix.hh:219
Matrix::size_type size_type
The size_type of the underlying matrix.
Definition: bcrsmatrix.hh:130
row_object operator[](size_type i) const
Returns a proxy for entries in row i.
Definition: bcrsmatrix.hh:207
size_type N() const
The number of rows in the matrix.
Definition: bcrsmatrix.hh:213
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:489
constexpr decltype(auto) elementAt(Container &&c, Index &&i)
Get element at given position from container.
Definition: hybridutilities.hh:131
constexpr T accumulate(Range &&range, T value, F &&f)
Accumulate values.
Definition: hybridutilities.hh:284
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.
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:90
size_type overflow_total
total number of elements written to the overflow area during construction.
Definition: bcrsmatrix.hh:96
size_type maximum
maximum number of non-zeroes per row.
Definition: bcrsmatrix.hh:94
double avg
average number of non-zeroes per row.
Definition: bcrsmatrix.hh:92
double mem_ratio
fraction of wasted memory resulting from non-used overflow area.
Definition: bcrsmatrix.hh:101
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 (Mar 8, 23:37, 2026)