dune-fem  2.4.1-rc
istlmatrix.hh
Go to the documentation of this file.
1 #ifndef DUNE_FEM_ISTLMATRIXWRAPPER_HH
2 #define DUNE_FEM_ISTLMATRIXWRAPPER_HH
3 
4 #if HAVE_DUNE_ISTL
5 
6 //- system includes
7 #include <vector>
8 #include <iostream>
9 #include <fstream>
10 #include <algorithm>
11 #include <set>
12 #include <map>
13 #include <string>
14 
15 //- Dune common includes
16 #include <dune/common/exceptions.hh>
17 
18 //- Dune istl includes
19 #include <dune/istl/bvector.hh>
20 #include <dune/istl/bcrsmatrix.hh>
21 #include <dune/istl/preconditioners.hh>
22 
23 //- Dune fem includes
31 #include <dune/fem/io/parameter.hh>
34 
37 
38 namespace Dune
39 {
40 
41  namespace Fem
42  {
43 
44  struct ISTLMatrixParameter
45  : public MatrixParameter
46  {
47 
48  ISTLMatrixParameter( const std::string keyPrefix = "istl." )
49  : keyPrefix_( keyPrefix )
50  {}
51 
52  virtual double overflowFraction () const
53  {
54  return Parameter::getValue< double >( keyPrefix_ + "matrix.overflowfraction", 1.0 );
55  }
56 
57  virtual int numIterations () const
58  {
59  return Parameter::getValue< int >( keyPrefix_ + "preconditioning.iterations", 5 );
60  }
61 
62  virtual double relaxation () const
63  {
64  return Parameter::getValue< double >( keyPrefix_ + "preconditioning.relaxation", 1.1 );
65  }
66 
67  virtual int method () const
68  {
69  static const std::string preConTable[]
70  = { "none", "ssor", "sor", "ilu-0", "ilu-n", "gauss-seidel", "jacobi", "amg-ilu-0", "amg-ilu-n", "amg-jacobi" };
71  return Parameter::getEnum( keyPrefix_ + "preconditioning.method", preConTable, 0 );
72  }
73 
74  virtual std::string preconditionName() const
75  {
76  static const std::string preConTable[]
77  = { "None", "SSOR", "SOR", "ILU-0", "ILU-n", "Gauss-Seidel", "Jacobi", "AMG-ILU-0", "AMG-ILU-n", "AMG-Jacobi" };
78  const int precond = method();
79  std::stringstream tmp;
80  tmp << preConTable[precond];
81 
82  if( precond != 3 )
83  tmp << " n=" << numIterations();
84  tmp << " relax=" << relaxation();
85  return tmp.str();
86  }
87 
88  private:
89  std::string keyPrefix_;
90 
91  };
92 
94  // --ISTLMatrixHandle
96  template <class LittleBlockType, class RowDiscreteFunctionImp, class ColDiscreteFunctionImp = RowDiscreteFunctionImp>
97  class ImprovedBCRSMatrix : public BCRSMatrix<LittleBlockType>
98  {
99  public:
100  typedef RowDiscreteFunctionImp RowDiscreteFunctionType;
101  typedef ColDiscreteFunctionImp ColDiscreteFunctionType;
102 
103  typedef BCRSMatrix<LittleBlockType> BaseType;
105  typedef BaseType MatrixBaseType;
106  typedef typename BaseType :: RowIterator RowIteratorType ;
107  typedef typename BaseType :: ColIterator ColIteratorType ;
108 
109  typedef ImprovedBCRSMatrix< LittleBlockType, RowDiscreteFunctionImp, ColDiscreteFunctionImp > ThisType;
110 
111  typedef typename BaseType :: size_type size_type;
112 
113  //===== type definitions and constants
114 
116  typedef typename BaseType::field_type field_type;
117 
119  typedef typename BaseType::block_type block_type;
120 
122  typedef typename BaseType:: allocator_type allocator_type;
123 
125  typedef typename BaseType :: row_type row_type;
126 
128  enum {
130  blocklevel = BaseType :: blocklevel
131  };
132 
134  typedef typename BaseType :: ColIterator ColIterator;
135 
137  typedef typename BaseType :: ConstColIterator ConstColIterator;
138 
140  typedef typename BaseType :: RowIterator RowIterator;
141 
143  typedef typename BaseType :: ConstRowIterator ConstRowIterator;
144 
146  typedef typename ColDiscreteFunctionType :: DiscreteFunctionSpaceType RangeSpaceType;
147 
149  typedef typename RowDiscreteFunctionType :: DofStorageType RowBlockVectorType;
150 
152  typedef typename ColDiscreteFunctionType :: DofStorageType ColBlockVectorType;
153 
155  typedef typename RangeSpaceType :: GridType :: Traits :: CollectiveCommunication CollectiveCommunictionType ;
156 
157  typedef typename BaseType :: BuildMode BuildMode ;
158 
159  public:
161  ImprovedBCRSMatrix(size_type rows, size_type cols, size_type nnz, double overflowFraction) :
162  BaseType (rows, cols, nnz, overflowFraction, BaseType::implicit)
163  {}
164 
166  ImprovedBCRSMatrix(size_type rows, size_type cols, size_type nz = 0 ) :
167  BaseType (rows, cols, BaseType :: row_wise)
168  {}
169 
171  ImprovedBCRSMatrix( ) :
172  BaseType ()
173  {}
174 
176  ImprovedBCRSMatrix(const ImprovedBCRSMatrix& org) :
177  BaseType(org)
178  {}
179 
180  template <class RowKeyType, class ColKeyType>
181  void createEntries(const std::map<RowKeyType , std::set<ColKeyType> >& indices)
182  {
183  // not insert map of indices into matrix
184  auto endcreate = this->createend();
185  for(auto create = this->createbegin(); create != endcreate; ++create)
186  {
187  const auto it = indices.find( create.index() );
188  if (it == indices.end() )
189  continue;
190  const auto& localIndices = it->second;
191  const auto end = localIndices.end();
192  // insert all indices for this row
193  for (auto it = localIndices.begin(); it != end; ++it)
194  create.insert( *it );
195  }
196  }
197 
199  void clear()
200  {
201  for (auto& row : *this)
202  for (auto& entry : row)
203  entry = 0;
204  }
205 
207  template <class HangingNodesType>
208  void setup(ThisType& oldMatrix, const HangingNodesType& hangingNodes)
209  {
210  // necessary because element traversal not necessaryly is in ascending order
211  typedef std::set< std::pair<int, block_type> > LocalEntryType;
212  typedef std::map< int , LocalEntryType > EntriesType;
213  EntriesType entries;
214 
215  // map of indices
216  std::map< int , std::set<int> > indices;
217  // not insert map of indices into matrix
218  auto rowend = oldMatrix.end();
219  for(auto it = oldMatrix.begin(); it != rowend; ++it)
220  {
221  const auto row = it.index();
222  auto& localIndices = indices[ row ];
223 
224  if( hangingNodes.isHangingNode( row ) )
225  {
226  // insert columns into other columns
227  const auto& cols = hangingNodes.associatedDofs( row );
228  const auto colSize = cols.size();
229  for(auto i=0; i<colSize; ++i)
230  {
231  assert( ! hangingNodes.isHangingNode( cols[i].first ) );
232 
233  // get local indices of col
234  auto& localColIndices = indices[ cols[i].first ];
235  auto& localEntry = entries[ cols[i].first ];
236 
237  // copy from old matrix
238  auto endj = (*it).end();
239  for (auto j= (*it).begin(); j!=endj; ++j)
240  {
241  localColIndices.insert( j.index () );
242  localEntry.insert( std::make_pair( j.index(), (cols[i].second * (*j)) ));
243  }
244  }
245 
246  // insert diagonal and hanging columns
247  localIndices.insert( row );
248  for(auto i=0; i<colSize; ++i)
249  localIndices.insert( cols[i].first );
250  }
251  else
252  {
253  // copy from old matrix
254  auto endj = (*it).end();
255  for (auto j= (*it).begin(); j!=endj; ++j)
256  localIndices.insert( j.index () );
257  }
258  }
259 
260  // create matrix from entry map
261  createEntries( indices );
262 
263  // not insert map of indices into matrix
264  auto rowit = oldMatrix.begin();
265 
266  auto endcreate = this->end();
267  for(auto create = this->begin(); create != endcreate; ++create, ++rowit )
268  {
269  assert( rowit != oldMatrix.end() );
270 
271  const auto row = create.index();
272  if( hangingNodes.isHangingNode( row ) )
273  {
274  const auto& cols = hangingNodes.associatedDofs( row );
275 
276  std::map< const int , block_type > colMap;
277  // only working for block size 1 ath the moment
278  assert( block_type :: rows == 1 );
279  // insert columns into map
280  const auto colSize = cols.size();
281  for( auto i=0; i<colSize; ++i)
282  colMap[ cols[i].first ] = -cols[i].second;
283  // insert diagonal into map
284  colMap[ row ] = 1;
285 
286  auto endj = (*create).end();
287  for (auto j= (*create).begin(); j!=endj; ++j)
288  {
289  assert( colMap.find( j.index() ) != colMap.end() );
290  (*j) = colMap[ j.index() ];
291  }
292  }
293  // if entries are equal, just copy
294  else if ( entries.find( row ) == entries.end() )
295  {
296  auto colit = (*rowit).begin();
297  auto endj = (*create).end();
298  for (auto j= (*create).begin(); j!=endj; ++j, ++colit )
299  {
300  assert( colit != (*rowit).end() );
301  (*j) = (*colit);
302  }
303  }
304  else
305  {
306  std::map< int , block_type > oldCols;
307 
308  {
309  auto colend = (*rowit).end();
310  for(auto colit = (*rowit).begin(); colit != colend; ++colit)
311  oldCols[ colit.index() ] = 0;
312  }
313 
314  auto entry = entries.find( row );
315  assert( entry != entries.end ());
316 
317  {
318  auto endcol = (*entry).second.end();
319  for( auto co = (*entry).second.begin(); co != endcol; ++co)
320  oldCols[ (*co).first ] = 0;
321  }
322 
323  {
324  auto colend = (*rowit).end();
325  for(auto colit = (*rowit).begin(); colit != colend; ++colit)
326  oldCols[ colit.index() ] += (*colit);
327  }
328 
329  {
330  auto endcol = (*entry).second.end();
331  for( auto co = (*entry).second.begin(); co != endcol; ++co)
332  oldCols[ (*co).first ] += (*co).second;
333  }
334 
335  auto endj = (*create).end();
336  for (auto j= (*create).begin(); j!=endj; ++j )
337  {
338  auto colEntry = oldCols.find( j.index() );
339  if( colEntry != oldCols.end() )
340  (*j) = (*colEntry).second;
341  else
342  abort();
343  }
344  }
345  }
346  }
347 
349  void extractDiagonal( ColBlockVectorType& diag ) const
350  {
351  const auto endi = this->end();
352  for (auto i = this->begin(); i!=endi; ++i)
353  {
354  // get diagonal entry of matrix
355  const auto row = i.index();
356  auto entry = (*i).find( row );
357  const LittleBlockType& block = (*entry);
358  enum { blockSize = LittleBlockType :: rows };
359  for( auto l=0; l<blockSize; ++l )
360  diag[ row ][ l ] = block[ l ][ l ];
361  }
362  }
363 
366  field_type operator()(const std::size_t row, const std::size_t col) const
367  {
368  const std::size_t blockRow(row/(LittleBlockType :: rows));
369  const std::size_t localRowIdx(row%(LittleBlockType :: rows));
370  const std::size_t blockCol(col/(LittleBlockType :: cols));
371  const std::size_t localColIdx(col%(LittleBlockType :: cols));
372 
373  const auto& matrixRow(this->operator[](blockRow));
374  auto entry = matrixRow.find( blockCol );
375  const LittleBlockType& block = (*entry);
376  return block[localRowIdx][localColIdx];
377  }
378 
381  void set(const std::size_t row, const std::size_t col, field_type value)
382  {
383  const std::size_t blockRow(row/(LittleBlockType :: rows));
384  const std::size_t localRowIdx(row%(LittleBlockType :: rows));
385  const std::size_t blockCol(col/(LittleBlockType :: cols));
386  const std::size_t localColIdx(col%(LittleBlockType :: cols));
387 
388  auto& matrixRow(this->operator[](blockRow));
389  auto entry = matrixRow.find( blockCol );
390  LittleBlockType& block = (*entry);
391  block[localRowIdx][localColIdx] = value;
392  }
393 
395  void print(std::ostream& s, unsigned int offset=0) const
396  {
397  s.precision( 6 );
398  const auto endi=this->end();
399  for (auto i=this->begin(); i!=endi; ++i)
400  {
401  const auto endj = (*i).end();
402  for (auto j=(*i).begin(); j!=endj; ++j)
403  if( (*j).infinity_norm() > 1.e-15)
404  s << i.index()+offset << " " << j.index()+offset << " " << *j << std::endl;
405  }
406  }
407  };
408 
409  template <class RowSpaceImp, class ColSpaceImp>
410  class ISTLMatrixObject;
411 
412  template <class RowSpaceImp, class ColSpaceImp = RowSpaceImp>
413  struct ISTLMatrixTraits
414  {
415  typedef RowSpaceImp RangeSpaceType;
416  typedef ColSpaceImp DomainSpaceType;
417  typedef ISTLMatrixTraits<DomainSpaceType,RangeSpaceType> ThisType;
418 
419  typedef ISTLMatrixObject<DomainSpaceType,RangeSpaceType> MatrixObjectType;
420  };
421 
423  template <class DomainSpaceImp, class RangeSpaceImp>
424  class ISTLMatrixObject
425  {
426  public:
428  typedef DomainSpaceImp DomainSpaceType;
430  typedef RangeSpaceImp RangeSpaceType;
431 
433  typedef ISTLMatrixObject<DomainSpaceType,RangeSpaceType> ThisType;
434 
435  typedef typename DomainSpaceType::GridType GridType;
436 
437  typedef typename RangeSpaceType :: EntityType RangeEntityType ;
438  typedef typename DomainSpaceType :: EntityType DomainEntityType ;
439 
440  enum { littleCols = DomainSpaceType :: localBlockSize };
441  enum { littleRows = RangeSpaceType :: localBlockSize };
442 
443  typedef FieldMatrix<typename DomainSpaceType :: RangeFieldType, littleRows, littleCols> LittleBlockType;
444 
445  typedef ISTLBlockVectorDiscreteFunction< RangeSpaceType > RowDiscreteFunctionType;
446  typedef ISTLBlockVectorDiscreteFunction< DomainSpaceType > ColumnDiscreteFunctionType;
447 
448  protected:
449  typedef typename RowDiscreteFunctionType :: DofStorageType RowBlockVectorType;
450  typedef typename ColumnDiscreteFunctionType :: DofStorageType ColumnBlockVectorType;
451 
452  typedef typename RangeSpaceType :: BlockMapperType RowMapperType;
453  typedef typename DomainSpaceType :: BlockMapperType ColMapperType;
454 
455  public:
457  typedef ImprovedBCRSMatrix< LittleBlockType , ColumnDiscreteFunctionType , RowDiscreteFunctionType > MatrixType;
458  typedef typename ISTLParallelMatrixAdapter<MatrixType,RangeSpaceType>::Type MatrixAdapterType;
459  typedef typename MatrixAdapterType :: ParallelScalarProductType ParallelScalarProductType;
460 
461  template <class MatrixObjectImp>
462  class LocalMatrix;
463 
464  struct LocalMatrixTraits
465  {
466  typedef DomainSpaceImp DomainSpaceType;
467  typedef RangeSpaceImp RangeSpaceType;
468  typedef typename DomainSpaceType :: RangeFieldType RangeFieldType;
469  typedef LocalMatrix<ThisType> LocalMatrixType;
470  typedef typename MatrixType:: block_type LittleBlockType;
471  };
472 
474  template <class MatrixObjectImp>
475  class LocalMatrix : public LocalMatrixDefault<LocalMatrixTraits>
476  {
477  public:
479  typedef LocalMatrixDefault<LocalMatrixTraits> BaseType;
480 
482  typedef MatrixObjectImp MatrixObjectType;
484  typedef typename MatrixObjectImp :: MatrixType MatrixType;
486  typedef typename MatrixType::block_type LittleBlockType;
487 
488  typedef typename MatrixObjectType::DomainSpaceType DomainSpaceType;
489  typedef typename MatrixObjectType::RangeSpaceType RangeSpaceType;
490 
491  typedef typename MatrixObjectType::DomainEntityType DomainEntityType;
492  typedef typename MatrixObjectType::RangeEntityType RangeEntityType;
493 
495  typedef typename DomainSpaceType::RangeFieldType DofType;
496  typedef typename MatrixType::row_type RowType;
497 
499  typedef typename DomainSpaceType::BlockMapperType ColMapperType;
501  typedef typename RangeSpaceType::BlockMapperType RowMapperType;
502 
503  static const int littleCols = MatrixObjectType::littleCols;
504  static const int littleRows = MatrixObjectType::littleRows;
505 
506  LocalMatrix ( const MatrixObjectType& mObj, const DomainSpaceType& domainSpace, const RangeSpaceType& rangeSpace )
507  : BaseType( domainSpace, rangeSpace ),
508  rowMapper_( rangeSpace.blockMapper() ),
509  colMapper_( domainSpace.blockMapper() ),
510  numRows_( rowMapper_.maxNumDofs() ),
511  numCols_( colMapper_.maxNumDofs() ),
512  matrixObj_( mObj )
513  {}
514 
515  LocalMatrix ( const LocalMatrix& org )
516  : BaseType( org ),
517  rowMapper_(org.rowMapper_),
518  colMapper_(org.colMapper_),
519  numRows_( org.numRows_ ),
520  numCols_( org.numCols_ ),
521  matrices_(org.matrices_),
522  matrixObj_(org.matrixObj_)
523  {}
524 
525 
527  void init ( const DomainEntityType &domainEntity, const RangeEntityType &rangeEntity )
528  {
529  // initialize base functions sets
530  BaseType :: init ( domainEntity, rangeEntity );
531 
532  numRows_ = rowMapper_.numDofs( rangeEntity );
533  numCols_ = colMapper_.numDofs( domainEntity );
534  matrices_.resize( numRows_, numCols_, nullptr );
535 
536  typedef typename MatrixType::size_type Index;
537  auto blockAccess = [ this ] ( const std::pair< Index, Index > &index ) -> LittleBlockType&
538  {
539  if( matrixObj_.implicitModeActive() )
540  return matrixObj_.matrix().entry( index.first, index.second );
541  else
542  return matrixObj_.matrix()[ index.first ][ index.second ];
543  };
544 
545  auto functor = [ this, blockAccess ] ( std::pair< int, int > local, const std::pair< Index, Index > &index )
546  {
547  matrices_[ local.first ][ local.second ] = &blockAccess( index );
548  };
549 
550  rowMapper_.mapEach( rangeEntity, makePairFunctor( colMapper_, domainEntity, functor ) );
551  }
552 
553  private:
554  // check whether given (row,col) pair is valid
555  void check(int localRow, int localCol) const
556  {
557  const std::size_t row = (int) localRow / littleRows;
558  const std::size_t col = (int) localCol / littleCols;
559  const int lRow = localRow%littleRows;
560  const int lCol = localCol%littleCols;
561  assert( row < matrices_.size() ) ;
562  assert( col < matrices_[row].size() );
563  assert( lRow < littleRows );
564  assert( lCol < littleCols );
565  }
566 
567  DofType& getValue(const int localRow, const int localCol)
568  {
569  const int row = (int) localRow / littleRows;
570  const int col = (int) localCol / littleCols;
571  const int lRow = localRow%littleRows;
572  const int lCol = localCol%littleCols;
573  return (*matrices_[row][col])[lRow][lCol];
574  }
575 
576  public:
577  const DofType get(const int localRow, const int localCol) const
578  {
579  const int row = (int) localRow / littleRows;
580  const int col = (int) localCol / littleCols;
581  const int lRow = localRow%littleRows;
582  const int lCol = localCol%littleCols;
583  return (*matrices_[row][col])[lRow][lCol];
584  }
585 
586  void scale (const DofType& scalar)
587  {
588  for(auto i=0; i<matrices_.size(); ++i)
589  for(auto j=0; j<matrices_[i].size(); ++j)
590  (*matrices_[i][j]) *= scalar;
591  }
592 
593  void add(const int localRow, const int localCol , const DofType value)
594  {
595 #ifndef NDEBUG
596  check(localRow,localCol);
597 #endif
598  getValue(localRow,localCol) += value;
599  }
600 
601  void set(const int localRow, const int localCol , const DofType value)
602  {
603 #ifndef NDEBUG
604  check(localRow,localCol);
605 #endif
606  getValue(localRow,localCol) = value;
607  }
608 
610  void unitRow(const int localRow)
611  {
612  const int row = (int) localRow / littleRows;
613  const int lRow = localRow%littleRows;
614 
615  // clear row
616  doClearRow( row, lRow );
617 
618  // set diagonal entry to 1
619  (*matrices_[row][row])[lRow][lRow] = 1;
620  }
621 
623  void clear ()
624  {
625  for(auto i=0; i<matrices_.size(); ++i)
626  for(auto j=0; j<matrices_[i].size(); ++j)
627  (*matrices_[i][j]) = (DofType) 0;
628  }
629 
631  void clearRow ( const int localRow )
632  {
633  const int row = (int) localRow / littleRows;
634  const int lRow = localRow%littleRows;
635 
636  // clear the row
637  doClearRow( row, lRow );
638  }
639 
641  void resort ()
642  {}
643 
644  protected:
646  void doClearRow ( const int row, const int lRow )
647  {
648  // get number of columns
649  const auto col = this->columns();
650  for(auto localCol=0; localCol<col; ++localCol)
651  {
652  const int col = (int) localCol / littleCols;
653  const int lCol = localCol%littleCols;
654  (*matrices_[row][col])[lRow][lCol] = 0;
655  }
656  }
657 
658  private:
659  // special mapper omiting block size
660  const RowMapperType& rowMapper_;
661  const ColMapperType& colMapper_;
662 
663  // number of local matrices
664  int numRows_;
665  int numCols_;
666 
667  // dynamic matrix with pointers to block matrices
668  Dune::DynamicMatrix< LittleBlockType* > matrices_;
669 
670  // matrix to build
671  const MatrixObjectType& matrixObj_;
672  }; // end of class LocalMatrix
673 
675  typedef LocalMatrix<ThisType> ObjectType;
676  typedef ThisType LocalMatrixFactoryType;
677  typedef ObjectStack< LocalMatrixFactoryType > LocalMatrixStackType;
679  typedef LocalMatrixWrapper< LocalMatrixStackType > LocalMatrixType;
680  typedef ColumnObject< ThisType > LocalColumnObjectType;
681 
682  protected:
683  const DomainSpaceType & domainSpace_;
684  const RangeSpaceType & rangeSpace_;
685 
686  // sepcial row mapper
687  RowMapperType& rowMapper_;
688  // special col mapper
689  ColMapperType& colMapper_;
690 
691  int size_;
692 
693  int sequence_;
694 
695  mutable MatrixType* matrix_;
696 
697  mutable LocalMatrixStackType localMatrixStack_;
698 
699  mutable MatrixAdapterType* matrixAdap_;
700  mutable ColumnBlockVectorType* Arg_;
701  mutable RowBlockVectorType* Dest_;
702  // overflow fraction for implicit build mode
703  const double overflowFraction_;
704 
705  // prohibit copy constructor
706  ISTLMatrixObject(const ISTLMatrixObject&);
707  public:
708 
712  ISTLMatrixObject ( const DomainSpaceType &domainSpace, const RangeSpaceType &rangeSpace, const std::string& paramfile )
713  DUNE_DEPRECATED_MSG("ISTLMatrixObject(...,string) is deprecated. Use ISTLMatrixObject(DomainSpace,RangeSpace,MatrixParameter) instead")
714  : domainSpace_(domainSpace)
715  , rangeSpace_(rangeSpace)
716  , rowMapper_( rangeSpace.blockMapper() )
717  , colMapper_( domainSpace.blockMapper() )
718  , size_(-1)
719  , sequence_(-1)
720  , matrix_(nullptr)
721  // , scp_(rangeSpace())
722  , localMatrixStack_( *this )
723  , matrixAdap_(nullptr)
724  , Arg_(nullptr)
725  , Dest_(nullptr)
726  , overflowFraction_( Parameter::getValue( "istl.matrix.overflowfraction", 1.0 ) )
727  {}
728 
732  ISTLMatrixObject ( const DomainSpaceType &domainSpace, const RangeSpaceType &rangeSpace, const MatrixParameter& param = ISTLMatrixParameter() ) :
733  domainSpace_(domainSpace)
734  , rangeSpace_(rangeSpace)
735  , rowMapper_( rangeSpace.blockMapper() )
736  , colMapper_( domainSpace.blockMapper() )
737  , size_(-1)
738  , sequence_(-1)
739  , matrix_(nullptr)
740  , localMatrixStack_( *this )
741  , matrixAdap_(nullptr)
742  , Arg_(nullptr)
743  , Dest_(nullptr)
744  , overflowFraction_( param.overflowFraction() )
745  {}
746 
748  static const bool assembled = true ;
749 
750  const ThisType& systemMatrix() const
751  {
752  return *this;
753  }
754  ThisType& systemMatrix()
755  {
756  return *this;
757  }
758 
760  ~ISTLMatrixObject()
761  {
762  removeObj( true );
763  }
764 
766  MatrixType & matrix() const
767  {
768  assert( matrix_ );
769  return *matrix_;
770  }
771 
772  void printTexInfo(std::ostream& out) const
773  {
774  out << "ISTL MatrixObj: ";
775  out << "\\\\ \n";
776  }
777 
779  std::string preconditionName() const
780  {
781  return "";
782  }
783 
784  void createMatrixAdapter () const
785  {
786  if( !matrixAdap_ )
787  matrixAdap_ = new MatrixAdapterType( matrixAdapterObject() );
788  assert( matrixAdap_ );
789  }
790 
792  const MatrixAdapterType& matrixAdapter() const
793  {
794  createMatrixAdapter();
795  return *matrixAdap_;
796  }
797  MatrixAdapterType& matrixAdapter()
798  {
799  createMatrixAdapter();
800  return *matrixAdap_;
801  }
802 
803  protected:
804  MatrixAdapterType matrixAdapterObject() const
805  {
806  // need some precondition object - empty here
807  typedef typename MatrixAdapterType :: PreconditionAdapterType PreConType;
808  return MatrixAdapterType(matrix(), domainSpace(), rangeSpace(), PreConType() );
809  }
810 
811  public:
812  bool implicitModeActive() const
813  {
814  // implicit build mode is only active when the
815  // build mode of the matrix is implicit and the
816  // matrix is currently being build
817 
818  if( matrix().buildMode() == MatrixType::implicit && matrix().buildStage() == MatrixType::building )
819  return true;
820  else
821  return false;
822  }
823 
824  // compress matrix if not already done before and only in implicit build mode
825  void communicate( )
826  {
827  if( implicitModeActive() )
828  matrix().compress();
829  }
830 
833  bool hasPreconditionMatrix() const
834  {
835  return true;
836  }
837 
839  void clear()
840  {
841  matrix().clear();
842  // clean matrix adapter and other helper classes
843  removeObj( false );
844  }
845 
847  template <class Stencil>
848  void reserve(const Stencil &stencil, const bool implicit = true )
849  {
850  // if grid sequence number changed, rebuild matrix
851  if(sequence_ != domainSpace().sequence())
852  {
853  removeObj( true );
854 
855  if( implicit )
856  {
857  auto nnz = stencil.maxNonZerosEstimate();
858  if( nnz == 0 )
859  {
860  Stencil tmpStencil( stencil );
861  tmpStencil.fill( *(domainSpace_.begin()), *(rangeSpace_.begin()) );
862  nnz = tmpStencil.maxNonZerosEstimate();
863  }
864  matrix_ = new MatrixType( rowMapper_.size(), colMapper_.size(), nnz, overflowFraction_ );
865  }
866  else
867  {
868  matrix_ = new MatrixType( rowMapper_.size(), colMapper_.size() );
869  matrix().createEntries( stencil.globalStencil() );
870  }
871 
872  sequence_ = domainSpace().sequence();
873  }
874  }
875 
877  template <class HangingNodesType>
878  void changeHangingNodes(const HangingNodesType& hangingNodes)
879  {
880  // create new matrix
881  MatrixType* newMatrix = new MatrixType(rowMapper_.size(), colMapper_.size());
882 
883  // setup with hanging rows
884  newMatrix->setup( *matrix_ , hangingNodes );
885 
886  // remove old matrix
887  removeObj( true );
888 
889  // store new matrix
890  matrix_ = newMatrix;
891  }
892 
895  void extractDiagonal( ColumnDiscreteFunctionType& diag ) const
896  {
897  // extract diagonal entries
898  matrix().extractDiagonal( diag.blockVector() );
899  }
900 
902  bool rightPrecondition() const
903  {
904  return false;
905  }
906 
908  void multOEM(const double* arg, double* dest) const
909  {
910  createBlockVectors();
911 
912  RowBlockVectorType& Arg = *Arg_;
913  ColumnBlockVectorType &Dest = *Dest_;
914 
915  std::copy_n( arg, Arg.size(), Arg.dbegin() );
916 
917  // call mult of matrix adapter
918  matrixAdapter().apply( Arg, Dest );
919 
920  std::copy( Dest.dbegin(), Dest.dend(), dest );
921  }
922 
924  void apply(const ColumnDiscreteFunctionType& arg, RowDiscreteFunctionType& dest) const
925  {
926  matrixAdapter().apply( arg.blockVector(), dest.blockVector() );
927  }
928 
930  template <class RowDFType, class ColDFType>
931  void apply(const RowDFType& arg, ColDFType& dest) const
932  {
933  multOEM( arg.dofVector(), dest.dofVector ());
934  }
935 
937  template <class ColumnLeakPointerType, class RowLeakPointerType>
938  void multOEM(const ColumnLeakPointerType& arg, RowLeakPointerType& dest) const
939  {
940  DUNE_THROW(NotImplemented,"Method has been removed");
941  }
942 
944  double ddotOEM(const double* v, const double* w) const
945  {
946  createBlockVectors();
947 
948  RowBlockVectorType& V = *Arg_;
949  ColumnBlockVectorType& W = *Dest_;
950 
951  std::copy_n( v, V.size(), V.dbegin() );
952  std::copy_n( w, W.size(), W.dbegin() );
953 
954 #if HAVE_MPI
955  // in parallel use scalar product of discrete functions
956  ISTLBlockVectorDiscreteFunction< DomainSpaceType > vF("ddotOEM:vF", domainSpace(), V );
957  ISTLBlockVectorDiscreteFunction< RangeSpaceType > wF("ddotOEM:wF", rangeSpace(), W );
958  return vF.scalarProductDofs( wF );
959 #else
960  return V * W;
961 #endif
962  }
963 
965  void resort()
966  {}
967 
970  void createPreconditionMatrix()
971  {}
972 
974  void print(std::ostream & s) const
975  {
976  matrix().print(s);
977  }
978 
979  const DomainSpaceType& domainSpace() const
980  {
981  return domainSpace_;
982  }
983  const RangeSpaceType& rangeSpace() const
984  {
985  return rangeSpace_;
986  }
987 
988  const RowMapperType& rowMapper() const
989  {
990  return rowMapper_;
991  }
992  const ColMapperType& colMapper() const
993  {
994  return colMapper_;
995  }
996 
998  ObjectType* newObject() const
999  {
1000  return new ObjectType(*this, domainSpace(), rangeSpace());
1001  }
1002 
1004  LocalMatrixType localMatrix( const DomainEntityType &domainEntity, const RangeEntityType &rangeEntity ) const
1005  {
1006  return LocalMatrixType( localMatrixStack_, domainEntity, rangeEntity );
1007  }
1008  LocalColumnObjectType localColumn( const DomainEntityType &domainEntity ) const
1009  {
1010  return LocalColumnObjectType ( *this, domainEntity );
1011  }
1012 
1013  template< class LocalMatrix >
1014  void addLocalMatrix ( const DomainEntityType &domainEntity, const RangeEntityType &rangeEntity, const LocalMatrix &localMat )
1015  {
1016  typedef typename MatrixType::size_type Index;
1017  auto blockAccess = [ this ] ( const std::pair< Index, Index > &index ) -> LittleBlockType&
1018  {
1019  if( implicitModeActive() )
1020  return matrix().entry( index.first, index.second );
1021  else
1022  return matrix()[ index.first ][ index.second ];
1023  };
1024 
1025  auto functor = [ &localMat, this, blockAccess ] ( std::pair< int, int > local, const std::pair< Index, Index > &index )
1026  {
1027  LittleBlockType& block = blockAccess( index );
1028  for( std::size_t i = 0; i < littleRows; ++i )
1029  for( std::size_t j = 0; j < littleCols; ++j )
1030  block[ i ][ j ] += localMat.get( local.first * littleRows + i, local.second *littleCols + j );
1031  };
1032 
1033  rowMapper_.mapEach( rangeEntity, makePairFunctor( colMapper_, domainEntity, functor ) );
1034  }
1035 
1036 
1037  template< class LocalMatrix >
1038  void setLocalMatrix ( const DomainEntityType &domainEntity, const RangeEntityType &rangeEntity, const LocalMatrix &localMat )
1039  {
1040  typedef typename MatrixType::size_type Index;
1041  auto blockAccess = [ this ] ( const std::pair< Index, Index > &index ) -> LittleBlockType&
1042  {
1043  if( implicitModeActive() )
1044  return matrix().entry( index.first, index.second );
1045  else
1046  return matrix()[ index.first ][ index.second ];
1047  };
1048 
1049  auto functor = [ &localMat, this, blockAccess ] ( std::pair< int, int > local, const std::pair< Index, Index > &index )
1050  {
1051  LittleBlockType& block = blockAccess( index );
1052  for( std::size_t i = 0; i < littleRows; ++i )
1053  for( std::size_t j = 0; j < littleCols; ++j )
1054  block[ i ][ j ] = localMat.get( local.first * littleRows + i, local.second *littleCols + j );
1055  };
1056 
1057  rowMapper_.mapEach( rangeEntity, makePairFunctor( colMapper_, domainEntity, functor ) );
1058  }
1059 
1060  template< class LocalMatrix, class Scalar >
1061  void addScaledLocalMatrix ( const DomainEntityType &domainEntity, const RangeEntityType &rangeEntity, const LocalMatrix &localMat, const Scalar &s )
1062  {
1063  typedef typename MatrixType::size_type Index;
1064  auto blockAccess = [ this ] ( const std::pair< Index, Index > &index ) -> LittleBlockType&
1065  {
1066  if( implicitModeActive() )
1067  return matrix().entry( index.first, index.second );
1068  else
1069  return matrix()[ index.first ][ index.second ];
1070  };
1071 
1072  auto functor = [ &localMat, &s, this, blockAccess ] ( std::pair< int, int > local, const std::pair< Index, Index > &index )
1073  {
1074  LittleBlockType& block = blockAccess( index );
1075  for( std::size_t i = 0; i < littleRows; ++i )
1076  for( std::size_t j = 0; j < littleCols; ++j )
1077  block[ i ][ j ] += s * localMat.get( local.first * littleRows + i, local.second *littleCols + j );
1078  };
1079 
1080  rowMapper_.mapEach( rangeEntity, makePairFunctor( colMapper_, domainEntity, functor ) );
1081  }
1082 
1083  template< class LocalMatrix >
1084  void getLocalMatrix ( const DomainEntityType &domainEntity, const RangeEntityType &rangeEntity, LocalMatrix &localMat ) const
1085  {
1086  typedef typename MatrixType::size_type Index;
1087  auto functor = [ &localMat, this ] ( std::pair< int, int > local, const std::pair< Index, Index > &global )
1088  {
1089  for( std::size_t i = 0; i < littleRows; ++i )
1090  for( std::size_t j = 0; j < littleCols; ++j )
1091  localMat.set( local.first * littleRows + i, local.second *littleCols + j,
1092  matrix()[ global.first ][ global.second ][i][j] );
1093  };
1094 
1095  rowMapper_.mapEach( rangeEntity, makePairFunctor( colMapper_, domainEntity, functor ) );
1096  }
1097 
1098  protected:
1099  void preConErrorMsg(int preCon) const
1100  {
1101  exit(1);
1102  }
1103 
1104  void removeObj( const bool alsoClearMatrix )
1105  {
1106  delete Dest_;
1107  Dest_ = nullptr;
1108  delete Arg_;
1109  Arg_ = nullptr;
1110  delete matrixAdap_;
1111  matrixAdap_ = nullptr;
1112 
1113  if( alsoClearMatrix )
1114  {
1115  delete matrix_;
1116  matrix_ = nullptr;
1117  }
1118  }
1119 
1120  void createBlockVectors() const
1121  {
1122  if( Arg_ == nullptr || Dest_ == nullptr )
1123  {
1124  delete Arg_;
1125  delete Dest_;
1126  Arg_ = new RowBlockVectorType( rowMapper_.size() );
1127  Dest_ = new ColumnBlockVectorType( colMapper_.size() );
1128  }
1129  createMatrixAdapter ();
1130  }
1131 
1132  };
1133 
1135  template <class SpaceImp>
1136  class ISTLMatrixObject<SpaceImp,SpaceImp>
1137  {
1138  public:
1139  typedef SpaceImp DomainSpaceImp;
1140  typedef SpaceImp RangeSpaceImp;
1142  typedef DomainSpaceImp DomainSpaceType;
1144  typedef RangeSpaceImp RangeSpaceType;
1145 
1147  typedef ISTLMatrixObject<DomainSpaceType,RangeSpaceType> ThisType;
1148 
1149  typedef typename DomainSpaceType::GridType GridType;
1150 
1151  typedef typename RangeSpaceType :: EntityType RangeEntityType ;
1152  typedef typename DomainSpaceType :: EntityType DomainEntityType;
1153 
1154  enum { littleRows = DomainSpaceType :: localBlockSize };
1155  enum { littleCols = RangeSpaceType :: localBlockSize };
1156 
1157  typedef FieldMatrix<typename DomainSpaceType :: RangeFieldType, littleRows, littleCols> LittleBlockType;
1158 
1159  typedef ISTLBlockVectorDiscreteFunction< DomainSpaceType > RowDiscreteFunctionType;
1160  typedef ISTLBlockVectorDiscreteFunction< RangeSpaceType > ColumnDiscreteFunctionType;
1161 
1162  protected:
1163  typedef typename RowDiscreteFunctionType :: DofContainerType RowBlockVectorType;
1164  typedef typename ColumnDiscreteFunctionType :: DofContainerType ColumnBlockVectorType;
1165 
1166  typedef typename DomainSpaceType :: BlockMapperType RowMapperType;
1167  typedef typename RangeSpaceType :: BlockMapperType ColMapperType;
1168 
1169  public:
1171  typedef ImprovedBCRSMatrix< LittleBlockType , RowDiscreteFunctionType , ColumnDiscreteFunctionType > MatrixType;
1172  typedef typename ISTLParallelMatrixAdapter<MatrixType,DomainSpaceType>::Type MatrixAdapterType;
1173  // get preconditioner type from MatrixAdapterType
1174  typedef ThisType PreconditionMatrixType;
1175  typedef typename MatrixAdapterType :: ParallelScalarProductType ParallelScalarProductType;
1176 
1177  template <class MatrixObjectImp>
1178  class LocalMatrix;
1179 
1180  struct LocalMatrixTraits
1181  {
1182  typedef DomainSpaceImp DomainSpaceType ;
1183  typedef RangeSpaceImp RangeSpaceType;
1184  typedef typename DomainSpaceType :: RangeFieldType RangeFieldType;
1185  typedef LocalMatrix<ThisType> LocalMatrixType;
1186  typedef typename MatrixType:: block_type LittleBlockType;
1187  };
1188 
1190  template <class MatrixObjectImp>
1191  class LocalMatrix : public LocalMatrixDefault<LocalMatrixTraits>
1192  {
1193  public:
1195  typedef LocalMatrixDefault<LocalMatrixTraits> BaseType;
1196 
1198  typedef MatrixObjectImp MatrixObjectType;
1200  typedef typename MatrixObjectImp :: MatrixType MatrixType;
1202  typedef typename MatrixType::block_type LittleBlockType;
1203 
1204  typedef typename MatrixObjectType::DomainSpaceType DomainSpaceType;
1205  typedef typename MatrixObjectType::RangeSpaceType RangeSpaceType;
1206 
1207  typedef typename MatrixObjectType::DomainEntityType DomainEntityType;
1208  typedef typename MatrixObjectType::RangeEntityType RangeEntityType;
1209 
1211  typedef typename DomainSpaceType::RangeFieldType DofType;
1212  typedef typename MatrixType::row_type RowType;
1213 
1215  typedef typename DomainSpaceType::BlockMapperType ColMapperType;
1217  typedef typename RangeSpaceType::BlockMapperType RowMapperType;
1218 
1219  static const int littleCols = MatrixObjectType::littleCols;
1220  static const int littleRows = MatrixObjectType::littleRows;
1221 
1222  LocalMatrix ( const MatrixObjectType& mObj, const DomainSpaceType& domainSpace, const RangeSpaceType& rangeSpace )
1223  : BaseType( domainSpace, rangeSpace ),
1224  rowMapper_( rangeSpace.blockMapper() ),
1225  colMapper_( domainSpace.blockMapper() ),
1226  numRows_( rowMapper_.maxNumDofs() ),
1227  numCols_( colMapper_.maxNumDofs() ),
1228  matrixObj_( mObj )
1229  {}
1230 
1231  LocalMatrix ( const LocalMatrix& org )
1232  : BaseType( org ),
1233  rowMapper_(org.rowMapper_),
1234  colMapper_(org.colMapper_),
1235  numRows_( org.numRows_ ),
1236  numCols_( org.numCols_ ),
1237  matrices_(org.matrices_),
1238  matrixObj_(org.matrixObj_)
1239  {}
1240 
1241 
1243  void init ( const DomainEntityType &domainEntity, const RangeEntityType &rangeEntity )
1244  {
1245  // initialize base functions sets
1246  BaseType :: init ( domainEntity, rangeEntity );
1247 
1248  numRows_ = rowMapper_.numDofs( rangeEntity );
1249  numCols_ = colMapper_.numDofs( domainEntity );
1250  matrices_.resize( numRows_, numCols_, nullptr );
1251 
1252  typedef typename MatrixType::size_type Index;
1253  auto blockAccess = [ this ] ( const std::pair< Index, Index > &index ) -> LittleBlockType&
1254  {
1255  if( matrixObj_.implicitModeActive() )
1256  return matrixObj_.matrix().entry( index.first, index.second );
1257  else
1258  return matrixObj_.matrix()[ index.first ][ index.second ];
1259  };
1260 
1261  auto functor = [ this, blockAccess ] ( std::pair< int, int > local, const std::pair< Index, Index > &index )
1262  {
1263  matrices_[ local.first ][ local.second ] = &blockAccess( index );
1264  };
1265 
1266  rowMapper_.mapEach( rangeEntity, makePairFunctor( colMapper_, domainEntity, functor ) );
1267  }
1268 
1269  private:
1270  // check whether given (row,col) pair is valid
1271  void check(int localRow, int localCol) const
1272  {
1273  const std::size_t row = (int) localRow / littleRows;
1274  const std::size_t col = (int) localCol / littleCols;
1275  const int lRow = localRow%littleRows;
1276  const int lCol = localCol%littleCols;
1277  assert( row < matrices_.size() ) ;
1278  assert( col < matrices_[row].size() );
1279  assert( lRow < littleRows );
1280  assert( lCol < littleCols );
1281  }
1282 
1283  DofType& getValue(const int localRow, const int localCol)
1284  {
1285  const int row = (int) localRow / littleRows;
1286  const int col = (int) localCol / littleCols;
1287  const int lRow = localRow%littleRows;
1288  const int lCol = localCol%littleCols;
1289  return (*matrices_[row][col])[lRow][lCol];
1290  }
1291 
1292  public:
1293  const DofType get(const int localRow, const int localCol) const
1294  {
1295  const int row = (int) localRow / littleRows;
1296  const int col = (int) localCol / littleCols;
1297  const int lRow = localRow%littleRows;
1298  const int lCol = localCol%littleCols;
1299  return (*matrices_[row][col])[lRow][lCol];
1300  }
1301 
1302  void scale (const DofType& scalar)
1303  {
1304  for(auto i=0; i<matrices_.size(); ++i)
1305  for(auto j=0; j<matrices_[i].size(); ++j)
1306  (*matrices_[i][j]) *= scalar;
1307  }
1308 
1309  void add(const int localRow, const int localCol , const DofType value)
1310  {
1311 #ifndef NDEBUG
1312  check(localRow,localCol);
1313 #endif
1314  getValue(localRow,localCol) += value;
1315  }
1316 
1317  void set(const int localRow, const int localCol , const DofType value)
1318  {
1319 #ifndef NDEBUG
1320  check(localRow,localCol);
1321 #endif
1322  getValue(localRow,localCol) = value;
1323  }
1324 
1326  void unitRow(const int localRow)
1327  {
1328  const int row = (int) localRow / littleRows;
1329  const int lRow = localRow%littleRows;
1330 
1331  // clear row
1332  doClearRow( row, lRow );
1333 
1334  // set diagonal entry to 1
1335  (*matrices_[row][row])[lRow][lRow] = 1;
1336  }
1337 
1339  void clear ()
1340  {
1341  for(auto i=0; i<matrices_.size(); ++i)
1342  for(auto j=0; j<matrices_[i].size(); ++j)
1343  (*matrices_[i][j]) = (DofType) 0;
1344  }
1345 
1347  void clearRow ( const int localRow )
1348  {
1349  const int row = (int) localRow / littleRows;
1350  const int lRow = localRow%littleRows;
1351 
1352  // clear the row
1353  doClearRow( row, lRow );
1354  }
1355 
1357  void resort ()
1358  {}
1359 
1360  protected:
1362  void doClearRow ( const int row, const int lRow )
1363  {
1364  // get number of columns
1365  const auto col = this->columns();
1366  for(auto localCol=0; localCol<col; ++localCol)
1367  {
1368  const int col = (int) localCol / littleCols;
1369  const int lCol = localCol%littleCols;
1370  (*matrices_[row][col])[lRow][lCol] = 0;
1371  }
1372  }
1373 
1374  private:
1375  // special mapper omiting block size
1376  const RowMapperType& rowMapper_;
1377  const ColMapperType& colMapper_;
1378 
1379  // number of local matrices
1380  int numRows_;
1381  int numCols_;
1382 
1383  // dynamic matrix with pointers to block matrices
1384  Dune::DynamicMatrix< LittleBlockType* > matrices_;
1385 
1386  // matrix to build
1387  const MatrixObjectType& matrixObj_;
1388  }; // end of class LocalMatrix
1389 
1391  typedef LocalMatrix<ThisType> ObjectType;
1392  typedef ThisType LocalMatrixFactoryType;
1393  typedef ObjectStack< LocalMatrixFactoryType > LocalMatrixStackType;
1395  typedef LocalMatrixWrapper< LocalMatrixStackType > LocalMatrixType;
1396  typedef ColumnObject< ThisType > LocalColumnObjectType;
1397 
1398  protected:
1399  const DomainSpaceType & domainSpace_;
1400  const RangeSpaceType & rangeSpace_;
1401 
1402  // sepcial row mapper
1403  RowMapperType& rowMapper_;
1404  // special col mapper
1405  ColMapperType& colMapper_;
1406 
1407  int size_;
1408 
1409  int sequence_;
1410 
1411  mutable MatrixType* matrix_;
1412 
1413  ParallelScalarProductType scp_;
1414 
1415  int numIterations_;
1416  double relaxFactor_;
1417 
1418  enum ISTLPreConder_Id { none = 0 , // no preconditioner
1419  ssor = 1 , // SSOR preconditioner
1420  sor = 2 , // SOR preconditioner
1421  ilu_0 = 3 , // ILU-0 preconditioner
1422  ilu_n = 4 , // ILU-n preconditioner
1423  gauss_seidel= 5, // Gauss-Seidel preconditioner
1424  jacobi = 6, // Jacobi preconditioner
1425  amg_ilu_0 = 7, // AMG with ILU-0 smoother
1426  amg_ilu_n = 8, // AMG with ILU-n smoother
1427  amg_jacobi = 9 // AMG with Jacobi smoother
1428  };
1429 
1430  ISTLPreConder_Id preconditioning_;
1431 
1432  mutable LocalMatrixStackType localMatrixStack_;
1433 
1434  mutable MatrixAdapterType* matrixAdap_;
1435  mutable RowBlockVectorType* Arg_;
1436  mutable ColumnBlockVectorType* Dest_;
1437  // overflow fraction for implicit build mode
1438  const double overflowFraction_;
1439  const MatrixParameter& param_;
1440 
1441  public:
1442  ISTLMatrixObject(const ISTLMatrixObject&) = delete;
1443 
1444  ISTLMatrixObject ( const DomainSpaceType &rowSpace, const RangeSpaceType &colSpace, const std :: string &paramfile )
1445  DUNE_DEPRECATED_MSG("ISTLMatrixObject(...,string) is deprecated. Use ISTLMatrixObject(DomainSpace,RangeSpace,ISTLMatrixParameter) instead")
1446  : domainSpace_(rowSpace)
1447  , rangeSpace_(colSpace)
1448  , colMapper_( colSpace.blockMapper() )
1449  , size_(-1)
1450  , sequence_(-1)
1451  , matrix_(nullptr)
1452  , scp_(rangeSpace())
1453  , numIterations_(5)
1454  , relaxFactor_(1.1)
1455  , preconditioning_(none)
1456  , localMatrixStack_( *this )
1457  , matrixAdap_(nullptr)
1458  , Arg_(nullptr)
1459  , Dest_(nullptr)
1460  , overflowFraction_( Parameter::getValue( "istl.matrix.overflowfraction", 1.0 ) )
1461  {
1462  int preCon = 0;
1463  if(paramfile != "")
1464  DUNE_THROW(InvalidStateException,"ISTLMatrixObject: old parameter method disabled");
1465  else
1466  {
1467  static const std::string preConTable[]
1468  = { "none", "ssor", "sor", "ilu-0", "ilu-n", "gauss-seidel", "jacobi", "amg-ilu-0", "amg-ilu-n", "amg-jacobi" };
1469  preCon = Parameter::getEnum( "istl.preconditioning.method", preConTable, preCon );
1470  numIterations_ = Parameter::getValue( "istl.preconditioning.iterations", numIterations_ );
1471  relaxFactor_ = Parameter::getValue( "istl.preconditioning.relaxation", relaxFactor_ );
1472  }
1473 
1474  if( preCon >= 0 && preCon <= 9)
1475  preconditioning_ = (ISTLPreConder_Id) preCon;
1476  else
1477  preConErrorMsg(preCon);
1478 
1479  assert( rowMapper_.size() == colMapper_.size() );
1480  }
1481 
1489  ISTLMatrixObject ( const DomainSpaceType &rowSpace, const RangeSpaceType &colSpace, const MatrixParameter& param = ISTLMatrixParameter() )
1490  : domainSpace_(rowSpace)
1491  , rangeSpace_(colSpace)
1492  // create scp to have at least one instance
1493  // otherwise instance will be deleted during setup
1494  // get new mappers with number of dofs without considerung block size
1495  , rowMapper_( rowSpace.blockMapper() )
1496  , colMapper_( colSpace.blockMapper() )
1497  , size_(-1)
1498  , sequence_(-1)
1499  , matrix_(nullptr)
1500  , scp_(rangeSpace())
1501  , numIterations_( param.numIterations() )
1502  , relaxFactor_( param.relaxation() )
1503  , preconditioning_( (ISTLPreConder_Id)param.method() )
1504  , localMatrixStack_( *this )
1505  , matrixAdap_(nullptr)
1506  , Arg_(nullptr)
1507  , Dest_(nullptr)
1508  , overflowFraction_( param.overflowFraction() )
1509  , param_( param )
1510  {
1511  assert( rowMapper_.size() == colMapper_.size() );
1512  }
1513 
1515  static const bool assembled = true ;
1516 
1517  const ThisType& systemMatrix() const
1518  {
1519  return *this;
1520  }
1521  ThisType& systemMatrix()
1522  {
1523  return *this;
1524  }
1525  public:
1527  ~ISTLMatrixObject()
1528  {
1529  removeObj( true );
1530  }
1531 
1533  MatrixType & matrix() const
1534  {
1535  assert( matrix_ );
1536  return *matrix_;
1537  }
1538 
1539  void printTexInfo(std::ostream& out) const
1540  {
1541  out << "ISTL MatrixObj: ";
1542  out << " preconditioner = " << preconditionName() ;
1543  out << "\\\\ \n";
1544  }
1545 
1547  std::string preconditionName() const
1548  {
1549  return param_.preconditionName();
1550  }
1551 
1552  template <class PreconditionerType>
1553  MatrixAdapterType createMatrixAdapter(const PreconditionerType* preconditioning, std::size_t numIterations) const
1554  {
1555  typedef typename MatrixAdapterType :: PreconditionAdapterType PreConType;
1556  PreConType preconAdapter(matrix(), numIterations, relaxFactor_, preconditioning );
1557  return MatrixAdapterType(matrix(), domainSpace(), rangeSpace(), preconAdapter );
1558  }
1559 
1560  template <class PreconditionerType>
1561  MatrixAdapterType createAMGMatrixAdapter(const PreconditionerType* preconditioning, std::size_t numIterations) const
1562  {
1563  typedef typename MatrixAdapterType :: PreconditionAdapterType PreConType;
1564  PreConType preconAdapter(matrix(), numIterations, relaxFactor_, preconditioning, domainSpace().gridPart().comm() );
1565  return MatrixAdapterType(matrix(), domainSpace(), rangeSpace(), preconAdapter );
1566  }
1567 
1569  const MatrixAdapterType& matrixAdapter() const
1570  {
1571  if( matrixAdap_ == nullptr )
1572  matrixAdap_ = new MatrixAdapterType( matrixAdapterObject() );
1573  return *matrixAdap_;
1574  }
1575  MatrixAdapterType& matrixAdapter()
1576  {
1577  if( matrixAdap_ == nullptr )
1578  matrixAdap_ = new MatrixAdapterType( matrixAdapterObject() );
1579  return *matrixAdap_;
1580  }
1581 
1582  protected:
1583  MatrixAdapterType matrixAdapterObject() const
1584  {
1585 #ifndef DISABLE_ISTL_PRECONDITIONING
1586  const auto procs = domainSpace().gridPart().comm().size();
1587 
1588  typedef typename MatrixType :: BaseType ISTLMatrixType;
1589  typedef typename MatrixAdapterType :: PreconditionAdapterType PreConType;
1590  // no preconditioner
1591  if( preconditioning_ == none )
1592  {
1593  return MatrixAdapterType(matrix(), domainSpace(),rangeSpace(), PreConType() );
1594  }
1595  // SSOR
1596  else if( preconditioning_ == ssor )
1597  {
1598  if( procs > 1 )
1599  DUNE_THROW(InvalidStateException,"ISTL::SeqSSOR not working in parallel computations");
1600 
1601  typedef SeqSSOR<ISTLMatrixType,RowBlockVectorType,ColumnBlockVectorType> PreconditionerType;
1602  return createMatrixAdapter( (PreconditionerType*)nullptr, numIterations_ );
1603  }
1604  // SOR
1605  else if(preconditioning_ == sor )
1606  {
1607  if( procs > 1 )
1608  DUNE_THROW(InvalidStateException,"ISTL::SeqSOR not working in parallel computations");
1609 
1610  typedef SeqSOR<ISTLMatrixType,RowBlockVectorType,ColumnBlockVectorType> PreconditionerType;
1611  return createMatrixAdapter( (PreconditionerType*)nullptr, numIterations_ );
1612  }
1613  // ILU-0
1614  else if(preconditioning_ == ilu_0)
1615  {
1616  if( procs > 1 )
1617  DUNE_THROW(InvalidStateException,"ISTL::SeqILU0 not working in parallel computations");
1618 
1619  typedef FemSeqILU0<ISTLMatrixType,RowBlockVectorType,ColumnBlockVectorType> PreconditionerType;
1620  return createMatrixAdapter( (PreconditionerType*)nullptr, numIterations_ );
1621  }
1622  // ILU-n
1623  else if(preconditioning_ == ilu_n)
1624  {
1625  if( procs > 1 )
1626  DUNE_THROW(InvalidStateException,"ISTL::SeqILUn not working in parallel computations");
1627 
1628  typedef SeqILUn<ISTLMatrixType,RowBlockVectorType,ColumnBlockVectorType> PreconditionerType;
1629  return createMatrixAdapter( (PreconditionerType*)nullptr, numIterations_ );
1630  }
1631  // Gauss-Seidel
1632  else if(preconditioning_ == gauss_seidel)
1633  {
1634  if( procs > 1 )
1635  DUNE_THROW(InvalidStateException,"ISTL::SeqGS not working in parallel computations");
1636 
1637  typedef SeqGS<ISTLMatrixType,RowBlockVectorType,ColumnBlockVectorType> PreconditionerType;
1638  return createMatrixAdapter( (PreconditionerType*)nullptr, numIterations_ );
1639  }
1640  // Jacobi
1641  else if(preconditioning_ == jacobi)
1642  {
1643  if( numIterations_ == 1 ) // diagonal preconditioning
1644  {
1645  typedef FemDiagonalPreconditioner< ThisType, RowBlockVectorType, ColumnBlockVectorType > PreconditionerType;
1646  typedef typename MatrixAdapterType :: PreconditionAdapterType PreConType;
1647  PreConType preconAdapter( matrix(), new PreconditionerType( *this ) );
1648  return MatrixAdapterType( matrix(), domainSpace(), rangeSpace(), preconAdapter );
1649  }
1650  else if ( procs == 1 )
1651  {
1652  typedef SeqJac<ISTLMatrixType,RowBlockVectorType,ColumnBlockVectorType> PreconditionerType;
1653  return createMatrixAdapter( (PreconditionerType*)nullptr, numIterations_ );
1654  }
1655  else
1656  {
1657  DUNE_THROW(InvalidStateException,"ISTL::SeqJac(Jacobi) only working with istl.preconditioning.iterations: 1 in parallel computations");
1658  }
1659  }
1660  // AMG ILU-0
1661  else if(preconditioning_ == amg_ilu_0)
1662  {
1663  // use original SeqILU0 because of some AMG traits classes.
1664  typedef SeqILU0<ISTLMatrixType,RowBlockVectorType,ColumnBlockVectorType> PreconditionerType;
1665  return createAMGMatrixAdapter( (PreconditionerType*)nullptr, numIterations_ );
1666  }
1667  // AMG ILU-n
1668  else if(preconditioning_ == amg_ilu_n)
1669  {
1670  typedef SeqILUn<ISTLMatrixType,RowBlockVectorType,ColumnBlockVectorType> PreconditionerType;
1671  return createAMGMatrixAdapter( (PreconditionerType*)nullptr, numIterations_ );
1672  }
1673  // AMG Jacobi
1674  else if(preconditioning_ == amg_jacobi)
1675  {
1676  typedef SeqJac<ISTLMatrixType,RowBlockVectorType,ColumnBlockVectorType,1> PreconditionerType;
1677  return createAMGMatrixAdapter( (PreconditionerType*)nullptr, numIterations_ );
1678  }
1679  else
1680  {
1681  preConErrorMsg(preconditioning_);
1682  }
1683 #endif
1684 
1685  return MatrixAdapterType(matrix(), domainSpace(), rangeSpace(), PreConType() );
1686  }
1687 
1688  public:
1689  bool implicitModeActive() const
1690  {
1691  // implicit build mode is only active when the
1692  // build mode of the matrix is implicit and the
1693  // matrix is currently being build
1694  if( matrix().buildMode() == MatrixType::implicit && matrix().buildStage() == MatrixType::building )
1695  return true;
1696  else
1697  return false;
1698  }
1699 
1700  // compress matrix if not already done before and only in implicit build mode
1701  void communicate( )
1702  {
1703  if( implicitModeActive() )
1704  matrix().compress();
1705  }
1706 
1709  bool hasPreconditionMatrix() const
1710  {
1711  return (preconditioning_ != none);
1712  }
1713 
1715  const PreconditionMatrixType& preconditionMatrix() const
1716  {
1717  return *this;
1718  }
1719 
1721  void clear()
1722  {
1723  matrix().clear();
1724  // clean matrix adapter and other helper classes
1725  removeObj( false );
1726  }
1727 
1729  template <class Stencil>
1730  void reserve( const Stencil &stencil, const bool implicit = true )
1731  {
1732  // if grid sequence number changed, rebuild matrix
1733  if(sequence_ != domainSpace().sequence())
1734  {
1735  removeObj( true );
1736 
1737  if( implicit )
1738  {
1739  auto nnz = stencil.maxNonZerosEstimate();
1740  if( nnz == 0 )
1741  {
1742  Stencil tmpStencil( stencil );
1743  tmpStencil.fill( *(domainSpace_.begin()), *(rangeSpace_.begin()) );
1744  nnz = tmpStencil.maxNonZerosEstimate();
1745  }
1746  matrix_ = new MatrixType( rowMapper_.size(), colMapper_.size(), nnz, overflowFraction_ );
1747  }
1748  else
1749  {
1750  matrix_ = new MatrixType( rowMapper_.size(), colMapper_.size() );
1751  matrix().createEntries( stencil.globalStencil() );
1752  }
1753 
1754  sequence_ = domainSpace().sequence();
1755  }
1756  }
1757 
1759  template <class HangingNodesType>
1760  void changeHangingNodes(const HangingNodesType& hangingNodes)
1761  {
1762  // create new matrix
1763  MatrixType* newMatrix = new MatrixType(rowMapper_.size(), colMapper_.size());
1764 
1765  // setup with hanging rows
1766  newMatrix->setup( *matrix_ , hangingNodes );
1767 
1768  // remove old matrix
1769  removeObj( true );
1770 
1771  // store new matrix
1772  matrix_ = newMatrix;
1773  }
1774 
1777  void extractDiagonal( ColumnDiscreteFunctionType& diag ) const
1778  {
1779  // extract diagonal entries
1780  matrix().extractDiagonal( diag.blockVector() );
1781  }
1782 
1784  bool rightPrecondition() const
1785  {
1786  return true;
1787  }
1788 
1792  void precondition(const double* arg, double* dest) const
1793  {
1794  createBlockVectors();
1795 
1796  assert( Arg_ );
1797  assert( Dest_ );
1798 
1799  RowBlockVectorType& Arg = *Arg_;
1800  ColumnBlockVectorType & Dest = *Dest_;
1801 
1802  std::copy_n( arg, Arg.size(), Arg.dbegin());
1803 
1804  // set Dest to zero
1805  Dest = 0;
1806 
1807  assert( matrixAdap_ );
1808  // not parameter swaped for preconditioner
1809  matrixAdap_->preconditionAdapter().apply(Dest , Arg);
1810 
1811  std::copy( Dest.dbegin(), Dest.dend(), dest );
1812  }
1813 
1815  void multOEM(const double* arg, double* dest) const
1816  {
1817  createBlockVectors();
1818 
1819  assert( Arg_ );
1820  assert( Dest_ );
1821 
1822  RowBlockVectorType& Arg = *Arg_;
1823  ColumnBlockVectorType & Dest = *Dest_;
1824 
1825  std::copy_n( arg, Arg.size(), Arg.dbegin() );
1826 
1827  // call mult of matrix adapter
1828  assert( matrixAdap_ );
1829  matrixAdap_->apply( Arg, Dest );
1830 
1831  std::copy( Dest.dbegin(), Dest.dend(), dest );
1832  }
1833 
1835  void apply(const RowDiscreteFunctionType& arg, ColumnDiscreteFunctionType& dest) const
1836  {
1837  createMatrixAdapter();
1838  assert( matrixAdap_ );
1839  matrixAdap_->apply( arg.blockVector(), dest.blockVector() );
1840  }
1841 
1843  template <class RowDFType, class ColDFType>
1844  void apply(const RowDFType& arg, ColDFType& dest) const
1845  {
1846  multOEM( arg.leakPointer(), dest.leakPointer ());
1847  }
1848 
1850  double ddotOEM(const double* v, const double* w) const
1851  {
1852  createBlockVectors();
1853 
1854  assert( Arg_ );
1855  assert( Dest_ );
1856 
1857  RowBlockVectorType& V = *Arg_;
1858  ColumnBlockVectorType& W = *Dest_;
1859 
1860  std::copy_n( v, V.size(), V.dbegin() );
1861  std::copy_n( w, W.size(), W.dbegin() );
1862 
1863 #if HAVE_MPI
1864  // in parallel use scalar product of discrete functions
1865  ISTLBlockVectorDiscreteFunction< DomainSpaceType > vF("ddotOEM:vF", domainSpace(), V );
1866  ISTLBlockVectorDiscreteFunction< RangeSpaceType > wF("ddotOEM:wF", rangeSpace(), W );
1867  return vF.scalarProductDofs( wF );
1868 #else
1869  return V * W;
1870 #endif
1871  }
1872 
1874  void resort()
1875  {}
1876 
1878  void createPreconditionMatrix()
1879  {}
1880 
1882  void print(std::ostream & s) const
1883  {
1884  matrix().print(s);
1885  }
1886 
1887  const DomainSpaceType& domainSpace() const { return domainSpace_; }
1888  const RangeSpaceType& rangeSpace() const { return rangeSpace_; }
1889 
1890  const RowMapperType& rowMapper() const { return rowMapper_; }
1891  const ColMapperType& colMapper() const { return colMapper_; }
1892 
1894  ObjectType* newObject() const
1895  {
1896  return new ObjectType(*this, domainSpace(), rangeSpace());
1897  }
1898 
1900  LocalMatrixType localMatrix(const DomainEntityType& domainEntity, const RangeEntityType &rangeEntity ) const
1901  {
1902  return LocalMatrixType( localMatrixStack_, domainEntity, rangeEntity );
1903  }
1904  LocalColumnObjectType localColumn( const DomainEntityType &domainEntity ) const
1905  {
1906  return LocalColumnObjectType ( *this, domainEntity );
1907  }
1908 
1909  template< class LocalMatrix >
1910  void addLocalMatrix ( const DomainEntityType &domainEntity, const RangeEntityType &rangeEntity, const LocalMatrix &localMat )
1911  {
1912  typedef typename MatrixType::size_type Index;
1913  auto blockAccess = [ this ] ( const std::pair< Index, Index > &index ) -> LittleBlockType&
1914  {
1915  if( implicitModeActive() )
1916  return matrix().entry( index.first, index.second );
1917  else
1918  return matrix()[ index.first ][ index.second ];
1919  };
1920 
1921  auto functor = [ &localMat, this, blockAccess ] ( std::pair< int, int > local, const std::pair< Index, Index > &index )
1922  {
1923  LittleBlockType& block = blockAccess( index );
1924  for( std::size_t i = 0; i < littleRows; ++i )
1925  for( std::size_t j = 0; j < littleCols; ++j )
1926  block[ i ][ j ] += localMat.get( local.first * littleRows + i, local.second *littleCols + j );
1927  };
1928 
1929  rowMapper_.mapEach( rangeEntity, makePairFunctor( colMapper_, domainEntity, functor ) );
1930  }
1931 
1932 
1933  template< class LocalMatrix >
1934  void setLocalMatrix ( const DomainEntityType &domainEntity, const RangeEntityType &rangeEntity, const LocalMatrix &localMat )
1935  {
1936  typedef typename MatrixType::size_type Index;
1937  auto blockAccess = [ this ] ( const std::pair< Index, Index > &index ) -> LittleBlockType&
1938  {
1939  if( implicitModeActive() )
1940  return matrix().entry( index.first, index.second );
1941  else
1942  return matrix()[ index.first ][ index.second ];
1943  };
1944 
1945  auto functor = [ &localMat, this, blockAccess ] ( std::pair< int, int > local, const std::pair< Index, Index > &index )
1946  {
1947  LittleBlockType& block = blockAccess( index );
1948  for( std::size_t i = 0; i < littleRows; ++i )
1949  for( std::size_t j = 0; j < littleCols; ++j )
1950  block[ i ][ j ] = localMat.get( local.first * littleRows + i, local.second *littleCols + j );
1951  };
1952 
1953  rowMapper_.mapEach( rangeEntity, makePairFunctor( colMapper_, domainEntity, functor ) );
1954  }
1955 
1956  template< class LocalMatrix, class Scalar >
1957  void addScaledLocalMatrix ( const DomainEntityType &domainEntity, const RangeEntityType &rangeEntity, const LocalMatrix &localMat, const Scalar &s )
1958  {
1959  typedef typename MatrixType::size_type Index;
1960  auto blockAccess = [ this ] ( const std::pair< Index, Index > &index ) -> LittleBlockType&
1961  {
1962  if( implicitModeActive() )
1963  return matrix().entry( index.first, index.second );
1964  else
1965  return matrix()[ index.first ][ index.second ];
1966  };
1967 
1968  auto functor = [ &localMat, &s, this, blockAccess ] ( std::pair< int, int > local, const std::pair< Index, Index > &index )
1969  {
1970  LittleBlockType& block = blockAccess( index );
1971  for( std::size_t i = 0; i < littleRows; ++i )
1972  for( std::size_t j = 0; j < littleCols; ++j )
1973  block[ i ][ j ] += s * localMat.get( local.first * littleRows + i, local.second *littleCols + j );
1974  };
1975 
1976  rowMapper_.mapEach( rangeEntity, makePairFunctor( colMapper_, domainEntity, functor ) );
1977  }
1978 
1979  template< class LocalMatrix >
1980  void getLocalMatrix ( const DomainEntityType &domainEntity, const RangeEntityType &rangeEntity, LocalMatrix &localMat ) const
1981  {
1982  typedef typename MatrixType::size_type Index;
1983  auto functor = [ &localMat, this ] ( std::pair< int, int > local, const std::pair< Index, Index > &global )
1984  {
1985  for( std::size_t i = 0; i < littleRows; ++i )
1986  for( std::size_t j = 0; j < littleCols; ++j )
1987  localMat.set( local.first * littleRows + i, local.second *littleCols + j,
1988  matrix()[ global.first ][ global.second ][i][j] );
1989  };
1990 
1991  rowMapper_.mapEach( rangeEntity, makePairFunctor( colMapper_, domainEntity, functor ) );
1992  }
1993 
1994  protected:
1995  void preConErrorMsg(int preCon) const
1996  {
1997  std::cerr << "ERROR: Wrong precoditioning number (p = " << preCon;
1998  std::cerr <<") in ISTLMatrixObject! \n";
1999  std::cerr <<"Valid values are: \n";
2000  std::cerr <<"0 == no \n";
2001  std::cerr <<"1 == SSOR \n";
2002  std::cerr <<"2 == SOR \n";
2003  std::cerr <<"3 == ILU-0 \n";
2004  std::cerr <<"4 == ILU-n \n";
2005  std::cerr <<"5 == Gauss-Seidel \n";
2006  std::cerr <<"6 == Jacobi \n";
2007  assert(false);
2008  exit(1);
2009  }
2010 
2011  void removeObj( const bool alsoClearMatrix )
2012  {
2013  delete Dest_;
2014  Dest_ = nullptr;
2015  delete Arg_;
2016  Arg_ = nullptr;
2017  delete matrixAdap_;
2018  matrixAdap_ = nullptr;
2019 
2020  if( alsoClearMatrix )
2021  {
2022  delete matrix_;
2023  matrix_ = nullptr;
2024  }
2025  }
2026 
2027  void createBlockVectors() const
2028  {
2029  if( Arg_ == nullptr || Dest_ == nullptr )
2030  {
2031  delete Arg_;
2032  delete Dest_;
2033  Arg_ = new RowBlockVectorType( rowMapper_.size() );
2034  Dest_ = new ColumnBlockVectorType( colMapper_.size() );
2035  }
2036 
2037  createMatrixAdapter ();
2038  }
2039 
2040  void createMatrixAdapter () const
2041  {
2042  if( matrixAdap_ == nullptr )
2043  matrixAdap_ = new MatrixAdapterType(matrixAdapterObject());
2044  }
2045 
2046  };
2047 
2048  } // namespace Fem
2049 
2050 } // namespace Dune
2051 
2052 #endif // #if HAVE_DUNE_ISTL
2053 
2054 #endif // #ifndef DUNE_FEM_ISTLMATRIXWRAPPER_HH
Definition: coordinate.hh:4
static int getEnum(const std::string &key, const std::string(&values)[n])
Definition: io/parameter.hh:388
PairFunctor< Mapper, Entity, Functor > makePairFunctor(const Mapper &mapper, const Entity &entity, Functor functor)
Definition: operator/matrix/functor.hh:65
static T getValue(const std::string &key)
get a mandatory parameter from the container
Definition: io/parameter.hh:332