dune-istl  2.2.0
matrixmatrix.hh
Go to the documentation of this file.
1 #ifndef DUNE_MATRIXMATRIX_HH
2 #define DUNE_MATRIXMATRIX_HH
4 #include <dune/common/fmatrix.hh>
5 #include <dune/common/tuples.hh>
6 #include <dune/common/timer.hh>
7 namespace Dune
8 {
9 
20  namespace
21  {
22 
31  template<int b>
32  struct NonzeroPatternTraverser
33  {};
34 
35 
36  template<>
37  struct NonzeroPatternTraverser<0>
38  {
39  template<class T,class A1, class A2, class F, int n, int m, int k>
40  static void traverse(const Dune::BCRSMatrix<Dune::FieldMatrix<T,n,k>,A1>& A,
41  const Dune::BCRSMatrix<Dune::FieldMatrix<T,k,m>,A2>& B,
42  F& func)
43  {
44  typedef typename Dune::BCRSMatrix<Dune::FieldMatrix<T,k,m>,A2>::size_type size_type;
45 
46  if(A.M()!=B.N())
47  DUNE_THROW(ISTLError, "The sizes of the matrices do not match: "<<A.M()<<"!="<<B.N());
48 
49  typedef typename Dune::BCRSMatrix<Dune::FieldMatrix<T,n,k>,A1>::ConstRowIterator Row;
50  typedef typename Dune::BCRSMatrix<Dune::FieldMatrix<T,n,k>,A1>::ConstColIterator Col;
51  typedef typename Dune::BCRSMatrix<Dune::FieldMatrix<T,k,m>,A2>::ConstRowIterator BRow;
52  typedef typename Dune::BCRSMatrix<Dune::FieldMatrix<T,k,m>,A2>::ConstColIterator BCol;
53  for(Row row= A.begin(); row != A.end(); ++row){
54  // Loop over all column entries
55  for(Col col = row->begin(); col != row->end(); ++col){
56  // entry at i,k
57  // search for all nonzeros in row k
58  for(BCol bcol = B[col.index()].begin(); bcol != B[col.index()].end(); ++bcol){
59  func(*col, *bcol, row.index(), bcol.index());
60  }
61  }
62  }
63  }
64 
65  };
66 
67  template<>
68  struct NonzeroPatternTraverser<1>
69  {
70  template<class T, class A1, class A2, class F, int n, int m, int k>
71  static void traverse(const Dune::BCRSMatrix<Dune::FieldMatrix<T,k,n>,A1>& A,
72  const Dune::BCRSMatrix<Dune::FieldMatrix<T,k,m>,A2>& B,
73  F& func)
74  {
75 
76  if(A.N()!=B.N())
77  DUNE_THROW(ISTLError, "The sizes of the matrices do not match: "<<A.N()<<"!="<<B.N());
78 
79  typedef typename Dune::BCRSMatrix<Dune::FieldMatrix<T,k,n>,A1>::ConstRowIterator Row;
80  typedef typename Dune::BCRSMatrix<Dune::FieldMatrix<T,k,n>,A1>::ConstColIterator Col;
81  typedef typename Dune::BCRSMatrix<Dune::FieldMatrix<T,k,m>,A2>::ConstColIterator BCol;
82  typedef typename Dune::BCRSMatrix<Dune::FieldMatrix<T,k,n>,A1>::size_type size_t1;
83  typedef typename Dune::BCRSMatrix<Dune::FieldMatrix<T,k,m>,A2>::size_type size_t2;
84 
85  for(Row row=A.begin(); row!=A.end(); ++row){
86  for(Col col=row->begin(); col!=row->end(); ++col){
87  for(BCol bcol = B[row.index()].begin(); bcol != B[row.index()].end(); ++bcol){
88  func(*col, *bcol, col.index(), bcol.index());
89  }
90  }
91  }
92  }
93  };
94 
95  template<>
96  struct NonzeroPatternTraverser<2>
97  {
98  template<class T, class A1, class A2, class F, int n, int m, int k>
99  static void traverse(const BCRSMatrix<FieldMatrix<T,n,m>,A1>& mat,
100  const BCRSMatrix<FieldMatrix<T,k,m>,A2>& matt,
101  F& func)
102  {
103  if(mat.M()!=matt.M())
104  DUNE_THROW(ISTLError, "The sizes of the matrices do not match: "<<mat.N()<<"!="<<matt.N());
105 
106  typedef typename BCRSMatrix<FieldMatrix<T,n,m>,A1>::ConstRowIterator row_iterator;
107  typedef typename BCRSMatrix<FieldMatrix<T,n,m>,A1>::ConstColIterator col_iterator;
108  typedef typename BCRSMatrix<FieldMatrix<T,k,m>,A2>::ConstRowIterator row_iterator_t;
109  typedef typename BCRSMatrix<FieldMatrix<T,k,m>,A2>::ConstColIterator col_iterator_t;
110 
111  for(row_iterator mrow=mat.begin(); mrow != mat.end(); ++mrow){
112  //iterate over the column entries
113  // mt is a transposed matrix crs therefore it is treated as a ccs matrix
114  // and the row_iterator iterates over the columns of the transposed matrix.
115  // search the row of the transposed matrix for an entry with the same index
116  // as the mcol iterator
117 
118  for(row_iterator_t mtcol=matt.begin(); mtcol != matt.end(); ++mtcol, func.nextCol()){
119  //Search for col entries in mat that have a corrsponding row index in matt
120  // (i.e. corresponding col index in the as this is the transposed matrix
121  col_iterator_t mtrow=mtcol->begin();
122 
123  for(col_iterator mcol=mrow->begin(); mcol != mrow->end(); ++mcol){
124  // search
125  // TODO: This should probably be substituted by a binary search
126  for( ;mtrow != mtcol->end(); ++mtrow)
127  if(mtrow.index()>=mcol.index())
128  break;
129  if(mtrow != mtcol->end() && mtrow.index()==mcol.index()){
130  func(*mcol, *mtrow, mtcol.index());
131  // In some cases we only search for one pair, then we break here
132  // and continue with the next column.
133  if(F::do_break)
134  break;
135  }
136  }
137  }
138  func.nextRow();
139  }
140  }
141  };
142 
143 
144 
145  template<class T, class A, int n, int m>
146  class SparsityPatternInitializer
147  {
148  public:
149  enum{do_break=true};
150  typedef typename BCRSMatrix<FieldMatrix<T,n,m>,A>::CreateIterator CreateIterator;
151  typedef typename BCRSMatrix<FieldMatrix<T,n,m>,A>::size_type size_type;
152 
153  SparsityPatternInitializer(CreateIterator iter)
154  : rowiter(iter)
155  {}
156 
157  template<class T1, class T2>
158  void operator()(const T1& t1, const T2& t2, size_type j)
159  {
160  rowiter.insert(j);
161  }
162 
163  void nextRow()
164  {
165  ++rowiter;
166  }
167  void nextCol()
168  {}
169 
170  private:
171  CreateIterator rowiter;
172  };
173 
174 
175  template<int transpose, class T, class TA, int n, int m>
176  class MatrixInitializer
177  {
178  public:
179  enum{do_break=true};
182  typedef typename Matrix::size_type size_type;
183 
184  MatrixInitializer(Matrix& A_, size_type rows)
185  : count(0), A(A_)
186  {}
187  template<class T1, class T2>
188  void operator()(const T1& t1, const T2& t2, int j)
189  {
190  ++count;
191  }
192 
193  void nextCol()
194  {}
195 
196  void nextRow()
197  {}
198 
199  std::size_t nonzeros()
200  {
201  return count;
202  }
203 
204  template<class A1, class A2, int n2, int m2, int n3, int m3>
205  void initPattern(const BCRSMatrix<FieldMatrix<T,n2,m2>,A1>& mat1,
206  const BCRSMatrix<FieldMatrix<T,n3,m3>,A2>& mat2)
207  {
208  SparsityPatternInitializer<T, TA, n, m> sparsity(A.createbegin());
209  NonzeroPatternTraverser<transpose>::traverse(mat1,mat2,sparsity);
210  }
211 
212  private:
213  std::size_t count;
214  Matrix& A;
215  };
216 
217  template<class T, class TA, int n, int m>
218  class MatrixInitializer<1,T,TA,n,m>
219  {
220  public:
221  enum{do_break=false};
224  typedef typename Matrix::size_type size_type;
225 
226  MatrixInitializer(Matrix& A_, size_type rows)
227  : A(A_), entries(rows)
228  {}
229 
230  template<class T1, class T2>
231  void operator()(const T1& t1, const T2& t2, size_type i, size_type j)
232  {
233  entries[i].insert(j);
234  }
235 
236  void nextCol()
237  {}
238 
239  size_type nonzeros()
240  {
241  size_type nnz=0;
242  typedef typename std::vector<std::set<size_t> >::const_iterator Iter;
243  for(Iter iter = entries.begin(); iter != entries.end(); ++iter)
244  nnz+=(*iter).size();
245  return nnz;
246  }
247  template<class A1, class A2, int n2, int m2, int n3, int m3>
248  void initPattern(const BCRSMatrix<FieldMatrix<T,n2,m2>,A1>& mat1,
249  const BCRSMatrix<FieldMatrix<T,n3,m3>,A2>& mat2)
250  {
251  typedef typename std::vector<std::set<size_t> >::const_iterator Iter;
252  CreateIterator citer = A.createbegin();
253  for(Iter iter = entries.begin(); iter != entries.end(); ++iter, ++citer){
254  typedef std::set<size_t>::const_iterator SetIter;
255  for(SetIter index=iter->begin(); index != iter->end(); ++index)
256  citer.insert(*index);
257  }
258  }
259 
260  private:
261  Matrix& A;
262  std::vector<std::set<size_t> > entries;
263  };
264 
265  template<class T, class TA, int n, int m>
266  struct MatrixInitializer<0,T,TA,n,m>
267  : public MatrixInitializer<1,T,TA,n,m>
268  {
269  MatrixInitializer(Dune::BCRSMatrix<Dune::FieldMatrix<T,n,m>,TA>& A_,
270  typename Dune::BCRSMatrix<Dune::FieldMatrix<T,n,m>,TA>::size_type rows)
271  :MatrixInitializer<1,T,TA,n,m>(A_,rows)
272  {}
273  };
274 
275 
276  template<class T, class T1, class T2, int n, int m, int k>
277  void addMatMultTransposeMat(FieldMatrix<T,n,k>& res, const FieldMatrix<T1,n,m>& mat,
278  const FieldMatrix<T2,k,m>& matt)
279  {
280  typedef typename FieldMatrix<T,n,k>::size_type size_type;
281 
282  for(size_type row=0; row<n; ++row)
283  for(size_type col=0; col<k;++col){
284  for(size_type i=0; i < m; ++i)
285  res[row][col]+=mat[row][i]*matt[col][i];
286  }
287  }
288 
289  template<class T, class T1, class T2, int n, int m, int k>
290  void addTransposeMatMultMat(FieldMatrix<T,n,k>& res, const FieldMatrix<T1,m,n>& mat,
291  const FieldMatrix<T2,m,k>& matt)
292  {
293  typedef typename FieldMatrix<T,n,k>::size_type size_type;
294  for(size_type i=0; i<m; ++i)
295  for(size_type row=0; row<n;++row){
296  for(size_type col=0; col < k; ++col)
297  res[row][col]+=mat[i][row]*matt[i][col];
298  }
299  }
300 
301  template<class T, class T1, class T2, int n, int m, int k>
302  void addMatMultMat(FieldMatrix<T,n,m>& res, const FieldMatrix<T1,n,k>& mat,
303  const FieldMatrix<T2,k,m>& matt)
304  {
305  typedef typename FieldMatrix<T,n,k>::size_type size_type;
306  for(size_type row=0; row<n; ++row)
307  for(size_type col=0; col<m;++col){
308  for(size_type i=0; i < k; ++i)
309  res[row][col]+=mat[row][i]*matt[i][col];
310  }
311  }
312 
313 
314  template<class T, class A, int n, int m>
315  class EntryAccumulatorFather
316  {
317  public:
318  enum{do_break=false};
319  typedef BCRSMatrix<FieldMatrix<T,n,m>,A> Matrix;
320  typedef typename Matrix::RowIterator Row;
321  typedef typename Matrix::ColIterator Col;
322 
323  EntryAccumulatorFather(Matrix& mat_)
324  :mat(mat_), row(mat.begin())
325  {
326  mat=0;
327  col=row->begin();
328  }
329  void nextRow()
330  {
331  ++row;
332  if(row!=mat.end())
333  col=row->begin();
334  }
335 
336  void nextCol()
337  {
338  ++this->col;
339  }
340  protected:
341  Matrix& mat;
342  private:
343  Row row;
344  protected:
345  Col col;
346  };
347 
348  template<class T, class A, int n, int m, int transpose>
349  class EntryAccumulator
350  : public EntryAccumulatorFather<T,A,n,m>
351  {
352  public:
353  typedef BCRSMatrix<FieldMatrix<T,n,m>,A> Matrix;
354  typedef typename Matrix::size_type size_type;
355 
356  EntryAccumulator(Matrix& mat_)
357  : EntryAccumulatorFather<T,A,n,m>(mat_)
358  {}
359 
360  template<class T1, class T2>
361  void operator()(const T1& t1, const T2& t2, size_type i)
362  {
363  assert(this->col.index()==i);
364  addMatMultMat(*(this->col),t1,t2);
365  }
366  };
367 
368  template<class T, class A, int n, int m>
369  class EntryAccumulator<T,A,n,m,0>
370  : public EntryAccumulatorFather<T,A,n,m>
371  {
372  public:
373  typedef BCRSMatrix<FieldMatrix<T,n,m>,A> Matrix;
374  typedef typename Matrix::size_type size_type;
375 
376  EntryAccumulator(Matrix& mat_)
377  : EntryAccumulatorFather<T,A,n,m>(mat_)
378  {}
379 
380  template<class T1, class T2>
381  void operator()(const T1& t1, const T2& t2, size_type i, size_type j)
382  {
383  addMatMultMat(this->mat[i][j], t1, t2);
384  }
385  };
386 
387  template<class T, class A, int n, int m>
388  class EntryAccumulator<T,A,n,m,1>
389  : public EntryAccumulatorFather<T,A,n,m>
390  {
391  public:
392  typedef BCRSMatrix<FieldMatrix<T,n,m>,A> Matrix;
393  typedef typename Matrix::size_type size_type;
394 
395  EntryAccumulator(Matrix& mat_)
396  : EntryAccumulatorFather<T,A,n,m>(mat_)
397  {}
398 
399  template<class T1, class T2>
400  void operator()(const T1& t1, const T2& t2, size_type i, size_type j)
401  {
402  addTransposeMatMultMat(this->mat[i][j], t1, t2);
403  }
404  };
405 
406  template<class T, class A, int n, int m>
407  class EntryAccumulator<T,A,n,m,2>
408  : public EntryAccumulatorFather<T,A,n,m>
409  {
410  public:
411  typedef BCRSMatrix<FieldMatrix<T,n,m>,A> Matrix;
412  typedef typename Matrix::size_type size_type;
413 
414  EntryAccumulator(Matrix& mat_)
415  : EntryAccumulatorFather<T,A,n,m>(mat_)
416  {}
417 
418  template<class T1, class T2>
419  void operator()(const T1& t1, const T2& t2, size_type i)
420  {
421  assert(this->col.index()==i);
422  addMatMultTransposeMat(*this->col,t1,t2);
423  }
424  };
425 
426 
427  template<int transpose>
428  struct SizeSelector
429  {
430  };
431 
432  template<>
433  struct SizeSelector<0>
434  {
435  template<class M1, class M2>
436  static tuple<typename M1::size_type, typename M2::size_type>
437  size(const M1& m1, const M2& m2)
438  {
439  return make_tuple(m1.N(), m2.M());
440  }
441  };
442 
443  template<>
444  struct SizeSelector<1>
445  {
446  template<class M1, class M2>
447  static tuple<typename M1::size_type, typename M2::size_type>
448  size(const M1& m1, const M2& m2)
449  {
450  return make_tuple(m1.M(), m2.M());
451  }
452  };
453 
454 
455  template<>
456  struct SizeSelector<2>
457  {
458  template<class M1, class M2>
459  static tuple<typename M1::size_type, typename M2::size_type>
460  size(const M1& m1, const M2& m2)
461  {
462  return make_tuple(m1.N(), m2.N());
463  }
464  };
465 
466  template<int transpose, class T, class A, class A1, class A2, int n1, int m1, int n2, int m2, int n3, int m3>
467  void matMultMat(BCRSMatrix<FieldMatrix<T,n1,m1>,A>& res, const BCRSMatrix<FieldMatrix<T,n2,m2>,A1>& mat1,
468  const BCRSMatrix<FieldMatrix<T,n3,m3>,A2>& mat2)
469  {
470  // First step is to count the number of nonzeros
471  typename BCRSMatrix<FieldMatrix<T,n1,m1>,A>::size_type rows, cols;
472  tie(rows,cols)=SizeSelector<transpose>::size(mat1, mat2);
473  MatrixInitializer<transpose,T,A,n1,m1> patternInit(res, rows);
474  Timer timer;
475  NonzeroPatternTraverser<transpose>::traverse(mat1,mat2,patternInit);
476  res.setSize(rows, cols, patternInit.nonzeros());
477  res.setBuildMode(BCRSMatrix<FieldMatrix<T,n1,m1>,A>::row_wise);
478 
479  std::cout<<"Counting nonzeros took "<<timer.elapsed()<<std::endl;
480  timer.reset();
481 
482  // Second step is to allocate the storage for the result and initialize the nonzero pattern
483  patternInit.initPattern(mat1, mat2);
484 
485  std::cout<<"Setting up sparsity pattern took "<<timer.elapsed()<<std::endl;
486  timer.reset();
487  // As a last step calculate the entries
488  EntryAccumulator<T,A,n1,m1, transpose> entriesAccu(res);
489  NonzeroPatternTraverser<transpose>::traverse(mat1,mat2,entriesAccu);
490  std::cout<<"Calculating entries took "<<timer.elapsed()<<std::endl;
491  }
492 
493  }
494 
502  template<typename M1, typename M2>
504  {
505  };
506 
507  template<typename T, int n, int k, int m>
508  struct MatMultMatResult<FieldMatrix<T,n,k>,FieldMatrix<T,k,m> >
509  {
510  typedef FieldMatrix<T,n,m> type;
511  };
512 
513  template<typename T, typename A, typename A1, int n, int k, int m>
514  struct MatMultMatResult<BCRSMatrix<FieldMatrix<T,n,k>,A >,BCRSMatrix<FieldMatrix<T,k,m>,A1 > >
515  {
517  FieldMatrix<T,k,m> >::type,A> type;
518  };
519 
528  template<class T, class A, class A1, class A2, int n, int m, int k>
529  void matMultTransposeMat(BCRSMatrix<FieldMatrix<T,n,k>,A>& res, const BCRSMatrix<FieldMatrix<T,n,m>,A1>& mat,
530  const BCRSMatrix<FieldMatrix<T,k,m>,A2>& matt, bool tryHard=false)
531  {
532  matMultMat<2>(res,mat, matt);
533  }
534 
543  template<class T, class A, class A1, class A2, int n, int m, int k>
544  void matMultMat(BCRSMatrix<FieldMatrix<T,n,m>,A>& res, const BCRSMatrix<FieldMatrix<T,n,k>,A1>& mat,
545  const BCRSMatrix<FieldMatrix<T,k,m>,A2>& matt, bool tryHard=false)
546  {
547  matMultMat<0>(res,mat, matt);
548  }
549 
558  template<class T, class A, class A1, class A2, int n, int m, int k>
559  void transposeMatMultMat(BCRSMatrix<FieldMatrix<T,n,m>,A>& res, const BCRSMatrix<FieldMatrix<T,k,n>,A1>& mat,
560  const BCRSMatrix<FieldMatrix<T,k,m>,A2>& matt, bool tryHard=false)
561  {
562  matMultMat<1>(res,mat, matt);
563  }
564 
565 }
566 #endif
567