dune-istl  2.3beta2
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_BCRSMATRIX_HH
5 #define DUNE_BCRSMATRIX_HH
6 
7 #include <cmath>
8 #include <complex>
9 #include <set>
10 #include <iostream>
11 #include <algorithm>
12 #include <numeric>
13 #include <vector>
14 #include <map>
15 
16 #include "istlexception.hh"
17 #include "bvector.hh"
18 #include "matrixutils.hh"
19 #include <dune/common/shared_ptr.hh>
20 #include <dune/common/stdstreams.hh>
21 #include <dune/common/iteratorfacades.hh>
22 #include <dune/common/typetraits.hh>
23 #include <dune/common/ftraits.hh>
24 #include <dune/common/static_assert.hh>
25 
30 namespace Dune {
31 
71  template<typename M>
72  struct MatrixDimension;
73 
75 
81  template<typename size_type>
83  {
85  double avg;
87  size_type maximum;
89  size_type overflow_total;
91 
94  double mem_ratio;
95  };
96 
98 
110  template<class M_>
112  {
113 
114  public:
115 
117  typedef M_ Matrix;
118 
120  typedef typename Matrix::block_type block_type;
121 
123  typedef typename Matrix::size_type size_type;
124 
126 
132  {
133 
134  public:
135 
138  {
139  return _m.entry(_i,j);
140  }
141 
142 #ifndef DOXYGEN
143 
145  : _m(m)
146  , _i(i)
147  {}
148 
149 #endif
150 
151  private:
152 
153  Matrix& _m;
154  size_type _i;
155 
156  };
157 
159 
166  : _m(m)
167  {
168  if (m.buildMode() != Matrix::implicit)
169  DUNE_THROW(BCRSMatrixError,"You can only create an ImplicitBuilder for a matrix in implicit build mode");
170  if (m.buildStage() != Matrix::building)
171  DUNE_THROW(BCRSMatrixError,"You can only create an ImplicitBuilder for a matrix with set size that has not been compressed() yet");
172  }
173 
175 
189  ImplicitMatrixBuilder(Matrix& m, size_type rows, size_type cols, size_type avg_cols_per_row, double overflow_fraction)
190  : _m(m)
191  {
192  if (m.buildStage() != Matrix::notAllocated)
193  DUNE_THROW(BCRSMatrixError,"You can only set up a matrix for this ImplicitBuilder if it has no memory allocated yet");
194  m.setBuildMode(Matrix::implicit);
195  m.setImplicitBuildModeParameters(avg_cols_per_row,overflow_fraction);
196  m.setSize(rows,cols);
197  }
198 
201  {
202  return row_object(_m,i);
203  }
204 
206  size_type N() const
207  {
208  return _m.N();
209  }
210 
212  size_type M() const
213  {
214  return _m.M();
215  }
216 
217  private:
218 
219  Matrix& _m;
220 
221  };
222 
412  template<class B, class A=std::allocator<B> >
414  {
415  friend struct MatrixDimension<BCRSMatrix>;
416  public:
417  enum BuildStage {
430  built=3,
431  };
432 
433  //===== type definitions and constants
434 
436  typedef typename B::field_type field_type;
437 
439  typedef B block_type;
440 
442  typedef A allocator_type;
443 
446 
448  typedef typename A::size_type size_type;
449 
451  typedef ::Dune::CompressionStatistics<size_type> CompressionStatistics;
452 
454  enum {
456  blocklevel = B::blocklevel+1
457  };
458 
460  enum BuildMode {
494  };
495 
496  //===== random access interface to rows of the matrix
497 
500  {
501 #ifdef DUNE_ISTL_WITH_CHECKING
502  if (build_mode == implicit && ready != built)
503  DUNE_THROW(BCRSMatrixError,"You cannot use operator[] in implicit build mode before calling compress()");
504  if (r==0) DUNE_THROW(BCRSMatrixError,"row not initialized yet");
505  if (i>=n) DUNE_THROW(BCRSMatrixError,"index out of range");
506  if (r[i].getptr()==0) DUNE_THROW(BCRSMatrixError,"row not initialized yet");
507 #endif
508  return r[i];
509  }
510 
513  {
514 #ifdef DUNE_ISTL_WITH_CHECKING
515  if (build_mode == implicit && ready != built)
516  DUNE_THROW(BCRSMatrixError,"You cannot use operator[] in implicit build mode before calling compress()");
517  if (built!=ready) DUNE_THROW(BCRSMatrixError,"row not initialized yet");
518  if (i>=n) DUNE_THROW(BCRSMatrixError,"index out of range");
519 #endif
520  return r[i];
521  }
522 
523 
524  //===== iterator interface to rows of the matrix
525 
527  template<class T>
529  : public RandomAccessIteratorFacade<RealRowIterator<T>, T>
530  {
531 
532  public:
535 
536  friend class RandomAccessIteratorFacade<RealRowIterator<const ValueType>, const ValueType>;
537  friend class RandomAccessIteratorFacade<RealRowIterator<ValueType>, ValueType>;
538  friend class RealRowIterator<const ValueType>;
539  friend class RealRowIterator<ValueType>;
540 
543  : p(_p), i(_i)
544  {}
545 
548  : p(0), i(0)
549  {}
550 
552  : p(it.p), i(it.i)
553  {}
554 
555 
557  size_type index () const
558  {
559  return i;
560  }
561 
562  std::ptrdiff_t distanceTo(const RealRowIterator<ValueType>& other) const
563  {
564  assert(other.p==p);
565  return (other.i-i);
566  }
567 
568  std::ptrdiff_t distanceTo(const RealRowIterator<const ValueType>& other) const
569  {
570  assert(other.p==p);
571  return (other.i-i);
572  }
573 
575  bool equals (const RealRowIterator<ValueType>& other) const
576  {
577  assert(other.p==p);
578  return i==other.i;
579  }
580 
582  bool equals (const RealRowIterator<const ValueType>& other) const
583  {
584  assert(other.p==p);
585  return i==other.i;
586  }
587 
588  private:
590  void increment()
591  {
592  ++i;
593  }
594 
596  void decrement()
597  {
598  --i;
599  }
600 
601  void advance(std::ptrdiff_t diff)
602  {
603  i+=diff;
604  }
605 
606  T& elementAt(std::ptrdiff_t diff) const
607  {
608  return p[i+diff];
609  }
610 
612  row_type& dereference () const
613  {
614  return p[i];
615  }
616 
617  row_type* p;
618  size_type i;
619  };
620 
624 
627  {
628  return Iterator(r,0);
629  }
630 
633  {
634  return Iterator(r,n);
635  }
636 
640  {
641  return Iterator(r,n-1);
642  }
643 
647  {
648  return Iterator(r,-1);
649  }
650 
653 
655  typedef typename row_type::Iterator ColIterator;
656 
660 
661 
664  {
665  return ConstIterator(r,0);
666  }
667 
670  {
671  return ConstIterator(r,n);
672  }
673 
677  {
678  return ConstIterator(r,n-1);
679  }
680 
684  {
685  return ConstIterator(r,-1);
686  }
687 
690 
693 
694  //===== constructors & resizers
695 
696  // we use a negative overflowsize to indicate that the implicit
697  // mode parameters have not been set yet
698 
701  : build_mode(unknown), ready(notAllocated), n(0), m(0), nnz(0),
702  allocationSize(0), r(0), a(0),
703  avg(0), overflowsize(-1.0)
704  {}
705 
708  : build_mode(bm), ready(notAllocated), n(0), m(0), nnz(0),
709  allocationSize(0), r(0), a(0),
710  avg(0), overflowsize(-1.0)
711  {
712  allocate(_n, _m, _nnz,true,false);
713  }
714 
717  : build_mode(bm), ready(notAllocated), n(0), m(0), nnz(0),
718  allocationSize(0), r(0), a(0),
719  avg(0), overflowsize(-1.0)
720  {
721  allocate(_n, _m,0,true,false);
722  }
723 
725 
735  BCRSMatrix (size_type _n, size_type _m, size_type _avg, double _overflowsize, BuildMode bm)
736  : build_mode(bm), ready(notAllocated), n(0), m(0), nnz(0),
737  allocationSize(0), r(0), a(0),
738  avg(_avg), overflowsize(_overflowsize)
739  {
740  if (bm != implicit)
741  DUNE_THROW(BCRSMatrixError,"Only call this constructor when using the implicit build mode");
742  // Prevent user from setting a negative overflowsize:
743  // 1) It doesn't make sense
744  // 2) We use a negative overflow value to indicate that the parameters
745  // have not been set yet
746  if (_overflowsize < 0.0)
747  DUNE_THROW(BCRSMatrixError,"You cannot set a negative overflow fraction");
748  implicit_allocate(_n,_m);
749  }
750 
756  BCRSMatrix (const BCRSMatrix& Mat)
757  : build_mode(Mat.build_mode), ready(notAllocated), n(0), m(0), nnz(0),
758  allocationSize(0), r(0), a(0),
759  avg(Mat.avg), overflowsize(Mat.overflowsize)
760  {
761  if (!(Mat.ready == notAllocated || Mat.ready == built))
762  DUNE_THROW(InvalidStateException,"BCRSMatrix can only be copy-constructed when source matrix is completely empty (size not set) or fully built)");
763 
764  // deep copy in global array
765  size_type _nnz = Mat.nnz;
766 
767  // in case of row-wise allocation
768  if (_nnz<=0)
769  {
770  _nnz = 0;
771  for (size_type i=0; i<Mat.n; i++)
772  _nnz += Mat.r[i].getsize();
773  }
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  deallocate(false);
880 
881  // reallocate the rows if required
882  if (n>0 && n!=Mat.n) {
883  // free rows
884  for(row_type *riter=r+(n-1), *rend=r-1; riter!=rend; --riter)
885  rowAllocator_.destroy(riter);
886  rowAllocator_.deallocate(r,n);
887  }
888 
889  nnz=Mat.nnz;
890  if (nnz<=0)
891  {
892  for (size_type i=0; i<Mat.n; i++)
893  nnz += Mat.r[i].getsize();
894  }
895 
896  // allocate a, share j
897  j = Mat.j;
898  allocate(Mat.n, Mat.m, nnz, n!=Mat.n, true);
899 
900  // build window structure
901  copyWindowStructure(Mat);
902  return *this;
903  }
904 
907  {
908 
909  if (!(ready == notAllocated || ready == built))
910  DUNE_THROW(InvalidStateException,"Scalar assignment only works on fully built BCRSMatrix)");
911 
912  for (size_type i=0; i<n; i++) r[i] = k;
913  return *this;
914  }
915 
916  //===== row-wise creation interface
917 
920  {
921  public:
924  : Mat(_Mat), i(_i), nnz(0), current_row(nullptr, Mat.j.get(), 0)
925  {
926  if (Mat.build_mode == unknown && Mat.ready == building)
927  {
928  Mat.build_mode = row_wise;
929  }
930  if (i==0 && Mat.ready != building)
931  DUNE_THROW(BCRSMatrixError,"creation only allowed for uninitialized matrix");
932  if(Mat.build_mode!=row_wise)
933  DUNE_THROW(BCRSMatrixError,"creation only allowed if row wise allocation was requested in the constructor");
934  }
935 
938  {
939  // this should only be called if matrix is in creation
940  if (Mat.ready != building)
941  DUNE_THROW(BCRSMatrixError,"matrix already built up");
942 
943  // row i is defined through the pattern
944  // get memory for the row and initialize the j array
945  // this depends on the allocation mode
946 
947  // compute size of the row
948  size_type s = pattern.size();
949 
950  if(s>0) {
951  // update number of nonzeroes including this row
952  nnz += s;
953 
954  // alloc memory / set window
955  if (Mat.nnz>0)
956  {
957  // memory is allocated in one long array
958 
959  // check if that memory is sufficient
960  if (nnz>Mat.nnz)
961  DUNE_THROW(BCRSMatrixError,"allocated nnz too small");
962 
963  // set row i
964  Mat.r[i].set(s,nullptr,current_row.getindexptr());
965  current_row.setindexptr(current_row.getindexptr()+s);
966  }else{
967  // memory is allocated individually per row
968  // allocate and set row i
969  B* a = Mat.allocator_.allocate(s);
970  // use placement new to call constructor that allocates
971  // additional memory.
972  new (a) B[s];
973  size_type* j = Mat.sizeAllocator_.allocate(s);
974  Mat.r[i].set(s,a,j);
975  }
976  }else
977  // setup empty row
978  Mat.r[i].set(0,0,0);
979 
980  // initialize the j array for row i from pattern
981  size_type k=0;
982  size_type *j = Mat.r[i].getindexptr();
983  for (typename PatternType::const_iterator it=pattern.begin(); it!=pattern.end(); ++it)
984  j[k++] = *it;
985 
986  // now go to next row
987  i++;
988  pattern.clear();
989 
990  // check if this was last row
991  if (i==Mat.n)
992  {
993  Mat.ready = built;
994  if(Mat.nnz>0)
995  {
996  // Set nnz to the exact number of nonzero blocks inserted
997  // as some methods rely on it
998  Mat.nnz=nnz;
999  // allocate data array
1000  Mat.allocateData();
1001  Mat.setDataPointers();
1002  }
1003  }
1004  // done
1005  return *this;
1006  }
1007 
1009  bool operator!= (const CreateIterator& it) const
1010  {
1011  return (i!=it.i) || (&Mat!=&it.Mat);
1012  }
1013 
1015  bool operator== (const CreateIterator& it) const
1016  {
1017  return (i==it.i) && (&Mat==&it.Mat);
1018  }
1019 
1021  size_type index () const
1022  {
1023  return i;
1024  }
1025 
1028  {
1029  pattern.insert(j);
1030  }
1031 
1034  {
1035  if (pattern.find(j)!=pattern.end())
1036  return true;
1037  else
1038  return false;
1039  }
1045  size_type size() const
1046  {
1047  return pattern.size();
1048  }
1049 
1050  private:
1051  BCRSMatrix& Mat; // the matrix we are defining
1052  size_type i; // current row to be defined
1053  size_type nnz; // count total number of nonzeros
1054  typedef std::set<size_type,std::less<size_type> > PatternType;
1055  PatternType pattern; // used to compile entries in a row
1056  row_type current_row; // row pointing to the current row to setup
1057  };
1058 
1060  friend class CreateIterator;
1061 
1064  {
1065  return CreateIterator(*this,0);
1066  }
1067 
1070  {
1071  return CreateIterator(*this,n);
1072  }
1073 
1074 
1075  //===== random creation interface
1076 
1079  {
1080  if (build_mode!=random)
1081  DUNE_THROW(BCRSMatrixError,"requires random build mode");
1082  if (ready != building)
1083  DUNE_THROW(BCRSMatrixError,"matrix row sizes already built up");
1084 
1085  r[i].setsize(s);
1086  }
1087 
1090  {
1091 #ifdef DUNE_ISTL_WITH_CHECKING
1092  if (r==0) DUNE_THROW(BCRSMatrixError,"row not initialized yet");
1093  if (i>=n) DUNE_THROW(BCRSMatrixError,"index out of range");
1094 #endif
1095  return r[i].getsize();
1096  }
1097 
1100  {
1101  if (build_mode!=random)
1102  DUNE_THROW(BCRSMatrixError,"requires random build mode");
1103  if (ready != building)
1104  DUNE_THROW(BCRSMatrixError,"matrix row sizes already built up");
1105 
1106  r[i].setsize(r[i].getsize()+s);
1107  }
1108 
1110  void endrowsizes ()
1111  {
1112  if (build_mode!=random)
1113  DUNE_THROW(BCRSMatrixError,"requires random build mode");
1114  if (ready != building)
1115  DUNE_THROW(BCRSMatrixError,"matrix row sizes already built up");
1116 
1117  // compute total size, check positivity
1118  size_type total=0;
1119  for (size_type i=0; i<n; i++)
1120  {
1121  total += r[i].getsize();
1122  }
1123 
1124  if(nnz==0)
1125  // allocate/check memory
1126  allocate(n,m,total,false,false);
1127  else if(nnz<total)
1128  DUNE_THROW(BCRSMatrixError,"Specified number of nonzeros ("<<nnz<<") not "
1129  <<"sufficient for calculated nonzeros ("<<total<<"! ");
1130 
1131  // set the window pointers correctly
1132  setColumnPointers(begin());
1133 
1134  // initialize j array with m (an invalid column index)
1135  // this indicates an unused entry
1136  for (size_type k=0; k<nnz; k++)
1137  j.get()[k] = m;
1138  ready = rowSizesBuilt;
1139  }
1140 
1142 
1153  {
1154  if (build_mode!=random)
1155  DUNE_THROW(BCRSMatrixError,"requires random build mode");
1156  if (ready==built)
1157  DUNE_THROW(BCRSMatrixError,"matrix already built up");
1158  if (ready==building)
1159  DUNE_THROW(BCRSMatrixError,"matrix row sizes not built up yet");
1160  if (ready==notAllocated)
1161  DUNE_THROW(BCRSMatrixError,"matrix size not set and no memory allocated yet");
1162 
1163  if (col >= m)
1164  DUNE_THROW(BCRSMatrixError,"column index exceeds matrix size");
1165 
1166  // get row range
1167  size_type* const first = r[row].getindexptr();
1168  size_type* const last = first + r[row].getsize();
1169 
1170  // find correct insertion position for new column index
1171  size_type* pos = std::lower_bound(first,last,col);
1172 
1173  // check if index is already in row
1174  if (pos!=last && *pos == col) return;
1175 
1176  // find end of already inserted column indices
1177  size_type* end = std::lower_bound(pos,last,m);
1178  if (end==last)
1179  DUNE_THROW(BCRSMatrixError,"row is too small");
1180 
1181  // insert new column index at correct position
1182  std::copy_backward(pos,end,end+1);
1183  *pos = col;
1184  }
1185 
1187 
1194  template<typename It>
1196  {
1197  size_type row_size = r[row].size();
1198  size_type* col_begin = r[row].getindexptr();
1199  size_type* col_end;
1200  // consistency check between allocated row size and number of passed column indices
1201  if ((col_end = std::copy(begin,end,r[row].getindexptr())) != col_begin + row_size)
1202  DUNE_THROW(BCRSMatrixError,"Given size of row " << row
1203  << " (" << row_size
1204  << ") does not match number of passed entries (" << (col_end - col_begin) << ")");
1205  std::sort(col_begin,col_end);
1206  }
1207 
1209  void endindices ()
1210  {
1211  if (build_mode!=random)
1212  DUNE_THROW(BCRSMatrixError,"requires random build mode");
1213  if (ready==built)
1214  DUNE_THROW(BCRSMatrixError,"matrix already built up");
1215  if (ready==building)
1216  DUNE_THROW(BCRSMatrixError,"row sizes are not built up yet");
1217  if (ready==notAllocated)
1218  DUNE_THROW(BCRSMatrixError,"matrix size not set and no memory allocated yet");
1219 
1220  // check if there are undefined indices
1221  RowIterator endi=end();
1222  for (RowIterator i=begin(); i!=endi; ++i)
1223  {
1224  ColIterator endj = (*i).end();
1225  for (ColIterator j=(*i).begin(); j!=endj; ++j) {
1226  if (j.index()>=m) {
1227  dwarn << "WARNING: size of row "<< i.index()<<" is "<<j.offset()<<". But was specified as being "<< (*i).end().offset()
1228  <<". This means you are wasting valuable space and creating additional cache misses!"<<std::endl;
1229  r[i.index()].setsize(j.offset());
1230  break;
1231  }
1232  }
1233  }
1234 
1235  allocateData();
1236  setDataPointers();
1237 
1238  // if not, set matrix to built
1239  ready = built;
1240  }
1241 
1242  //===== implicit creation interface
1243 
1245 
1257  {
1258 #ifdef DUNE_ISTL_WITH_CHECKING
1259  if (build_mode!=implicit)
1260  DUNE_THROW(BCRSMatrixError,"requires implicit build mode");
1261  if (ready==built)
1262  DUNE_THROW(BCRSMatrixError,"matrix already built up, use operator[] for entry access now");
1263  if (ready==notAllocated)
1264  DUNE_THROW(BCRSMatrixError,"matrix size not set and no memory allocated yet");
1265  if (ready!=building)
1266  DUNE_THROW(InvalidStateException,"You may only use entry() during the 'building' stage");
1267 
1268  if (row >= n)
1269  DUNE_THROW(BCRSMatrixError,"row index exceeds matrix size");
1270  if (col >= m)
1271  DUNE_THROW(BCRSMatrixError,"column index exceeds matrix size");
1272 #endif
1273 
1274  size_type* begin = r[row].getindexptr();
1275  size_type* end = begin + r[row].getsize();
1276 
1277  size_type* pos = std::find(begin, end, col);
1278 
1279  //treat the case that there was a match in the array
1280  if (pos != end)
1281  if (*pos == col)
1282  {
1283  std::ptrdiff_t offset = pos - r[row].getindexptr();
1284  B* aptr = r[row].getptr() + offset;
1285 
1286  return *aptr;
1287  }
1288 
1289  //determine whether overflow has to be taken into account or not
1290  if (r[row].getsize() == avg)
1291  return overflow[std::make_pair(row,col)];
1292  else
1293  {
1294  //modify index array
1295  *end = col;
1296 
1297  //do simulatenous operations on data array a
1298  std::ptrdiff_t offset = end - r[row].getindexptr();
1299  B* apos = r[row].getptr() + offset;
1300 
1301  //increase rowsize
1302  r[row].setsize(r[row].getsize()+1);
1303 
1304  //return reference to the newly created entry
1305  return *apos;
1306  }
1307  }
1308 
1310 
1321  {
1322  if (build_mode!=implicit)
1323  DUNE_THROW(BCRSMatrixError,"requires implicit build mode");
1324  if (ready==built)
1325  DUNE_THROW(BCRSMatrixError,"matrix already built up, no more need for compression");
1326  if (ready==notAllocated)
1327  DUNE_THROW(BCRSMatrixError,"matrix size not set and no memory allocated yet");
1328  if (ready!=building)
1329  DUNE_THROW(InvalidStateException,"You may only call compress() at the end of the 'building' stage");
1330 
1331  //calculate statistics
1332  CompressionStatistics stats;
1333  stats.overflow_total = overflow.size();
1334  stats.maximum = 0;
1335 
1336  //get insertion iterators pointing to one before start (for later use of ++it)
1337  size_type* jiit = j.get();
1338  B* aiit = a;
1339 
1340  //get iterator to the smallest overflow element
1341  typename OverflowType::iterator oit = overflow.begin();
1342 
1343  //store a copy of index pointers on which to perform sortation
1344  std::vector<size_type*> perm;
1345 
1346  //iterate over all rows and copy elements into their position in the compressed array
1347  for (size_type i=0; i<n; i++)
1348  {
1349  //get old pointers into a and j and size without overflow changes
1350  size_type* begin = r[i].getindexptr();
1351  //B* apos = r[i].getptr();
1352  size_type size = r[i].getsize();
1353 
1354  perm.resize(size);
1355 
1356  typename std::vector<size_type*>::iterator it = perm.begin();
1357  for (size_type* iit = begin; iit < begin + size; ++iit, ++it)
1358  *it = iit;
1359 
1360  //sort permutation array
1361  std::sort(perm.begin(),perm.end(),PointerCompare<size_type>());
1362 
1363  //change row window pointer to their new positions
1364  r[i].setindexptr(jiit);
1365  r[i].setptr(aiit);
1366 
1367  for (it = perm.begin(); it != perm.end(); ++it)
1368  {
1369  //check whether there are elements in the overflow area which take precedence
1370  while ((oit!=overflow.end()) && (oit->first < std::make_pair(i,**it)))
1371  {
1372  //check whether there is enough memory to write to
1373  if (jiit > begin)
1375  "Allocated memory for BCRSMatrix exhausted during compress()!"
1376  "Please increase either the average number of entries per row or the overflow fraction."
1377  );
1378  //copy and element from the overflow area to the insertion position in a and j
1379  *jiit = oit->first.second;
1380  ++jiit;
1381  *aiit = oit->second;
1382  ++aiit;
1383  ++oit;
1384  r[i].setsize(r[i].getsize()+1);
1385  }
1386 
1387  //check whether there is enough memory to write to
1388  if (jiit > begin)
1390  "Allocated memory for BCRSMatrix exhausted during compress()!"
1391  "Please increase either the average number of entries per row or the overflow fraction."
1392  );
1393 
1394  //copy element from array
1395  *jiit = **it;
1396  ++jiit;
1397  B* apos = *it-j.get()+a;
1398  *aiit = *apos;
1399  ++aiit;
1400  }
1401 
1402  //copy remaining elements from the overflow area
1403  while ((oit!=overflow.end()) && (oit->first.first == i))
1404  {
1405  //check whether there is enough memory to write to
1406  if (jiit > begin)
1408  "Allocated memory for BCRSMatrix exhausted during compress()!"
1409  "Please increase either the average number of entries per row or the overflow fraction."
1410  );
1411 
1412  //copy and element from the overflow area to the insertion position in a and j
1413  *jiit = oit->first.second;
1414  ++jiit;
1415  *aiit = oit->second;
1416  ++aiit;
1417  ++oit;
1418  r[i].setsize(r[i].getsize()+1);
1419  }
1420 
1421  // update maximum row size
1422  if (r[i].getsize()>stats.maximum)
1423  stats.maximum = r[i].getsize();
1424  }
1425 
1426  // overflow area may be cleared
1427  overflow.clear();
1428 
1429  //determine average number of entries and memory usage
1430  std::ptrdiff_t diff = (r[n-1].getindexptr() + r[n-1].getsize() - j.get());
1431  nnz = diff;
1432  stats.avg = (double) (nnz) / (double) n;
1433  stats.mem_ratio = (double) (nnz)/(double) allocationSize;
1434 
1435  //matrix is now built
1436  ready = built;
1437 
1438  return stats;
1439  }
1440 
1441  //===== vector space arithmetic
1442 
1445  {
1446 #ifdef DUNE_ISTL_WITH_CHECKING
1447  if (ready != built)
1448  DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1449 #endif
1450 
1451  if (nnz>0)
1452  {
1453  // process 1D array
1454  for (size_type i=0; i<nnz; i++)
1455  a[i] *= k;
1456  }
1457  else
1458  {
1459  RowIterator endi=end();
1460  for (RowIterator i=begin(); i!=endi; ++i)
1461  {
1462  ColIterator endj = (*i).end();
1463  for (ColIterator j=(*i).begin(); j!=endj; ++j)
1464  (*j) *= k;
1465  }
1466  }
1467 
1468  return *this;
1469  }
1470 
1473  {
1474 #ifdef DUNE_ISTL_WITH_CHECKING
1475  if (ready != built)
1476  DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1477 #endif
1478 
1479  if (nnz>0)
1480  {
1481  // process 1D array
1482  for (size_type i=0; i<nnz; i++)
1483  a[i] /= k;
1484  }
1485  else
1486  {
1487  RowIterator endi=end();
1488  for (RowIterator i=begin(); i!=endi; ++i)
1489  {
1490  ColIterator endj = (*i).end();
1491  for (ColIterator j=(*i).begin(); j!=endj; ++j)
1492  (*j) /= k;
1493  }
1494  }
1495 
1496  return *this;
1497  }
1498 
1499 
1506  {
1507 #ifdef DUNE_ISTL_WITH_CHECKING
1508  if (ready != built || b.ready != built)
1509  DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1510  if(N()!=b.N() || M() != b.M())
1511  DUNE_THROW(RangeError, "Matrix sizes do not match!");
1512 #endif
1513  RowIterator endi=end();
1514  ConstRowIterator j=b.begin();
1515  for (RowIterator i=begin(); i!=endi; ++i, ++j) {
1516  i->operator+=(*j);
1517  }
1518 
1519  return *this;
1520  }
1521 
1528  {
1529 #ifdef DUNE_ISTL_WITH_CHECKING
1530  if (ready != built || b.ready != built)
1531  DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1532  if(N()!=b.N() || M() != b.M())
1533  DUNE_THROW(RangeError, "Matrix sizes do not match!");
1534 #endif
1535  RowIterator endi=end();
1536  ConstRowIterator j=b.begin();
1537  for (RowIterator i=begin(); i!=endi; ++i, ++j) {
1538  i->operator-=(*j);
1539  }
1540 
1541  return *this;
1542  }
1543 
1553  {
1554 #ifdef DUNE_ISTL_WITH_CHECKING
1555  if (ready != built || b.ready != built)
1556  DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1557  if(N()!=b.N() || M() != b.M())
1558  DUNE_THROW(RangeError, "Matrix sizes do not match!");
1559 #endif
1560  RowIterator endi=end();
1561  ConstRowIterator j=b.begin();
1562  for(RowIterator i=begin(); i!=endi; ++i, ++j)
1563  i->axpy(alpha, *j);
1564 
1565  return *this;
1566  }
1567 
1568  //===== linear maps
1569 
1571  template<class X, class Y>
1572  void mv (const X& x, Y& y) const
1573  {
1574 #ifdef DUNE_ISTL_WITH_CHECKING
1575  if (ready != built)
1576  DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1577  if (x.N()!=M()) DUNE_THROW(BCRSMatrixError,
1578  "Size mismatch: M: " << N() << "x" << M() << " x: " << x.N());
1579  if (y.N()!=N()) DUNE_THROW(BCRSMatrixError,
1580  "Size mismatch: M: " << N() << "x" << M() << " y: " << y.N());
1581 #endif
1582  ConstRowIterator endi=end();
1583  for (ConstRowIterator i=begin(); i!=endi; ++i)
1584  {
1585  y[i.index()]=0;
1586  ConstColIterator endj = (*i).end();
1587  for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1588  (*j).umv(x[j.index()],y[i.index()]);
1589  }
1590  }
1591 
1593  template<class X, class Y>
1594  void umv (const X& x, Y& y) const
1595  {
1596 #ifdef DUNE_ISTL_WITH_CHECKING
1597  if (ready != built)
1598  DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1599  if (x.N()!=M()) DUNE_THROW(BCRSMatrixError,"index out of range");
1600  if (y.N()!=N()) DUNE_THROW(BCRSMatrixError,"index out of range");
1601 #endif
1602  ConstRowIterator endi=end();
1603  for (ConstRowIterator i=begin(); i!=endi; ++i)
1604  {
1605  ConstColIterator endj = (*i).end();
1606  for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1607  (*j).umv(x[j.index()],y[i.index()]);
1608  }
1609  }
1610 
1612  template<class X, class Y>
1613  void mmv (const X& x, Y& y) const
1614  {
1615 #ifdef DUNE_ISTL_WITH_CHECKING
1616  if (ready != built)
1617  DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1618  if (x.N()!=M()) DUNE_THROW(BCRSMatrixError,"index out of range");
1619  if (y.N()!=N()) DUNE_THROW(BCRSMatrixError,"index out of range");
1620 #endif
1621  ConstRowIterator endi=end();
1622  for (ConstRowIterator i=begin(); i!=endi; ++i)
1623  {
1624  ConstColIterator endj = (*i).end();
1625  for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1626  (*j).mmv(x[j.index()],y[i.index()]);
1627  }
1628  }
1629 
1631  template<class X, class Y>
1632  void usmv (const field_type& alpha, const X& x, Y& y) const
1633  {
1634 #ifdef DUNE_ISTL_WITH_CHECKING
1635  if (ready != built)
1636  DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1637  if (x.N()!=M()) DUNE_THROW(BCRSMatrixError,"index out of range");
1638  if (y.N()!=N()) DUNE_THROW(BCRSMatrixError,"index out of range");
1639 #endif
1640  ConstRowIterator endi=end();
1641  for (ConstRowIterator i=begin(); i!=endi; ++i)
1642  {
1643  ConstColIterator endj = (*i).end();
1644  for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1645  (*j).usmv(alpha,x[j.index()],y[i.index()]);
1646  }
1647  }
1648 
1650  template<class X, class Y>
1651  void mtv (const X& x, Y& y) const
1652  {
1653 #ifdef DUNE_ISTL_WITH_CHECKING
1654  if (ready != built)
1655  DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1656  if (x.N()!=N()) DUNE_THROW(BCRSMatrixError,"index out of range");
1657  if (y.N()!=M()) DUNE_THROW(BCRSMatrixError,"index out of range");
1658 #endif
1659  for(size_type i=0; i<y.N(); ++i)
1660  y[i]=0;
1661  umtv(x,y);
1662  }
1663 
1665  template<class X, class Y>
1666  void umtv (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  ConstRowIterator endi=end();
1675  for (ConstRowIterator i=begin(); i!=endi; ++i)
1676  {
1677  ConstColIterator endj = (*i).end();
1678  for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1679  (*j).umtv(x[i.index()],y[j.index()]);
1680  }
1681  }
1682 
1684  template<class X, class Y>
1685  void mmtv (const X& x, Y& y) const
1686  {
1687 #ifdef DUNE_ISTL_WITH_CHECKING
1688  if (x.N()!=N()) DUNE_THROW(BCRSMatrixError,"index out of range");
1689  if (y.N()!=M()) DUNE_THROW(BCRSMatrixError,"index out of range");
1690 #endif
1691  ConstRowIterator endi=end();
1692  for (ConstRowIterator i=begin(); i!=endi; ++i)
1693  {
1694  ConstColIterator endj = (*i).end();
1695  for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1696  (*j).mmtv(x[i.index()],y[j.index()]);
1697  }
1698  }
1699 
1701  template<class X, class Y>
1702  void usmtv (const field_type& alpha, const X& x, Y& y) const
1703  {
1704 #ifdef DUNE_ISTL_WITH_CHECKING
1705  if (ready != built)
1706  DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
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=begin(); i!=endi; ++i)
1712  {
1713  ConstColIterator endj = (*i).end();
1714  for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1715  (*j).usmtv(alpha,x[i.index()],y[j.index()]);
1716  }
1717  }
1718 
1720  template<class X, class Y>
1721  void umhv (const X& x, Y& y) const
1722  {
1723 #ifdef DUNE_ISTL_WITH_CHECKING
1724  if (ready != built)
1725  DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1726  if (x.N()!=N()) DUNE_THROW(BCRSMatrixError,"index out of range");
1727  if (y.N()!=M()) DUNE_THROW(BCRSMatrixError,"index out of range");
1728 #endif
1729  ConstRowIterator endi=end();
1730  for (ConstRowIterator i=begin(); i!=endi; ++i)
1731  {
1732  ConstColIterator endj = (*i).end();
1733  for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1734  (*j).umhv(x[i.index()],y[j.index()]);
1735  }
1736  }
1737 
1739  template<class X, class Y>
1740  void mmhv (const X& x, Y& y) const
1741  {
1742 #ifdef DUNE_ISTL_WITH_CHECKING
1743  if (ready != built)
1744  DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1745  if (x.N()!=N()) DUNE_THROW(BCRSMatrixError,"index out of range");
1746  if (y.N()!=M()) DUNE_THROW(BCRSMatrixError,"index out of range");
1747 #endif
1748  ConstRowIterator endi=end();
1749  for (ConstRowIterator i=begin(); i!=endi; ++i)
1750  {
1751  ConstColIterator endj = (*i).end();
1752  for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1753  (*j).mmhv(x[i.index()],y[j.index()]);
1754  }
1755  }
1756 
1758  template<class X, class Y>
1759  void usmhv (const field_type& alpha, const X& x, Y& y) const
1760  {
1761 #ifdef DUNE_ISTL_WITH_CHECKING
1762  if (ready != built)
1763  DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1764  if (x.N()!=N()) DUNE_THROW(BCRSMatrixError,"index out of range");
1765  if (y.N()!=M()) DUNE_THROW(BCRSMatrixError,"index out of range");
1766 #endif
1767  ConstRowIterator endi=end();
1768  for (ConstRowIterator i=begin(); i!=endi; ++i)
1769  {
1770  ConstColIterator endj = (*i).end();
1771  for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1772  (*j).usmhv(alpha,x[i.index()],y[j.index()]);
1773  }
1774  }
1775 
1776 
1777  //===== norms
1778 
1780  typename FieldTraits<field_type>::real_type frobenius_norm2 () const
1781  {
1782 #ifdef DUNE_ISTL_WITH_CHECKING
1783  if (ready != built)
1784  DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1785 #endif
1786 
1787  double sum=0;
1788 
1789  ConstRowIterator endi=end();
1790  for (ConstRowIterator i=begin(); i!=endi; ++i)
1791  {
1792  ConstColIterator endj = (*i).end();
1793  for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1794  sum += (*j).frobenius_norm2();
1795  }
1796 
1797  return sum;
1798  }
1799 
1801  typename FieldTraits<field_type>::real_type frobenius_norm () const
1802  {
1803  return sqrt(frobenius_norm2());
1804  }
1805 
1807  typename FieldTraits<field_type>::real_type infinity_norm () const
1808  {
1809  if (ready != built)
1810  DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1811 
1812  double max=0;
1813  ConstRowIterator endi=end();
1814  for (ConstRowIterator i=begin(); i!=endi; ++i)
1815  {
1816  double sum=0;
1817  ConstColIterator endj = (*i).end();
1818  for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1819  sum += (*j).infinity_norm();
1820  max = std::max(max,sum);
1821  }
1822  return max;
1823  }
1824 
1826  typename FieldTraits<field_type>::real_type infinity_norm_real () const
1827  {
1828 #ifdef DUNE_ISTL_WITH_CHECKING
1829  if (ready != built)
1830  DUNE_THROW(BCRSMatrixError,"You can only call arithmetic operations on fully built BCRSMatrix instances");
1831 #endif
1832 
1833  double max=0;
1834  ConstRowIterator endi=end();
1835  for (ConstRowIterator i=begin(); i!=endi; ++i)
1836  {
1837  double sum=0;
1838  ConstColIterator endj = (*i).end();
1839  for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1840  sum += (*j).infinity_norm_real();
1841  max = std::max(max,sum);
1842  }
1843  return max;
1844  }
1845 
1846 
1847  //===== sizes
1848 
1850  size_type N () const
1851  {
1852  return n;
1853  }
1854 
1856  size_type M () const
1857  {
1858  return m;
1859  }
1860 
1863  {
1864  return nnz;
1865  }
1866 
1869  {
1870  return ready;
1871  }
1872 
1875  {
1876  return build_mode;
1877  }
1878 
1879  //===== query
1880 
1882  bool exists (size_type i, size_type j) const
1883  {
1884 #ifdef DUNE_ISTL_WITH_CHECKING
1885  if (i<0 || i>=n) DUNE_THROW(BCRSMatrixError,"row index out of range");
1886  if (j<0 || j>=m) DUNE_THROW(BCRSMatrixError,"column index out of range");
1887 #endif
1888  if (r[i].size() && r[i].find(j)!=r[i].end())
1889  return true;
1890  else
1891  return false;
1892  }
1893 
1894 
1895  private:
1896  // state information
1897  BuildMode build_mode; // row wise or whole matrix
1898  BuildStage ready; // indicate the stage the matrix building is in
1899 
1900  // The allocator used for memory management
1901  typename A::template rebind<B>::other allocator_;
1902 
1903  typename A::template rebind<row_type>::other rowAllocator_;
1904 
1905  typename A::template rebind<size_type>::other sizeAllocator_;
1906 
1907  // size of the matrix
1908  size_type n; // number of rows
1909  size_type m; // number of columns
1910  size_type nnz; // number of nonzeroes contained in the matrix
1911  size_type allocationSize; //allocated size of a and j arrays, except in implicit mode: nnz==allocationsSize
1912  // zero means that memory is allocated separately for each row.
1913 
1914  // the rows are dynamically allocated
1915  row_type* r; // [n] the individual rows having pointers into a,j arrays
1916 
1917  // dynamically allocated memory
1918  B* a; // [allocationSize] non-zero entries of the matrix in row-wise ordering
1919  // If a single array of column indices is used, it can be shared
1920  // between different matrices with the same sparsity pattern
1921  Dune::shared_ptr<size_type> j; // [allocationSize] column indices of entries
1922 
1923  // additional data is needed in implicit buildmode
1924  size_type avg;
1925  double overflowsize;
1926 
1927  typedef std::map<std::pair<size_type,size_type>, B> OverflowType;
1928  OverflowType overflow;
1929 
1930  void setWindowPointers(ConstRowIterator row)
1931  {
1932  row_type current_row(a,j.get(),0); // Pointers to current row data
1933  for (size_type i=0; i<n; i++, ++row) {
1934  // set row i
1935  size_type s = row->getsize();
1936 
1937  if (s>0) {
1938  // setup pointers and size
1939  r[i].set(s,current_row.getptr(), current_row.getindexptr());
1940  // update pointer for next row
1941  current_row.setptr(current_row.getptr()+s);
1942  current_row.setindexptr(current_row.getindexptr()+s);
1943  } else{
1944  // empty row
1945  r[i].set(0,0,0);
1946  }
1947  }
1948  }
1949 
1951 
1955  void setColumnPointers(ConstRowIterator row)
1956  {
1957  size_type* jptr = j.get();
1958  for (size_type i=0; i<n; ++i, ++row) {
1959  // set row i
1960  size_type s = row->getsize();
1961 
1962  if (s>0) {
1963  // setup pointers and size
1964  r[i].setsize(s);
1965  r[i].setindexptr(jptr);
1966  } else{
1967  // empty row
1968  r[i].set(0,0,0);
1969  }
1970 
1971  // advance position in global array
1972  jptr += s;
1973  }
1974  }
1975 
1977 
1981  void setDataPointers()
1982  {
1983  B* aptr = a;
1984  for (size_type i=0; i<n; ++i) {
1985  // set row i
1986  if (r[i].getsize() > 0) {
1987  // setup pointers and size
1988  r[i].setptr(aptr);
1989  } else{
1990  // empty row
1991  r[i].set(0,0,0);
1992  }
1993 
1994  // advance position in global array
1995  aptr += r[i].getsize();
1996  }
1997  }
1998 
2000  void copyWindowStructure(const BCRSMatrix& Mat)
2001  {
2002  setWindowPointers(Mat.begin());
2003 
2004  // copy data
2005  for (size_type i=0; i<n; i++) r[i] = Mat.r[i];
2006 
2007  // finish off
2008  build_mode = row_wise; // dummy
2009  ready = built;
2010  }
2011 
2017  void deallocate(bool deallocateRows=true)
2018  {
2019 
2020  if (notAllocated)
2021  return;
2022 
2023  if (allocationSize>0)
2024  {
2025  // a,j have been allocated as one long vector
2026  j.reset();
2027  if (a)
2028  {
2029  for(B *aiter=a+(allocationSize-1), *aend=a-1; aiter!=aend; --aiter)
2030  allocator_.destroy(aiter);
2031  allocator_.deallocate(a,allocationSize);
2032  a = nullptr;
2033  }
2034  }
2035  else if (r)
2036  {
2037  // check if memory for rows have been allocated individually
2038  for (size_type i=0; i<n; i++)
2039  if (r[i].getsize()>0)
2040  {
2041  for (B *col=r[i].getptr()+(r[i].getsize()-1),
2042  *colend = r[i].getptr()-1; col!=colend; --col) {
2043  allocator_.destroy(col);
2044  }
2045  sizeAllocator_.deallocate(r[i].getindexptr(),1);
2046  allocator_.deallocate(r[i].getptr(),1);
2047  // clear out row data in case we don't want to deallocate the rows
2048  // otherwise we might run into a double free problem here later
2049  r[i].set(0,nullptr,nullptr);
2050  }
2051  }
2052 
2053  // deallocate the rows
2054  if (n>0 && deallocateRows && r) {
2055  for(row_type *riter=r+(n-1), *rend=r-1; riter!=rend; --riter)
2056  rowAllocator_.destroy(riter);
2057  rowAllocator_.deallocate(r,n);
2058  r = nullptr;
2059  }
2060 
2061  // Mark matrix as not built at all.
2062  ready=notAllocated;
2063 
2064  }
2065 
2067  class Deallocator
2068  {
2069  typename A::template rebind<size_type>::other& sizeAllocator_;
2070 
2071  public:
2072  Deallocator(typename A::template rebind<size_type>::other& sizeAllocator)
2073  : sizeAllocator_(sizeAllocator)
2074  {}
2075 
2076  void operator()(size_type* p) { sizeAllocator_.deallocate(p,1); }
2077  };
2078 
2079 
2097  void allocate(size_type rows, size_type columns, size_type allocationSize_, bool allocateRows, bool allocate_data)
2098  {
2099  // Store size
2100  n = rows;
2101  m = columns;
2102  nnz = allocationSize_;
2103  allocationSize = allocationSize_;
2104 
2105  // allocate rows
2106  if(allocateRows) {
2107  if (n>0) {
2108  if (r)
2109  DUNE_THROW(InvalidStateException,"Rows have already been allocated, cannot allocate a second time");
2110  r = rowAllocator_.allocate(rows);
2111  }else{
2112  r = 0;
2113  }
2114  }
2115 
2116  // allocate a and j array
2117  if (allocate_data)
2118  allocateData();
2119  if (allocationSize>0) {
2120  // allocate column indices only if not yet present (enable sharing)
2121  if (!j.get())
2122  j.reset(sizeAllocator_.allocate(allocationSize),Deallocator(sizeAllocator_));
2123  }else{
2124  j.reset();
2125  for(row_type* ri=r; ri!=r+rows; ++ri)
2126  rowAllocator_.construct(ri, row_type());
2127  }
2128 
2129  // Mark the matrix as not built.
2130  ready = building;
2131  }
2132 
2133  void allocateData()
2134  {
2135  if (a)
2136  DUNE_THROW(InvalidStateException,"Cannot allocate data array (already allocated)");
2137  if (allocationSize>0) {
2138  a = allocator_.allocate(allocationSize);
2139  // use placement new to call constructor that allocates
2140  // additional memory.
2141  new (a) B[allocationSize];
2142  } else {
2143  a = nullptr;
2144  }
2145  }
2146 
2152  void implicit_allocate(size_type _n, size_type _m)
2153  {
2154  if (build_mode != implicit)
2155  DUNE_THROW(InvalidStateException,"implicit_allocate() may only be called in implicit build mode");
2156  if (ready != notAllocated)
2157  DUNE_THROW(InvalidStateException,"memory has already been allocated");
2158 
2159  // check to make sure the user has actually set the parameters
2160  if (overflowsize < 0)
2161  DUNE_THROW(InvalidStateException,"You have to set the implicit build mode parameters before starting to build the matrix");
2162  //calculate size of overflow region, add buffer for row sort!
2163  size_type osize = (size_type) (_n*avg)*overflowsize + 4*avg;
2164  allocationSize = _n*avg + osize;
2165 
2166  allocate(_n, _m, allocationSize,true,true);
2167 
2168  //set row pointers correctly
2169  size_type* jptr = j.get() + osize;
2170  B* aptr = a + osize;
2171  for (size_type i=0; i<n; i++)
2172  {
2173  r[i].set(0,aptr,jptr);
2174  jptr = jptr + avg;
2175  aptr = aptr + avg;
2176  }
2177 
2178  ready = building;
2179  }
2180  };
2181 
2182 
2185 } // end namespace
2186 
2187 #endif