Dune Core Modules (2.7.1)

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