Dune Core Modules (2.9.0)

bcrsmatrix.hh
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (C) 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"
28
30
35namespace Dune {
36
76 template<typename M>
77 struct MatrixDimension;
78
80
86 template<typename size_type>
88 {
90 double avg;
92 size_type maximum;
94 size_type overflow_total;
96
99 double mem_ratio;
100 };
101
103
115 template<class M_>
117 {
118
119 public:
120
122 typedef M_ Matrix;
123
126
128 typedef typename Matrix::size_type size_type;
129
131
137 {
138
139 public:
140
143 {
144 return _m.entry(_i,j);
145 }
146
147#ifndef DOXYGEN
148
150 : _m(m)
151 , _i(i)
152 {}
153
154#endif
155
156 private:
157
158 Matrix& _m;
159 size_type _i;
160
161 };
162
164
171 : _m(m)
172 {
173 if (m.buildMode() != Matrix::implicit)
174 DUNE_THROW(BCRSMatrixError,"You can only create an ImplicitBuilder for a matrix in implicit build mode");
175 if (m.buildStage() != Matrix::building)
176 DUNE_THROW(BCRSMatrixError,"You can only create an ImplicitBuilder for a matrix with set size that has not been compressed() yet");
177 }
178
180
194 ImplicitMatrixBuilder(Matrix& m, size_type rows, size_type cols, size_type avg_cols_per_row, double overflow_fraction)
195 : _m(m)
196 {
197 if (m.buildStage() != Matrix::notAllocated)
198 DUNE_THROW(BCRSMatrixError,"You can only set up a matrix for this ImplicitBuilder if it has no memory allocated yet");
199 m.setBuildMode(Matrix::implicit);
200 m.setImplicitBuildModeParameters(avg_cols_per_row,overflow_fraction);
201 m.setSize(rows,cols);
202 }
203
206 {
207 return row_object(_m,i);
208 }
209
211 size_type N() const
212 {
213 return _m.N();
214 }
215
217 size_type M() const
218 {
219 return _m.M();
220 }
221
222 private:
223
224 Matrix& _m;
225
226 };
227
464 template<class B, class A=std::allocator<B> >
466 {
467 friend struct MatrixDimension<BCRSMatrix>;
468 public:
482 built=3
483 };
484
485 //===== type definitions and constants
486
488 using field_type = typename Imp::BlockTraits<B>::field_type;
489
491 typedef B block_type;
492
494 typedef A allocator_type;
495
497 typedef Imp::CompressedBlockVectorWindow<B,A> row_type;
498
500 typedef typename A::size_type size_type;
501
503 typedef ::Dune::CompressionStatistics<size_type> CompressionStatistics;
504
506 [[deprecated("Use free function blockLevel(). Will be removed after 2.8.")]]
507 static constexpr unsigned int blocklevel = blockLevel<B>()+1;
508
543 unknown
544 };
545
546 //===== random access interface to rows of the matrix
547
550 {
551#ifdef DUNE_ISTL_WITH_CHECKING
552 if (build_mode == implicit && ready != built)
553 DUNE_THROW(BCRSMatrixError,"You cannot use operator[] in implicit build mode before calling compress()");
554 if (r==0) DUNE_THROW(BCRSMatrixError,"row not initialized yet");
555 if (i>=n) DUNE_THROW(BCRSMatrixError,"index out of range");
556#endif
557 return r[i];
558 }
559
562 {
563#ifdef DUNE_ISTL_WITH_CHECKING
564 if (build_mode == implicit && ready != built)
565 DUNE_THROW(BCRSMatrixError,"You cannot use operator[] in implicit build mode before calling compress()");
566 if (built!=ready) DUNE_THROW(BCRSMatrixError,"row not initialized yet");
567 if (i>=n) DUNE_THROW(BCRSMatrixError,"index out of range");
568#endif
569 return r[i];
570 }
571
572
573 //===== iterator interface to rows of the matrix
574
576 template<class T>
578 : public RandomAccessIteratorFacade<RealRowIterator<T>, T>
579 {
580
581 public:
583 typedef typename std::remove_const<T>::type ValueType;
584
587 friend class RealRowIterator<const ValueType>;
588 friend class RealRowIterator<ValueType>;
589
592 : p(_p), i(_i)
593 {}
594
597 : p(0), i(0)
598 {}
599
601 : p(it.p), i(it.i)
602 {}
603
604
607 {
608 return i;
609 }
610
611 std::ptrdiff_t distanceTo(const RealRowIterator<ValueType>& other) const
612 {
613 assert(other.p==p);
614 return (other.i-i);
615 }
616
617 std::ptrdiff_t distanceTo(const RealRowIterator<const ValueType>& other) const
618 {
619 assert(other.p==p);
620 return (other.i-i);
621 }
622
624 bool equals (const RealRowIterator<ValueType>& other) const
625 {
626 assert(other.p==p);
627 return i==other.i;
628 }
629
631 bool equals (const RealRowIterator<const ValueType>& other) const
632 {
633 assert(other.p==p);
634 return i==other.i;
635 }
636
637 private:
639 void increment()
640 {
641 ++i;
642 }
643
645 void decrement()
646 {
647 --i;
648 }
649
650 void advance(std::ptrdiff_t diff)
651 {
652 i+=diff;
653 }
654
655 T& elementAt(std::ptrdiff_t diff) const
656 {
657 return p[i+diff];
658 }
659
661 row_type& dereference () const
662 {
663 return p[i];
664 }
665
666 row_type* p;
667 size_type i;
668 };
669
673
676 {
677 return Iterator(r,0);
678 }
679
682 {
683 return Iterator(r,n);
684 }
685
689 {
690 return Iterator(r,n-1);
691 }
692
696 {
697 return Iterator(r,-1);
698 }
699
702
704 typedef typename row_type::Iterator ColIterator;
705
709
710
713 {
714 return ConstIterator(r,0);
715 }
716
719 {
720 return ConstIterator(r,n);
721 }
722
726 {
727 return ConstIterator(r,n-1);
728 }
729
733 {
734 return ConstIterator(r,-1);
735 }
736
739
741 typedef typename row_type::ConstIterator ConstColIterator;
742
743 //===== constructors & resizers
744
745 // we use a negative compressionBufferSize to indicate that the implicit
746 // mode parameters have not been set yet
747
750 : build_mode(unknown), ready(notAllocated), n(0), m(0), nnz_(0),
751 allocationSize_(0), r(0), a(0),
752 avg(0), compressionBufferSize_(-1.0)
753 {}
754
757 : build_mode(bm), ready(notAllocated), n(0), m(0), nnz_(0),
758 allocationSize_(0), r(0), a(0),
759 avg(0), compressionBufferSize_(-1.0)
760 {
761 allocate(_n, _m, _nnz,true,false);
762 }
763
766 : build_mode(bm), ready(notAllocated), n(0), m(0), nnz_(0),
767 allocationSize_(0), r(0), a(0),
768 avg(0), compressionBufferSize_(-1.0)
769 {
770 allocate(_n, _m,0,true,false);
771 }
772
774
784 BCRSMatrix (size_type _n, size_type _m, size_type _avg, double compressionBufferSize, BuildMode bm)
785 : build_mode(bm), ready(notAllocated), n(0), m(0), nnz_(0),
786 allocationSize_(0), r(0), a(0),
787 avg(_avg), compressionBufferSize_(compressionBufferSize)
788 {
789 if (bm != implicit)
790 DUNE_THROW(BCRSMatrixError,"Only call this constructor when using the implicit build mode");
791 // Prevent user from setting a negative compression buffer size:
792 // 1) It doesn't make sense
793 // 2) We use a negative value to indicate that the parameters
794 // have not been set yet
795 if (compressionBufferSize_ < 0.0)
796 DUNE_THROW(BCRSMatrixError,"You cannot set a negative overflow fraction");
797 implicit_allocate(_n,_m);
798 }
799
806 : build_mode(Mat.build_mode), ready(notAllocated), n(0), m(0), nnz_(0),
807 allocationSize_(0), r(0), a(0),
808 avg(Mat.avg), compressionBufferSize_(Mat.compressionBufferSize_)
809 {
810 if (!(Mat.ready == notAllocated || Mat.ready == built))
811 DUNE_THROW(InvalidStateException,"BCRSMatrix can only be copy-constructed when source matrix is completely empty (size not set) or fully built)");
812
813 // deep copy in global array
814 size_type _nnz = Mat.nonzeroes();
815
816 j_ = Mat.j_; // enable column index sharing, release array in case of row-wise allocation
817 allocate(Mat.n, Mat.m, _nnz, true, true);
818
819 // build window structure
821 }
822
825 {
826 deallocate();
827 }
828
834 {
835 if (ready == notAllocated)
836 {
837 build_mode = bm;
838 return;
839 }
840 if (ready == building && (build_mode == unknown || build_mode == random || build_mode == row_wise) && (bm == row_wise || bm == random))
841 build_mode = bm;
842 else
843 DUNE_THROW(InvalidStateException, "Matrix structure cannot be changed at this stage anymore (ready == "<<ready<<").");
844 }
845
861 void setSize(size_type rows, size_type columns, size_type nnz=0)
862 {
863 // deallocate already setup memory
864 deallocate();
865
866 if (build_mode == implicit)
867 {
868 if (nnz>0)
869 DUNE_THROW(Dune::BCRSMatrixError,"number of non-zeroes may not be set in implicit mode, use setImplicitBuildModeParameters() instead");
870
871 // implicit allocates differently
872 implicit_allocate(rows,columns);
873 }
874 else
875 {
876 // allocate matrix memory
877 allocate(rows, columns, nnz, true, false);
878 }
879 }
880
889 void setImplicitBuildModeParameters(size_type _avg, double compressionBufferSize)
890 {
891 // Prevent user from setting a negative compression buffer size:
892 // 1) It doesn't make sense
893 // 2) We use a negative value to indicate that the parameters
894 // have not been set yet
895 if (compressionBufferSize < 0.0)
896 DUNE_THROW(BCRSMatrixError,"You cannot set a negative compressionBufferSize value");
897
898 // make sure the parameters aren't changed after memory has been allocated
899 if (ready != notAllocated)
900 DUNE_THROW(InvalidStateException,"You cannot modify build mode parameters at this stage anymore");
901 avg = _avg;
902 compressionBufferSize_ = compressionBufferSize;
903 }
904
912 {
913 // return immediately when self-assignment
914 if (&Mat==this) return *this;
915
916 if (!((ready == notAllocated || ready == built) && (Mat.ready == notAllocated || Mat.ready == built)))
917 DUNE_THROW(InvalidStateException,"BCRSMatrix can only be copied when both target and source are empty or fully built)");
918
919 // make it simple: ALWAYS throw away memory for a and j_
920 // and deallocate rows only if n != Mat.n
921 deallocate(n!=Mat.n);
922
923 // reallocate the rows if required
924 if (n>0 && n!=Mat.n) {
925 // free rows
926 for(row_type *riter=r+(n-1), *rend=r-1; riter!=rend; --riter)
927 std::allocator_traits<decltype(rowAllocator_)>::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(nullptr, Mat.j_.get(), 0)
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.getindexptr());
1005 current_row.setindexptr(current_row.getindexptr()+s);
1006 }else{
1007 // memory is allocated individually per row
1008 // allocate and set row i
1009 B* b = Mat.allocator_.allocate(s);
1010 // use placement new to call constructor that allocates
1011 // additional memory.
1012 new (b) B[s];
1013 size_type* j = Mat.sizeAllocator_.allocate(s);
1014 Mat.r[i].set(s,b,j);
1015 }
1016 }else
1017 // setup empty row
1018 Mat.r[i].set(0,nullptr,nullptr);
1019
1020 // initialize the j array for row i from pattern
1021 std::copy(pattern.cbegin(), pattern.cend(), Mat.r[i].getindexptr());
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.find(j) != pattern.end();
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 typedef std::set<size_type,std::less<size_type> > PatternType;
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
1233 template<typename It>
1234 void setIndices(size_type row, It begin, It end)
1235 {
1236 size_type row_size = r[row].size();
1237 size_type* col_begin = r[row].getindexptr();
1238 size_type* col_end;
1239 // consistency check between allocated row size and number of passed column indices
1240 if ((col_end = std::copy(begin,end,r[row].getindexptr())) != col_begin + row_size)
1241 DUNE_THROW(BCRSMatrixError,"Given size of row " << row
1242 << " (" << row_size
1243 << ") does not match number of passed entries (" << (col_end - col_begin) << ")");
1244 std::sort(col_begin,col_end);
1245 }
1246
1249 {
1250 if (build_mode!=random)
1251 DUNE_THROW(BCRSMatrixError,"requires random build mode");
1252 if (ready==built)
1253 DUNE_THROW(BCRSMatrixError,"matrix already built up");
1254 if (ready==building)
1255 DUNE_THROW(BCRSMatrixError,"row sizes are not built up yet");
1256 if (ready==notAllocated)
1257 DUNE_THROW(BCRSMatrixError,"matrix size not set and no memory allocated yet");
1258
1259 // check if there are undefined indices
1260 RowIterator endi=end();
1261 for (RowIterator i=begin(); i!=endi; ++i)
1262 {
1263 ColIterator endj = (*i).end();
1264 for (ColIterator j=(*i).begin(); j!=endj; ++j) {
1265 if (j.index() >= m) {
1266 dwarn << "WARNING: size of row "<< i.index()<<" is "<<j.offset()<<". But was specified as being "<< (*i).end().offset()
1267 <<". This means you are wasting valuable space and creating additional cache misses!"<<std::endl;
1268 nnz_ -= ((*i).end().offset() - j.offset());
1269 r[i.index()].setsize(j.offset());
1270 break;
1271 }
1272 }
1273 }
1274
1275 allocateData();
1277
1278 // if not, set matrix to built
1279 ready = built;
1280 }
1281
1282 //===== implicit creation interface
1283
1285
1297 {
1298#ifdef DUNE_ISTL_WITH_CHECKING
1299 if (build_mode!=implicit)
1300 DUNE_THROW(BCRSMatrixError,"requires implicit build mode");
1301 if (ready==built)
1302 DUNE_THROW(BCRSMatrixError,"matrix already built up, use operator[] for entry access now");
1303 if (ready==notAllocated)
1304 DUNE_THROW(BCRSMatrixError,"matrix size not set and no memory allocated yet");
1305 if (ready!=building)
1306 DUNE_THROW(InvalidStateException,"You may only use entry() during the 'building' stage");
1307
1308 if (row >= n)
1309 DUNE_THROW(BCRSMatrixError,"row index exceeds matrix size");
1310 if (col >= m)
1311 DUNE_THROW(BCRSMatrixError,"column index exceeds matrix size");
1312#endif
1313
1314 size_type* begin = r[row].getindexptr();
1315 size_type* end = begin + r[row].getsize();
1316
1317 size_type* pos = std::find(begin, end, col);
1318
1319 //treat the case that there was a match in the array
1320 if (pos != end)
1321 if (*pos == col)
1322 {
1323 std::ptrdiff_t offset = pos - r[row].getindexptr();
1324 B* aptr = r[row].getptr() + offset;
1325
1326 return *aptr;
1327 }
1328
1329 //determine whether overflow has to be taken into account or not
1330 if (r[row].getsize() == avg)
1331 return overflow[std::make_pair(row,col)];
1332 else
1333 {
1334 //modify index array
1335 *end = col;
1336
1337 //do simultaneous operations on data array a
1338 std::ptrdiff_t offset = end - r[row].getindexptr();
1339 B* apos = r[row].getptr() + offset;
1340
1341 //increase rowsize
1342 r[row].setsize(r[row].getsize()+1);
1343
1344 //return reference to the newly created entry
1345 return *apos;
1346 }
1347 }
1348
1350
1361 {
1362 if (build_mode!=implicit)
1363 DUNE_THROW(BCRSMatrixError,"requires implicit build mode");
1364 if (ready==built)
1365 DUNE_THROW(BCRSMatrixError,"matrix already built up, no more need for compression");
1366 if (ready==notAllocated)
1367 DUNE_THROW(BCRSMatrixError,"matrix size not set and no memory allocated yet");
1368 if (ready!=building)
1369 DUNE_THROW(InvalidStateException,"You may only call compress() at the end of the 'building' stage");
1370
1371 //calculate statistics
1373 stats.overflow_total = overflow.size();
1374 stats.maximum = 0;
1375
1376 //get insertion iterators pointing to one before start (for later use of ++it)
1377 size_type* jiit = j_.get();
1378 B* aiit = a;
1379
1380 //get iterator to the smallest overflow element
1381 typename OverflowType::iterator oit = overflow.begin();
1382
1383 //store a copy of index pointers on which to perform sorting
1384 std::vector<size_type*> perm;
1385
1386 //iterate over all rows and copy elements into their position in the compressed array
1387 for (size_type i=0; i<n; i++)
1388 {
1389 //get old pointers into a and j and size without overflow changes
1390 size_type* begin = r[i].getindexptr();
1391 //B* apos = r[i].getptr();
1392 size_type size = r[i].getsize();
1393
1394 perm.resize(size);
1395
1396 typename std::vector<size_type*>::iterator it = perm.begin();
1397 for (size_type* iit = begin; iit < begin + size; ++iit, ++it)
1398 *it = iit;
1399
1400 //sort permutation array
1401 std::sort(perm.begin(),perm.end(),PointerCompare<size_type>());
1402
1403 //change row window pointer to their new positions
1404 r[i].setindexptr(jiit);
1405 r[i].setptr(aiit);
1406
1407 for (it = perm.begin(); it != perm.end(); ++it)
1408 {
1409 //check whether there are elements in the overflow area which take precedence
1410 while ((oit!=overflow.end()) && (oit->first < std::make_pair(i,**it)))
1411 {
1412 //check whether there is enough memory to write to
1413 if (jiit > begin)
1415 "Allocated memory for BCRSMatrix exhausted during compress()!"
1416 "Please increase either the average number of entries per row or the compressionBufferSize value."
1417 );
1418 //copy an element from the overflow area to the insertion position in a and j
1419 *jiit = oit->first.second;
1420 ++jiit;
1421 *aiit = oit->second;
1422 ++aiit;
1423 ++oit;
1424 r[i].setsize(r[i].getsize()+1);
1425 }
1426
1427 //check whether there is enough memory to write to
1428 if (jiit > begin)
1430 "Allocated memory for BCRSMatrix exhausted during compress()!"
1431 "Please increase either the average number of entries per row or the compressionBufferSize value."
1432 );
1433
1434 //copy element from array
1435 *jiit = **it;
1436 ++jiit;
1437 B* apos = *it - j_.get() + a;
1438 *aiit = *apos;
1439 ++aiit;
1440 }
1441
1442 //copy remaining elements from the overflow area
1443 while ((oit!=overflow.end()) && (oit->first.first == i))
1444 {
1445 //check whether there is enough memory to write to
1446 if (jiit > begin)
1448 "Allocated memory for BCRSMatrix exhausted during compress()!"
1449 "Please increase either the average number of entries per row or the compressionBufferSize value."
1450 );
1451
1452 //copy and element from the overflow area to the insertion position in a and j
1453 *jiit = oit->first.second;
1454 ++jiit;
1455 *aiit = oit->second;
1456 ++aiit;
1457 ++oit;
1458 r[i].setsize(r[i].getsize()+1);
1459 }
1460
1461 // update maximum row size
1462 if (r[i].getsize()>stats.maximum)
1463 stats.maximum = r[i].getsize();
1464 }
1465
1466 // overflow area may be cleared
1467 overflow.clear();
1468
1469 //determine average number of entries and memory usage
1470 std::ptrdiff_t diff = (r[n-1].getindexptr() + r[n-1].getsize() - j_.get());
1471 nnz_ = diff;
1472 stats.avg = (double) (nnz_) / (double) n;
1473 stats.mem_ratio = (double) (nnz_) / (double) allocationSize_;
1474
1475 //matrix is now built
1476 ready = built;
1477
1478 return stats;
1479 }
1480
1481 //===== vector space arithmetic
1482
1485 {
1486#ifdef DUNE_ISTL_WITH_CHECKING
1487 if (ready != built)
1488 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1489#endif
1490
1491 if (nnz_ > 0)
1492 {
1493 // process 1D array
1494 for (size_type i=0; i<nnz_; i++)
1495 a[i] *= k;
1496 }
1497 else
1498 {
1499 RowIterator endi=end();
1500 for (RowIterator i=begin(); i!=endi; ++i)
1501 {
1502 ColIterator endj = (*i).end();
1503 for (ColIterator j=(*i).begin(); j!=endj; ++j)
1504 (*j) *= k;
1505 }
1506 }
1507
1508 return *this;
1509 }
1510
1513 {
1514#ifdef DUNE_ISTL_WITH_CHECKING
1515 if (ready != built)
1516 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1517#endif
1518
1519 if (nnz_ > 0)
1520 {
1521 // process 1D array
1522 for (size_type i=0; i<nnz_; i++)
1523 a[i] /= k;
1524 }
1525 else
1526 {
1527 RowIterator endi=end();
1528 for (RowIterator i=begin(); i!=endi; ++i)
1529 {
1530 ColIterator endj = (*i).end();
1531 for (ColIterator j=(*i).begin(); j!=endj; ++j)
1532 (*j) /= k;
1533 }
1534 }
1535
1536 return *this;
1537 }
1538
1539
1546 {
1547#ifdef DUNE_ISTL_WITH_CHECKING
1548 if (ready != built || b.ready != built)
1549 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1550 if(N()!=b.N() || M() != b.M())
1551 DUNE_THROW(RangeError, "Matrix sizes do not match!");
1552#endif
1553 RowIterator endi=end();
1554 ConstRowIterator j=b.begin();
1555 for (RowIterator i=begin(); i!=endi; ++i, ++j) {
1556 i->operator+=(*j);
1557 }
1558
1559 return *this;
1560 }
1561
1568 {
1569#ifdef DUNE_ISTL_WITH_CHECKING
1570 if (ready != built || b.ready != built)
1571 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1572 if(N()!=b.N() || M() != b.M())
1573 DUNE_THROW(RangeError, "Matrix sizes do not match!");
1574#endif
1575 RowIterator endi=end();
1576 ConstRowIterator j=b.begin();
1577 for (RowIterator i=begin(); i!=endi; ++i, ++j) {
1578 i->operator-=(*j);
1579 }
1580
1581 return *this;
1582 }
1583
1593 {
1594#ifdef DUNE_ISTL_WITH_CHECKING
1595 if (ready != built || b.ready != built)
1596 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1597 if(N()!=b.N() || M() != b.M())
1598 DUNE_THROW(RangeError, "Matrix sizes do not match!");
1599#endif
1600 RowIterator endi=end();
1601 ConstRowIterator j=b.begin();
1602 for(RowIterator i=begin(); i!=endi; ++i, ++j)
1603 i->axpy(alpha, *j);
1604
1605 return *this;
1606 }
1607
1608 //===== linear maps
1609
1611 template<class X, class Y>
1612 void mv (const X& x, Y& y) const
1613 {
1614#ifdef DUNE_ISTL_WITH_CHECKING
1615 if (ready != built)
1616 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1617 if (x.N()!=M()) DUNE_THROW(BCRSMatrixError,
1618 "Size mismatch: M: " << N() << "x" << M() << " x: " << x.N());
1619 if (y.N()!=N()) DUNE_THROW(BCRSMatrixError,
1620 "Size mismatch: M: " << N() << "x" << M() << " y: " << y.N());
1621#endif
1622 ConstRowIterator endi=end();
1623 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1624 {
1625 y[i.index()]=0;
1626 ConstColIterator endj = (*i).end();
1627 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1628 {
1629 auto&& xj = Impl::asVector(x[j.index()]);
1630 auto&& yi = Impl::asVector(y[i.index()]);
1631 Impl::asMatrix(*j).umv(xj, yi);
1632 }
1633 }
1634 }
1635
1637 template<class X, class Y>
1638 void umv (const X& x, Y& y) const
1639 {
1640#ifdef DUNE_ISTL_WITH_CHECKING
1641 if (ready != built)
1642 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1643 if (x.N()!=M()) DUNE_THROW(BCRSMatrixError,"index out of range");
1644 if (y.N()!=N()) DUNE_THROW(BCRSMatrixError,"index out of range");
1645#endif
1646 ConstRowIterator endi=end();
1647 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1648 {
1649 ConstColIterator endj = (*i).end();
1650 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1651 {
1652 auto&& xj = Impl::asVector(x[j.index()]);
1653 auto&& yi = Impl::asVector(y[i.index()]);
1654 Impl::asMatrix(*j).umv(xj,yi);
1655 }
1656 }
1657 }
1658
1660 template<class X, class Y>
1661 void mmv (const X& x, Y& y) const
1662 {
1663#ifdef DUNE_ISTL_WITH_CHECKING
1664 if (ready != built)
1665 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1666 if (x.N()!=M()) DUNE_THROW(BCRSMatrixError,"index out of range");
1667 if (y.N()!=N()) DUNE_THROW(BCRSMatrixError,"index out of range");
1668#endif
1669 ConstRowIterator endi=end();
1670 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1671 {
1672 ConstColIterator endj = (*i).end();
1673 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1674 {
1675 auto&& xj = Impl::asVector(x[j.index()]);
1676 auto&& yi = Impl::asVector(y[i.index()]);
1677 Impl::asMatrix(*j).mmv(xj,yi);
1678 }
1679 }
1680 }
1681
1683 template<class X, class Y, class F>
1684 void usmv (F&& alpha, const X& x, Y& y) const
1685 {
1686#ifdef DUNE_ISTL_WITH_CHECKING
1687 if (ready != built)
1688 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1689 if (x.N()!=M()) DUNE_THROW(BCRSMatrixError,"index out of range");
1690 if (y.N()!=N()) DUNE_THROW(BCRSMatrixError,"index out of range");
1691#endif
1692 ConstRowIterator endi=end();
1693 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1694 {
1695 ConstColIterator endj = (*i).end();
1696 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1697 {
1698 auto&& xj = Impl::asVector(x[j.index()]);
1699 auto&& yi = Impl::asVector(y[i.index()]);
1700 Impl::asMatrix(*j).usmv(alpha,xj,yi);
1701 }
1702 }
1703 }
1704
1706 template<class X, class Y>
1707 void mtv (const X& x, Y& y) const
1708 {
1709#ifdef DUNE_ISTL_WITH_CHECKING
1710 if (ready != built)
1711 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1712 if (x.N()!=N()) DUNE_THROW(BCRSMatrixError,"index out of range");
1713 if (y.N()!=M()) DUNE_THROW(BCRSMatrixError,"index out of range");
1714#endif
1715 for(size_type i=0; i<y.N(); ++i)
1716 y[i]=0;
1717 umtv(x,y);
1718 }
1719
1721 template<class X, class Y>
1722 void umtv (const X& x, Y& y) const
1723 {
1724#ifdef DUNE_ISTL_WITH_CHECKING
1725 if (ready != built)
1726 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1727 if (x.N()!=N()) DUNE_THROW(BCRSMatrixError,"index out of range");
1728 if (y.N()!=M()) DUNE_THROW(BCRSMatrixError,"index out of range");
1729#endif
1730 ConstRowIterator endi=end();
1731 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1732 {
1733 ConstColIterator endj = (*i).end();
1734 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1735 {
1736 auto&& xi = Impl::asVector(x[i.index()]);
1737 auto&& yj = Impl::asVector(y[j.index()]);
1738 Impl::asMatrix(*j).umtv(xi,yj);
1739 }
1740 }
1741 }
1742
1744 template<class X, class Y>
1745 void mmtv (const X& x, Y& y) const
1746 {
1747#ifdef DUNE_ISTL_WITH_CHECKING
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 ConstRowIterator endi=end();
1752 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1753 {
1754 ConstColIterator endj = (*i).end();
1755 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1756 {
1757 auto&& xi = Impl::asVector(x[i.index()]);
1758 auto&& yj = Impl::asVector(y[j.index()]);
1759 Impl::asMatrix(*j).mmtv(xi,yj);
1760 }
1761 }
1762 }
1763
1765 template<class X, class Y>
1766 void usmtv (const field_type& alpha, const X& x, Y& y) const
1767 {
1768#ifdef DUNE_ISTL_WITH_CHECKING
1769 if (ready != built)
1770 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1771 if (x.N()!=N()) DUNE_THROW(BCRSMatrixError,"index out of range");
1772 if (y.N()!=M()) DUNE_THROW(BCRSMatrixError,"index out of range");
1773#endif
1774 ConstRowIterator endi=end();
1775 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1776 {
1777 ConstColIterator endj = (*i).end();
1778 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1779 {
1780 auto&& xi = Impl::asVector(x[i.index()]);
1781 auto&& yj = Impl::asVector(y[j.index()]);
1782 Impl::asMatrix(*j).usmtv(alpha,xi,yj);
1783 }
1784 }
1785 }
1786
1788 template<class X, class Y>
1789 void umhv (const X& x, Y& y) const
1790 {
1791#ifdef DUNE_ISTL_WITH_CHECKING
1792 if (ready != built)
1793 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1794 if (x.N()!=N()) DUNE_THROW(BCRSMatrixError,"index out of range");
1795 if (y.N()!=M()) DUNE_THROW(BCRSMatrixError,"index out of range");
1796#endif
1797 ConstRowIterator endi=end();
1798 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1799 {
1800 ConstColIterator endj = (*i).end();
1801 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1802 {
1803 auto&& xi = Impl::asVector(x[i.index()]);
1804 auto&& yj = Impl::asVector(y[j.index()]);
1805 Impl::asMatrix(*j).umhv(xi,yj);
1806 }
1807 }
1808 }
1809
1811 template<class X, class Y>
1812 void mmhv (const X& x, Y& y) const
1813 {
1814#ifdef DUNE_ISTL_WITH_CHECKING
1815 if (ready != built)
1816 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1817 if (x.N()!=N()) DUNE_THROW(BCRSMatrixError,"index out of range");
1818 if (y.N()!=M()) DUNE_THROW(BCRSMatrixError,"index out of range");
1819#endif
1820 ConstRowIterator endi=end();
1821 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1822 {
1823 ConstColIterator endj = (*i).end();
1824 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1825 {
1826 auto&& xi = Impl::asVector(x[i.index()]);
1827 auto&& yj = Impl::asVector(y[j.index()]);
1828 Impl::asMatrix(*j).mmhv(xi,yj);
1829 }
1830 }
1831 }
1832
1834 template<class X, class Y>
1835 void usmhv (const field_type& alpha, const X& x, Y& y) const
1836 {
1837#ifdef DUNE_ISTL_WITH_CHECKING
1838 if (ready != built)
1839 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1840 if (x.N()!=N()) DUNE_THROW(BCRSMatrixError,"index out of range");
1841 if (y.N()!=M()) DUNE_THROW(BCRSMatrixError,"index out of range");
1842#endif
1843 ConstRowIterator endi=end();
1844 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1845 {
1846 ConstColIterator endj = (*i).end();
1847 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1848 {
1849 auto&& xi = Impl::asVector(x[i.index()]);
1850 auto&& yj = Impl::asVector(y[j.index()]);
1851 Impl::asMatrix(*j).usmhv(alpha,xi,yj);
1852 }
1853 }
1854 }
1855
1856
1857 //===== norms
1858
1860 typename FieldTraits<field_type>::real_type frobenius_norm2 () const
1861 {
1862#ifdef DUNE_ISTL_WITH_CHECKING
1863 if (ready != built)
1864 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1865#endif
1866
1867 typename FieldTraits<field_type>::real_type sum=0;
1868
1869 for (auto&& row : (*this))
1870 for (auto&& entry : row)
1871 sum += Impl::asMatrix(entry).frobenius_norm2();
1872
1873 return sum;
1874 }
1875
1877 typename FieldTraits<field_type>::real_type frobenius_norm () const
1878 {
1879 return sqrt(frobenius_norm2());
1880 }
1881
1883 template <typename ft = field_type,
1884 typename std::enable_if<!HasNaN<ft>::value, int>::type = 0>
1885 typename FieldTraits<ft>::real_type infinity_norm() const {
1886 if (ready != built)
1887 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1888
1889 using real_type = typename FieldTraits<ft>::real_type;
1890 using std::max;
1891
1892 real_type norm = 0;
1893 for (auto const &x : *this) {
1894 real_type sum = 0;
1895 for (auto const &y : x)
1896 sum += Impl::asMatrix(y).infinity_norm();
1897 norm = max(sum, norm);
1898 }
1899 return norm;
1900 }
1901
1903 template <typename ft = field_type,
1904 typename std::enable_if<!HasNaN<ft>::value, int>::type = 0>
1905 typename FieldTraits<ft>::real_type infinity_norm_real() const {
1906 if (ready != built)
1907 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1908
1909 using real_type = typename FieldTraits<ft>::real_type;
1910 using std::max;
1911
1912 real_type norm = 0;
1913 for (auto const &x : *this) {
1914 real_type sum = 0;
1915 for (auto const &y : x)
1916 sum += Impl::asMatrix(y).infinity_norm_real();
1917 norm = max(sum, norm);
1918 }
1919 return norm;
1920 }
1921
1923 template <typename ft = field_type,
1924 typename std::enable_if<HasNaN<ft>::value, int>::type = 0>
1925 typename FieldTraits<ft>::real_type infinity_norm() const {
1926 if (ready != built)
1927 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1928
1929 using real_type = typename FieldTraits<ft>::real_type;
1930 using std::max;
1931
1932 real_type norm = 0;
1933 real_type isNaN = 1;
1934 for (auto const &x : *this) {
1935 real_type sum = 0;
1936 for (auto const &y : x)
1937 sum += Impl::asMatrix(y).infinity_norm();
1938 norm = max(sum, norm);
1939 isNaN += sum;
1940 }
1941
1942 return norm * (isNaN / isNaN);
1943 }
1944
1946 template <typename ft = field_type,
1947 typename std::enable_if<HasNaN<ft>::value, int>::type = 0>
1948 typename FieldTraits<ft>::real_type infinity_norm_real() const {
1949 if (ready != built)
1950 DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1951
1952 using real_type = typename FieldTraits<ft>::real_type;
1953 using std::max;
1954
1955 real_type norm = 0;
1956 real_type isNaN = 1;
1957
1958 for (auto const &x : *this) {
1959 real_type sum = 0;
1960 for (auto const &y : x)
1961 sum += Impl::asMatrix(y).infinity_norm_real();
1962 norm = max(sum, norm);
1963 isNaN += sum;
1964 }
1965
1966 return norm * (isNaN / isNaN);
1967 }
1968
1969 //===== sizes
1970
1972 size_type N () const
1973 {
1974 return n;
1975 }
1976
1978 size_type M () const
1979 {
1980 return m;
1981 }
1982
1985 {
1986 // in case of row-wise allocation
1987 if( nnz_ <= 0 )
1988 nnz_ = std::accumulate( begin(), end(), size_type( 0 ), [] ( size_type s, const row_type &row ) { return s+row.getsize(); } );
1989 return nnz_;
1990 }
1991
1994 {
1995 return ready;
1996 }
1997
2000 {
2001 return build_mode;
2002 }
2003
2004 //===== query
2005
2007 bool exists (size_type i, size_type j) const
2008 {
2009#ifdef DUNE_ISTL_WITH_CHECKING
2010 if (i<0 || i>=n) DUNE_THROW(BCRSMatrixError,"row index out of range");
2011 if (j<0 || j>=m) DUNE_THROW(BCRSMatrixError,"column index out of range");
2012#endif
2013 return (r[i].size() && r[i].find(j) != r[i].end());
2014 }
2015
2016
2017 protected:
2018 // state information
2019 BuildMode build_mode; // row wise or whole matrix
2020 BuildStage ready; // indicate the stage the matrix building is in
2021
2022 // The allocator used for memory management
2023 typename std::allocator_traits<A>::template rebind_alloc<B> allocator_;
2024
2025 typename std::allocator_traits<A>::template rebind_alloc<row_type> rowAllocator_;
2026
2027 typename std::allocator_traits<A>::template rebind_alloc<size_type> sizeAllocator_;
2028
2029 // size of the matrix
2030 size_type n; // number of rows
2031 size_type m; // number of columns
2032 mutable size_type nnz_; // number of nonzeroes contained in the matrix
2033 size_type allocationSize_; //allocated size of a and j arrays, except in implicit mode: nnz_==allocationsSize_
2034 // zero means that memory is allocated separately for each row.
2035
2036 // the rows are dynamically allocated
2037 row_type* r; // [n] the individual rows having pointers into a,j arrays
2038
2039 // dynamically allocated memory
2040 B* a; // [allocationSize] non-zero entries of the matrix in row-wise ordering
2041 // If a single array of column indices is used, it can be shared
2042 // between different matrices with the same sparsity pattern
2043 std::shared_ptr<size_type> j_; // [allocationSize] column indices of entries
2044
2045 // additional data is needed in implicit buildmode
2046 size_type avg;
2047 double compressionBufferSize_;
2048
2049 typedef std::map<std::pair<size_type,size_type>, B> OverflowType;
2050 OverflowType overflow;
2051
2052 void setWindowPointers(ConstRowIterator row)
2053 {
2054 row_type current_row(a,j_.get(),0); // Pointers to current row data
2055 for (size_type i=0; i<n; i++, ++row) {
2056 // set row i
2057 size_type s = row->getsize();
2058
2059 if (s>0) {
2060 // setup pointers and size
2061 r[i].set(s,current_row.getptr(), current_row.getindexptr());
2062 // update pointer for next row
2063 current_row.setptr(current_row.getptr()+s);
2064 current_row.setindexptr(current_row.getindexptr()+s);
2065 } else{
2066 // empty row
2067 r[i].set(0,nullptr,nullptr);
2068 }
2069 }
2070 }
2071
2073
2078 {
2079 size_type* jptr = j_.get();
2080 for (size_type i=0; i<n; ++i, ++row) {
2081 // set row i
2082 size_type s = row->getsize();
2083
2084 if (s>0) {
2085 // setup pointers and size
2086 r[i].setsize(s);
2087 r[i].setindexptr(jptr);
2088 } else{
2089 // empty row
2090 r[i].set(0,nullptr,nullptr);
2091 }
2092
2093 // advance position in global array
2094 jptr += s;
2095 }
2096 }
2097
2099
2104 {
2105 B* aptr = a;
2106 for (size_type i=0; i<n; ++i) {
2107 // set row i
2108 if (r[i].getsize() > 0) {
2109 // setup pointers and size
2110 r[i].setptr(aptr);
2111 } else{
2112 // empty row
2113 r[i].set(0,nullptr,nullptr);
2114 }
2115
2116 // advance position in global array
2117 aptr += r[i].getsize();
2118 }
2119 }
2120
2123 {
2124 setWindowPointers(Mat.begin());
2125
2126 // copy data
2127 for (size_type i=0; i<n; i++) r[i] = Mat.r[i];
2128
2129 // finish off
2130 build_mode = row_wise; // dummy
2131 ready = built;
2132 }
2133
2139 void deallocate(bool deallocateRows=true)
2140 {
2141
2142 if (notAllocated)
2143 return;
2144
2145 if (allocationSize_>0)
2146 {
2147 // a,j_ have been allocated as one long vector
2148 j_.reset();
2149 if (a)
2150 {
2151 for(B *aiter=a+(allocationSize_-1), *aend=a-1; aiter!=aend; --aiter)
2152 std::allocator_traits<decltype(allocator_)>::destroy(allocator_, aiter);
2153 allocator_.deallocate(a,allocationSize_);
2154 a = nullptr;
2155 }
2156 }
2157 else if (r)
2158 {
2159 // check if memory for rows have been allocated individually
2160 for (size_type i=0; i<n; i++)
2161 if (r[i].getsize()>0)
2162 {
2163 for (B *col=r[i].getptr()+(r[i].getsize()-1),
2164 *colend = r[i].getptr()-1; col!=colend; --col) {
2165 std::allocator_traits<decltype(allocator_)>::destroy(allocator_, col);
2166 }
2167 sizeAllocator_.deallocate(r[i].getindexptr(),1);
2168 allocator_.deallocate(r[i].getptr(),1);
2169 // clear out row data in case we don't want to deallocate the rows
2170 // otherwise we might run into a double free problem here later
2171 r[i].set(0,nullptr,nullptr);
2172 }
2173 }
2174
2175 // deallocate the rows
2176 if (n>0 && deallocateRows && r) {
2177 for(row_type *riter=r+(n-1), *rend=r-1; riter!=rend; --riter)
2178 std::allocator_traits<decltype(rowAllocator_)>::destroy(rowAllocator_, riter);
2179 rowAllocator_.deallocate(r,n);
2180 r = nullptr;
2181 }
2182
2183 // Mark matrix as not built at all.
2184 ready=notAllocated;
2185
2186 }
2187
2205 void allocate(size_type rows, size_type columns, size_type allocationSize, bool allocateRows, bool allocate_data)
2206 {
2207 // Store size
2208 n = rows;
2209 m = columns;
2210 nnz_ = allocationSize;
2211 allocationSize_ = allocationSize;
2212
2213 // allocate rows
2214 if(allocateRows) {
2215 if (n>0) {
2216 if (r)
2217 DUNE_THROW(InvalidStateException,"Rows have already been allocated, cannot allocate a second time");
2218 r = rowAllocator_.allocate(rows);
2219 // initialize row entries
2220 for(row_type* ri=r; ri!=r+rows; ++ri)
2221 std::allocator_traits<decltype(rowAllocator_)>::construct(rowAllocator_, ri, row_type());
2222 }else{
2223 r = 0;
2224 }
2225 }
2226
2227 // allocate a and j_ array
2228 if (allocate_data)
2229 allocateData();
2230 // allocate column indices only if not yet present (enable sharing)
2231 if (allocationSize_>0) {
2232 // we copy allocator and size to the deleter since _j may outlive this class
2233 if (!j_.get())
2234 j_.reset(sizeAllocator_.allocate(allocationSize_),
2235 [alloc = sizeAllocator_, size = allocationSize_](auto ptr) mutable {
2236 alloc.deallocate(ptr, size);
2237 });
2238 }else{
2239 j_.reset();
2240 }
2241
2242 // Mark the matrix as not built.
2243 ready = building;
2244 }
2245
2246 void allocateData()
2247 {
2248 if (a)
2249 DUNE_THROW(InvalidStateException,"Cannot allocate data array (already allocated)");
2250 if (allocationSize_>0) {
2251 a = allocator_.allocate(allocationSize_);
2252 // use placement new to call constructor that allocates
2253 // additional memory.
2254 new (a) B[allocationSize_];
2255 } else {
2256 a = nullptr;
2257 }
2258 }
2259
2266 {
2267 if (build_mode != implicit)
2268 DUNE_THROW(InvalidStateException,"implicit_allocate() may only be called in implicit build mode");
2269 if (ready != notAllocated)
2270 DUNE_THROW(InvalidStateException,"memory has already been allocated");
2271
2272 // check to make sure the user has actually set the parameters
2273 if (compressionBufferSize_ < 0)
2274 DUNE_THROW(InvalidStateException,"You have to set the implicit build mode parameters before starting to build the matrix");
2275 //calculate size of overflow region, add buffer for row sort!
2276 size_type osize = (size_type) (_n*avg)*compressionBufferSize_ + 4*avg;
2277 allocationSize_ = _n*avg + osize;
2278
2279 allocate(_n, _m, allocationSize_,true,true);
2280
2281 //set row pointers correctly
2282 size_type* jptr = j_.get() + osize;
2283 B* aptr = a + osize;
2284 for (size_type i=0; i<n; i++)
2285 {
2286 r[i].set(0,aptr,jptr);
2287 jptr = jptr + avg;
2288 aptr = aptr + avg;
2289 }
2290
2291 ready = building;
2292 }
2293 };
2294
2295
2296 template<class B, class A>
2297 struct FieldTraits< BCRSMatrix<B, A> >
2298 {
2299 using field_type = typename BCRSMatrix<B, A>::field_type;
2300 using real_type = typename FieldTraits<field_type>::real_type;
2301 };
2302
2305} // end namespace
2306
2307#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:1052
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:1058
bool operator!=(const CreateIterator &it) const
inequality
Definition: bcrsmatrix.hh:1046
CreateIterator(BCRSMatrix &_Mat, size_type _i)
constructor
Definition: bcrsmatrix.hh:960
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:579
RealRowIterator()
empty constructor, use with care!
Definition: bcrsmatrix.hh:596
bool equals(const RealRowIterator< ValueType > &other) const
equality
Definition: bcrsmatrix.hh:624
std::remove_const< T >::type ValueType
The unqualified value type.
Definition: bcrsmatrix.hh:583
bool equals(const RealRowIterator< const ValueType > &other) const
equality
Definition: bcrsmatrix.hh:631
RealRowIterator(row_type *_p, size_type _i)
constructor
Definition: bcrsmatrix.hh:591
size_type index() const
return index
Definition: bcrsmatrix.hh:606
A sparse block matrix with compressed row storage.
Definition: bcrsmatrix.hh:466
typename Imp::BlockTraits< B >::field_type field_type
export the type representing the field
Definition: bcrsmatrix.hh:488
bool exists(size_type i, size_type j) const
return true if (i,j) is in pattern
Definition: bcrsmatrix.hh:2007
BuildStage buildStage() const
The current build stage of the matrix.
Definition: bcrsmatrix.hh:1993
Iterator begin()
Get iterator to first row.
Definition: bcrsmatrix.hh:675
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:2122
B & entry(size_type row, size_type col)
Returns reference to entry (row,col) of the matrix.
Definition: bcrsmatrix.hh:1296
void mmhv(const X &x, Y &y) const
y -= A^H x
Definition: bcrsmatrix.hh:1812
void usmhv(const field_type &alpha, const X &x, Y &y) const
y += alpha A^H x
Definition: bcrsmatrix.hh:1835
void umtv(const X &x, Y &y) const
y += A^T x
Definition: bcrsmatrix.hh:1722
RealRowIterator< const row_type > const_iterator
The const iterator over the matrix rows.
Definition: bcrsmatrix.hh:707
static constexpr unsigned int blocklevel
increment block level counter
Definition: bcrsmatrix.hh:507
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:2205
BCRSMatrix & axpy(field_type alpha, const BCRSMatrix &b)
Add the scaled entries of another matrix to this one.
Definition: bcrsmatrix.hh:1592
FieldTraits< ft >::real_type infinity_norm_real() const
simplified infinity norm (uses Manhattan norm for complex values)
Definition: bcrsmatrix.hh:1905
~BCRSMatrix()
destructor
Definition: bcrsmatrix.hh:824
void deallocate(bool deallocateRows=true)
deallocate memory of the matrix.
Definition: bcrsmatrix.hh:2139
Iterator RowIterator
rename the iterators for easier access
Definition: bcrsmatrix.hh:701
row_type & operator[](size_type i)
random access to the rows
Definition: bcrsmatrix.hh:549
BCRSMatrix()
an empty matrix
Definition: bcrsmatrix.hh:749
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:1707
void umhv(const X &x, Y &y) const
y += A^H x
Definition: bcrsmatrix.hh:1789
size_type nonzeroes() const
number of blocks that are stored (the number of blocks that possibly are nonzero)
Definition: bcrsmatrix.hh:1984
ConstIterator ConstRowIterator
rename the const row iterator for easier access
Definition: bcrsmatrix.hh:738
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:1484
RealRowIterator< row_type > iterator
The iterator over the (mutable matrix rows.
Definition: bcrsmatrix.hh:671
void usmtv(const field_type &alpha, const X &x, Y &y) const
y += alpha A^T x
Definition: bcrsmatrix.hh:1766
ConstIterator beforeBegin() const
Definition: bcrsmatrix.hh:732
Iterator beforeBegin()
Definition: bcrsmatrix.hh:695
BuildMode
we support two modes
Definition: bcrsmatrix.hh:510
@ implicit
Build entries randomly with an educated guess for the number of entries per row.
Definition: bcrsmatrix.hh:539
@ unknown
Build mode not set!
Definition: bcrsmatrix.hh:543
@ random
Build entries randomly.
Definition: bcrsmatrix.hh:530
@ row_wise
Build in a row-wise manner.
Definition: bcrsmatrix.hh:521
BCRSMatrix(size_type _n, size_type _m, size_type _nnz, BuildMode bm)
matrix with known number of nonzeroes
Definition: bcrsmatrix.hh:756
::Dune::CompressionStatistics< size_type > CompressionStatistics
The type for the statistics object returned by compress()
Definition: bcrsmatrix.hh:503
BCRSMatrix & operator-=(const BCRSMatrix &b)
Subtract the entries of another matrix from this one.
Definition: bcrsmatrix.hh:1567
BCRSMatrix(const BCRSMatrix &Mat)
copy constructor
Definition: bcrsmatrix.hh:805
Iterator end()
Get iterator to one beyond last row.
Definition: bcrsmatrix.hh:681
void setIndices(size_type row, It begin, It end)
Set all column indices for row from the given iterator range.
Definition: bcrsmatrix.hh:1234
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:704
FieldTraits< field_type >::real_type frobenius_norm() const
frobenius norm: sqrt(sum over squared values of entries)
Definition: bcrsmatrix.hh:1877
A::size_type size_type
The type for the index access and the size.
Definition: bcrsmatrix.hh:500
BCRSMatrix & operator/=(const field_type &k)
vector space division by scalar
Definition: bcrsmatrix.hh:1512
BCRSMatrix & operator+=(const BCRSMatrix &b)
Add the entries of another matrix to this one.
Definition: bcrsmatrix.hh:1545
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:784
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:1860
Iterator beforeEnd()
Definition: bcrsmatrix.hh:688
row_type::ConstIterator ConstColIterator
Const iterator to the entries of a row.
Definition: bcrsmatrix.hh:741
void usmv(F &&alpha, const X &x, Y &y) const
y += alpha A x
Definition: bcrsmatrix.hh:1684
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:1978
Imp::CompressedBlockVectorWindow< B, A > row_type
implement row_type with compressed vector
Definition: bcrsmatrix.hh:497
CreateIterator createbegin()
get initial create iterator
Definition: bcrsmatrix.hh:1097
BuildStage
Definition: bcrsmatrix.hh:469
@ rowSizesBuilt
The row sizes of the matrix are known.
Definition: bcrsmatrix.hh:480
@ built
The matrix structure is fully built.
Definition: bcrsmatrix.hh:482
@ notbuilt
Matrix is not built at all, no memory has been allocated, build mode and size can still be set.
Definition: bcrsmatrix.hh:471
@ notAllocated
Matrix is not built at all, no memory has been allocated, build mode and size can still be set.
Definition: bcrsmatrix.hh:473
@ building
Matrix is currently being built, some memory has been allocated, build mode and size are fixed.
Definition: bcrsmatrix.hh:475
BuildMode buildMode() const
The currently selected build mode of the matrix.
Definition: bcrsmatrix.hh:1999
void mmv(const X &x, Y &y) const
y -= A x
Definition: bcrsmatrix.hh:1661
FieldTraits< ft >::real_type infinity_norm() const
infinity norm (row sum norm, how to generalize for blocks?)
Definition: bcrsmatrix.hh:1885
void mv(const X &x, Y &y) const
y = A x
Definition: bcrsmatrix.hh:1612
B block_type
export the type representing the components
Definition: bcrsmatrix.hh:491
void mmtv(const X &x, Y &y) const
y -= A^T x
Definition: bcrsmatrix.hh:1745
void umv(const X &x, Y &y) const
y += A x
Definition: bcrsmatrix.hh:1638
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:2265
void setImplicitBuildModeParameters(size_type _avg, double compressionBufferSize)
Set parameters needed for creation in implicit build mode.
Definition: bcrsmatrix.hh:889
BCRSMatrix(size_type _n, size_type _m, BuildMode bm)
matrix with unknown number of nonzeroes
Definition: bcrsmatrix.hh:765
void endindices()
indicate that all indices are defined, check consistency
Definition: bcrsmatrix.hh:1248
CompressionStatistics compress()
Finishes the buildstage in implicit mode.
Definition: bcrsmatrix.hh:1360
void setDataPointers()
Set data pointers for all rows.
Definition: bcrsmatrix.hh:2103
size_type N() const
number of rows (counted in blocks)
Definition: bcrsmatrix.hh:1972
void setBuildMode(BuildMode bm)
Sets the build mode of the matrix.
Definition: bcrsmatrix.hh:833
void setSize(size_type rows, size_type columns, size_type nnz=0)
Set the size of the matrix.
Definition: bcrsmatrix.hh:861
BCRSMatrix & operator=(const BCRSMatrix &Mat)
assignment
Definition: bcrsmatrix.hh:911
void setColumnPointers(ConstRowIterator row)
Copy row sizes from iterator range starting at row and set column index pointers for all rows.
Definition: bcrsmatrix.hh:2077
ConstIterator end() const
Get const iterator to one beyond last row.
Definition: bcrsmatrix.hh:718
ConstIterator begin() const
Get const iterator to first row.
Definition: bcrsmatrix.hh:712
A allocator_type
export the allocator type
Definition: bcrsmatrix.hh:494
ConstIterator beforeEnd() const
Definition: bcrsmatrix.hh:725
Proxy row object for entry access.
Definition: bcrsmatrix.hh:137
block_type & operator[](size_type j) const
Returns entry in column j.
Definition: bcrsmatrix.hh:142
A wrapper for uniform access to the BCRSMatrix during and after the build stage in implicit build mod...
Definition: bcrsmatrix.hh:117
Matrix::block_type block_type
The block_type of the underlying matrix.
Definition: bcrsmatrix.hh:125
ImplicitMatrixBuilder(Matrix &m)
Creates an ImplicitMatrixBuilder for matrix m.
Definition: bcrsmatrix.hh:170
M_ Matrix
The underlying matrix.
Definition: bcrsmatrix.hh:122
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:194
size_type M() const
The number of columns in the matrix.
Definition: bcrsmatrix.hh:217
Matrix::size_type size_type
The size_type of the underlying matrix.
Definition: bcrsmatrix.hh:128
row_object operator[](size_type i) const
Returns a proxy for entries in row i.
Definition: bcrsmatrix.hh:205
size_type N() const
The number of rows in the matrix.
Definition: bcrsmatrix.hh:211
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:281
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:434
Default exception class for range errors.
Definition: exceptions.hh:254
Traits for type conversions and type information.
Type traits to determine the type of reals (when working with complex numbers)
#define DUNE_THROW(E, m)
Definition: exceptions.hh:218
constexpr decltype(auto) elementAt(Container &&c, Index &&i)
Get element at given position from container.
Definition: hybridutilities.hh:135
constexpr T accumulate(Range &&range, T value, F &&f)
Accumulate values.
Definition: hybridutilities.hh:291
auto max(ADLTag< 0 >, const V &v1, const V &v2)
implements binary Simd::max()
Definition: defaults.hh:81
DWarnType dwarn(std::cerr)
Stream for warnings indicating problems.
Definition: stdstreams.hh:161
This file implements iterator facade classes for writing stl conformant iterators.
Some handy generic functions for ISTL matrices.
Dune namespace.
Definition: alignedallocator.hh:13
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:88
size_type overflow_total
total number of elements written to the overflow area during construction.
Definition: bcrsmatrix.hh:94
size_type maximum
maximum number of non-zeroes per row.
Definition: bcrsmatrix.hh:92
double avg
average number of non-zeroes per row.
Definition: bcrsmatrix.hh:90
double mem_ratio
fraction of wasted memory resulting from non-used overflow area.
Definition: bcrsmatrix.hh:99
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Dec 21, 23:30, 2024)