dune-fem  2.4.1-rc
spmatrix.hh
Go to the documentation of this file.
1 #ifndef DUNE_FEM_SPMATRIX_HH
2 #define DUNE_FEM_SPMATRIX_HH
3 
4 //- system includes
5 #include <vector>
6 #include <set>
7 #include <algorithm>
8 #include <iostream>
9 #include <fstream>
10 
11 //- local includes
13 #include <dune/fem/misc/functor.hh>
17 #include <dune/fem/io/parameter.hh>
23 
25 
26 namespace Dune
27 {
28 
29  namespace Fem
30  {
31 
34  namespace
35  {
41  const int checkNonConstMethods = 0;
42  }
43 
46  : public Dune::Fem::LocalParameter< MatrixParameter, MatrixParameter >
47  {
48 
49  MatrixParameter( const std::string keyPrefix = "" )
50  {}
51 
52  virtual double overflowFraction () const
53  {
54  return 1.0;
55  }
56 
57  virtual int numIterations () const
58  {
59  return 0;
60  }
61 
62  virtual double relaxation () const
63  {
64  return 1.0;
65  }
66 
67  virtual int method () const
68  {
69  return 0;
70  }
71 
72  virtual std::string preconditionName() const
73  {
74  return "None";
75  }
76  };
77 
78 
80  : public MatrixParameter
81  {
83 
84  SparseRowMatrixParameter( const std::string keyPrefix = "spmatrix." )
85  : BaseType( keyPrefix )
86  {}
87  };
88 
89 
90  //*****************************************************************
91  //
92  // --SparseRowMatrix
93  //
100  //*****************************************************************
101 
103  template <class T>
105  {
106  public:
107  typedef T Ttype;
108 
111  typedef ThisType MatrixBaseType;
112 
113  enum { defaultCol = -1 };
114  enum { firstCol = defaultCol + 1 };
115 
116  protected:
117  T* values_ ;
118  int* col_;
119  int* nonZeros_;
120  int dim_[2];
121  int nz_;
122 
123  int memSize_;
124  bool sorted_;
125 
126  // temporary mem for resort
127  std::vector<int> newIndices_;
128  std::vector<T> newValues_;
129 
130  // omega for ssor preconditioning
131  const double omega_;
132 
133  public:
134  SparseRowMatrix(const SparseRowMatrix<T> &S) = delete;
135 
137  SparseRowMatrix(double omega = 1.1);
138 
142  SparseRowMatrix(int rows, int cols, int nz, const T& val = 0,
143  double omega = 1.1);
144 
150  void reserve(int rows, int cols, int nz, const T& dummy);
151 
153  void resize ( int newSize );
154 
156  void resize ( int newRow, int newCol , int newNz = -1 );
157 
159  ~SparseRowMatrix();
160 
161  /*******************************/
162  /* Access and info functions */
163  /*******************************/
165  int numberOfValues () const { return dim_[0]*nz_; }
167  T& val(int i)
168  {
169  assert( i >= 0 );
170  assert( i < numberOfValues () );
171  return values_[i];
172  }
173 
175  const T& val(int i) const {
176  assert( i >= 0 );
177  assert( i < numberOfValues () );
178  return values_[i];
179  }
180 
182  T popValue (int i)
183  {
184  if (checkNonConstMethods) assert(checkConsistency());
185  T v = val(i);
186  values_[i] = 0;
187  col_[i] = -1;
188  if (checkNonConstMethods) assert(checkConsistency());
189  return v;
190  }
191 
193  int realCol (int row, int fakeCol) const
194  {
195  assert(fakeCol<dim_[1]);
196  assert( row < dim_[0] );
197  int pos = row*nz_ + fakeCol;
198  return col_[pos];
199  }
200 
202  std::pair < T , int > realValue(int row, int fakeCol)
203  {
204  assert( row < dim_[0] );
205  assert( fakeCol < nz_ );
206  int pos = row*nz_ + fakeCol;
207  return realValue(pos);
208  }
209 
211  std::pair < const T , int > realValue(int row, int fakeCol) const
212  {
213  assert( row < dim_[0] );
214  assert( fakeCol < nz_ );
215  int pos = row*nz_ + fakeCol;
216  return realValue(pos);
217  }
218 
220  std::pair < T , int > realValue(int index)
221  {
222  return std::pair< T, int > (values_[index], col_[index]);
223  }
224 
226  std::pair < const T , int > realValue(int index) const
227  {
228  return std::pair< const T, int > (values_[index], col_[index]);
229  }
230 
232  int colIndex(int row, int col);
233 
235  bool find (int row, int col) const;
236 
238  int dim(int i) const {return dim_[i];}
240  int size(int i) const {return dim_[i];}
241 
243  int rows() const {return dim_[0];}
244 
246  int cols() const {return dim_[1];}
247 
249  int numNonZeros() const {return nz_;}
250 
252  int numNonZeros(int i) const
253  {
254  assert( nonZeros_ );
255  return nonZeros_[i];
256  }
257 
259  T operator() ( const int row, const int col ) const;
260  T operator() ( const unsigned int row, const unsigned int col ) const;
261  T operator() ( const long unsigned int row, const long unsigned int col ) const
262  {
263  return this->operator()((unsigned int)(row), (unsigned int)(col) );
264  }
265 
270  void set(int row, int col, T val);
271 
273  void clearRow (int row);
274 
276  void clearCol ( int col );
277 
279  void scaleRow (int row, const T& val );
280 
283  void clear();
284 
286  void add(int row, int col, T val);
287 
289  void multScalar(int row, int col, T val);
290 
292  void kroneckerKill(int row, int col);
293 
295  template <class VECtype>
296  void mult(const VECtype *x, VECtype * ret) const;
297 
299  template <class VECtype>
300  T multOEMRow (const VECtype *x , const int row ) const;
301 
303  template <class VECtype>
304  void multOEM(const VECtype *x, VECtype * ret) const;
305 
307  template <class VECtype>
308  void multOEMAdd(const VECtype *x, VECtype * ret) const;
309 
311  template <class VECtype>
312  void multOEM_t(const VECtype *x, VECtype * ret) const;
313 
315  template <class DiscFType, class DiscFuncType>
316  void apply(const DiscFType &f, DiscFuncType &ret) const;
317 
319  template <class ArgDFType, class DestDFType>
320  void apply_t(const ArgDFType& f, DestDFType &ret) const;
321 
323  template <class DiscFuncType>
324  void operator () (const DiscFuncType &f, DiscFuncType &ret) const
325  {
326  apply(f,ret);
327  }
328 
330  template <class DiscFuncType>
331  void getDiag(const ThisType&A, const ThisType &B, DiscFuncType &rhs) const;
332 
334  template <class DiscFuncType>
335  void getDiag(const ThisType &A, DiscFuncType &rhs) const;
336 
338  template <class DiscFuncType>
339  void getDiag(DiscFuncType &rhs) const;
340 
342  template <class DiscFuncType>
343  void addDiag(DiscFuncType &rhs) const;
344 
346  void print (std::ostream& s, unsigned int offset=0) const;
347 
349  void printReal (std::ostream& s) const;
350 
352  void printColumns(std::ostream& s) const;
353 
358  void printNonZeros(std::ostream& s) const;
359 
369  bool checkConsistency() const;
370 
372  void unitRow(int row);
373 
375  void unitCol(int col);
376 
378  void checkSym ();
379 
380  // res = this * B
381  void multiply(const ThisType & B, ThisType & res) const;
382 
384  void scale(const T& factor);
385 
387  void add(const ThisType & B );
388 
390  void resort();
391 
393  void resortRow(const int row);
394 
396  void ssorPrecondition (const T*, T*) const;
397 
399  bool rightPrecondition() const { return true; }
400 
402  void precondition (const T*u , T*x) const { ssorPrecondition(u,x); }
403 
404  void DUNE_DEPRECATED_MSG("Use directly Dune::Fem::UMFPACKOp instead of solveUMF") solveUMF(const T* b, T* x);
405  void DUNE_DEPRECATED_MSG("Use directly Dune::Fem::UMFPACKOp instead of solveUMFNonSymmetric") solveUMFNonSymmetric(const T* b, T* x);
406 
407  private:
409  void removeObj();
410  };
411 
412 
413  // SparseRowMatrixObject
414  // ---------------------
415  template <class DomainSpace, class RangeSpace,
418 
419  template <class RowSpaceImp, class ColSpaceImp = RowSpaceImp>
421  {
422  typedef RowSpaceImp RowSpaceType;
423  typedef ColSpaceImp ColumnSpaceType;
426 
427  typedef RowSpaceImp RangeSpaceType;
428  typedef ColSpaceImp DomainSpaceType;
429 
430  };
431 
432  template< class DomainSpace, class RangeSpace, class Matrix >
434  : public Fem :: OEMMatrix
435  {
436  public:
437  typedef DomainSpace DomainSpaceType;
438  typedef RangeSpace RangeSpaceType;
439 
440  /*******************************************************************
441  * Rows belong to the DomainSpace and Columns to the RangeSpace *
442  *******************************************************************/
443  typedef typename DomainSpaceType :: EntityType DomainEntityType ;
444  typedef typename RangeSpaceType :: EntityType RangeEntityType ;
445  typedef typename DomainSpaceType :: EntityType ColumnEntityType ;
446  typedef typename RangeSpaceType :: EntityType RowEntityType ;
447 
448  typedef typename DomainSpaceType :: BlockMapperType DomainBlockMapperType ;
449  typedef NonBlockMapper< DomainBlockMapperType,
450  DomainSpaceType :: localBlockSize > DomainMapperType;
451  typedef typename RangeSpaceType :: BlockMapperType RangeBlockMapperType ;
452  typedef NonBlockMapper< RangeBlockMapperType,
453  RangeSpaceType :: localBlockSize > RangeMapperType;
454 
455  typedef Matrix MatrixType ;
456  private:
458 
459  protected:
460  typedef typename DomainSpaceType :: GridType GridType;
461 
462 
463  template< class MatrixObject >
465 
466  template< class MatrixObject >
467  class LocalMatrix;
468 
469  public:
470  typedef MatrixType PreconditionMatrixType;
471 
473  typedef LocalMatrix<ThisType> ObjectType;
474  typedef ThisType LocalMatrixFactoryType;
478 
480 
481  protected:
482  const DomainSpaceType &domainSpace_;
483  const RangeSpaceType &rangeSpace_;
484 
485  DomainMapperType domainMapper_ ;
486  RangeMapperType rangeMapper_ ;
487 
489 
490  mutable MatrixType matrix_;
492 
493  mutable LocalMatrixStackType localMatrixStack_;
494 
495  public:
496 
498  inline SparseRowMatrixObject( const DomainSpaceType &domainSpace,
499  const RangeSpaceType &rangeSpace,
500  const std::string& paramfile )
501  DUNE_DEPRECATED_MSG("SparseRowMatrixObject(...,string) is deprecated. Use SparseRowMatrixObject(string,DomainSpace,RangeSpace,MatrixParameter) instead")
502  : domainSpace_( domainSpace ),
503  rangeSpace_( rangeSpace ),
504  domainMapper_( domainSpace_.blockMapper() ),
505  rangeMapper_( rangeSpace_.blockMapper() ),
506  sequence_( -1 ),
507  matrix_(),
508  preconditioning_( false ),
509  localMatrixStack_( *this )
510  {
511  int precon = 0;
512  if( paramfile != "" )
513  readParameter( paramfile, "Preconditioning", precon );
514  else
515  precon = Parameter :: getValue("Preconditioning", precon );
516  preconditioning_ = (precon > 0) ? true : false;
517  }
518 
519 
521  inline SparseRowMatrixObject( const DomainSpaceType &domainSpace,
522  const RangeSpaceType &rangeSpace,
523  const MatrixParameter& param = SparseRowMatrixParameter() )
524  : domainSpace_( domainSpace ),
525  rangeSpace_( rangeSpace ),
526  domainMapper_( domainSpace_.blockMapper() ),
527  rangeMapper_( rangeSpace_.blockMapper() ),
528  sequence_( -1 ),
529  matrix_(),
530  preconditioning_( param.method() != 0 ),
531  localMatrixStack_( *this )
532  {}
533 
535  const DomainSpaceType& domainSpace() const { return domainSpace_; }
536 
538  const RangeSpaceType& rangeSpace() const { return rangeSpace_; }
539 
541  inline MatrixType &matrix () const
542  {
543  return matrix_;
544  }
545 
547  inline ObjectType *newObject () const
548  {
549  return new ObjectType( *this, domainSpace_, rangeSpace_, domainMapper_, rangeMapper_ );
550  }
551 
553  inline LocalMatrixType localMatrix( const DomainEntityType &domainEntity,
554  const RangeEntityType &rangeEntity ) const
555  {
556  /*******************************************************************
557  * Rows belong to the DomainSpace and Columns to the RangeSpace *
558  *******************************************************************/
559  return LocalMatrixType( localMatrixStack_, domainEntity, rangeEntity );
560  }
561 
562  LocalColumnObjectType localColumn( const DomainEntityType &domainEntity ) const
563  {
564  return LocalColumnObjectType ( *this, domainEntity );
565  }
566 
567  template< class LocalMatrix >
568  void addLocalMatrix ( const DomainEntityType &domainEntity, const RangeEntityType &rangeEntity, const LocalMatrix &localMat )
569  {
570  typedef int Index;
571  auto functor = [ &localMat, this ] ( std::pair< int, int > local, const std::pair< Index, Index >& global )
572  {
573  matrix_.add( global.first, global.second, localMat.get( local.first, local.second ) );
574  };
575 
576  rangeMapper_.mapEach( rangeEntity, makePairFunctor( domainMapper_, domainEntity, functor ) );
577  }
578 
579 
580  template< class LocalMatrix, class Scalar >
581  void addScaledLocalMatrix ( const DomainEntityType &domainEntity, const RangeEntityType &rangeEntity, const LocalMatrix &localMat, const Scalar &s )
582  {
583  typedef int Index;
584  auto functor = [ &localMat, &s, this ] ( std::pair< int, int > local, const std::pair< Index, Index >& global )
585  {
586  matrix_.add( global.first, global.second, s * localMat.get( local.first, local.second ) );
587  };
588 
589  rangeMapper_.mapEach( rangeEntity, makePairFunctor( domainMapper_, domainEntity, functor ) );
590  }
591 
592 
593  template< class LocalMatrix >
594  void setLocalMatrix ( const DomainEntityType &domainEntity, const RangeEntityType &rangeEntity, const LocalMatrix &localMat )
595  {
596  typedef int Index;
597  auto functor = [ &localMat, this ] ( std::pair< int, int > local, const std::pair< Index, Index >& global )
598  {
599  matrix_.set( global.first, global.second, localMat.get( local.first, local.second ) );
600  };
601 
602  rangeMapper_.mapEach( rangeEntity, makePairFunctor( domainMapper_, domainEntity, functor ) );
603  }
604 
605 
606  template< class LocalMatrix >
607  void getLocalMatrix ( const DomainEntityType &domainEntity, const RangeEntityType &rangeEntity, LocalMatrix &localMat ) const
608  {
609  typedef int Index;
610  auto functor = [ &localMat, this ] ( std::pair< int, int > local, const std::pair< Index, Index >& global )
611  {
612  localMat.set( local.first, local.second, matrix_( global.first, global.second ) );
613  };
614 
615  rangeMapper_.mapEach( rangeEntity, makePairFunctor( domainMapper_, domainEntity, functor ) );
616  }
617 
619  inline void clear ()
620  {
621  matrix_.clear();
622  }
623 
625  bool hasPreconditionMatrix () const
626  {
627  return preconditioning_;
628  }
629 
631  const PreconditionMatrixType &preconditionMatrix () const
632  {
633  return matrix_;
634  }
635 
637  template <class Stencil>
638  inline void reserve(const Stencil &stencil, bool verbose = false )
639  {
640  if( sequence_ != domainSpace_.sequence() )
641  {
642 #ifndef DUNE_FEM_DONT_CHECKITERATORS_OF_SPACE
643  // if empty grid do nothing (can appear in parallel runs)
644  if( (domainSpace_.begin() != domainSpace_.end())
645  && (rangeSpace_.begin() != rangeSpace_.end()) )
646 #endif
647  {
648 
649  if( verbose )
650  {
651  const int rowMaxNumbers = rangeMapper_.maxNumDofs();
652  const int colMaxNumbers = domainMapper_.maxNumDofs();
653 
654  std::cout << "Reserve Matrix with (" << rangeSpace_.size() << "," << domainSpace_.size()<< ")\n";
655  std::cout << "Max number of base functions = (" << rowMaxNumbers << "," << colMaxNumbers << ")\n";
656  }
657 
658  // upper estimate for number of non-zeros
659  const static size_t domainLocalBlockSize = DomainSpaceType::localBlockSize;
660  const int nonZeros = std::max( (int)(stencil.maxNonZerosEstimate()*domainLocalBlockSize),
661  matrix_.numNonZeros() );
662  matrix_.reserve( rangeSpace_.size(), domainSpace_.size(), nonZeros, 0.0 );
663  }
664  sequence_ = domainSpace_.sequence();
665  }
666  }
667 
668  template< class DomainFunction, class RangeFunction >
669  void DUNE_DEPRECATED_MSG("Use directly Dune::Fem::UMFPACKOp instead of solveUMF")
670  solveUMF( const DomainFunction &arg, RangeFunction &dest ) const
671  {
672  matrix_.solveUMF( arg.leakPointer(), dest.leakPointer() );
673  }
674 
675  template< class DomainFunction, class RangeFunction >
676  void DUNE_DEPRECATED_MSG("Use directly Dune::Fem::UMFPACKOp instead of solveUMFNonSymmetric")
677  solveUMFNonSymmetric( const DomainFunction &arg, RangeFunction &dest ) const
678  {
679  matrix_.solveUMFNonSymmetric( arg.leakPointer(), dest.leakPointer() );
680  }
681 
683  template< class DomainFunction, class RangeFunction >
684  void apply ( const DomainFunction &arg, RangeFunction &dest ) const
685  {
686  // do matrix vector multiplication
687  matrix_.apply( arg, dest );
688 
689  // communicate data
690  dest.communicate();
691  }
692 
696  {
697  // do matrix vector multiplication
698  matrix_.multOEM( arg.leakPointer(), dest.leakPointer() );
699 
700  // communicate data
701  dest.communicate();
702  }
703 
705  template< class DomainFunction, class RangeFunction >
706  void apply_t ( const RangeFunction &arg, DomainFunction &dest ) const
707  {
708  // do matrix vector multiplication
709  matrix_.apply_t( arg, dest );
710 
711  // communicate data
712  dest.communicate();
713  }
714 
718  {
719  // do matrix vector multiplication
720  matrix_.multOEM_t( arg.leakPointer(), dest.leakPointer() );
721 
722  // communicate data
723  dest.communicate();
724  }
725 
726 
728  double ddotOEM( const double *v, const double *w ) const
729  {
730  typedef AdaptiveDiscreteFunction< DomainSpaceType > DomainFunctionType;
731  DomainFunctionType V( "ddot V", domainSpace_, v );
732  DomainFunctionType W( "ddot W", domainSpace_, w );
733  return V.scalarProductDofs( W );
734  }
735 
737  void multOEM( const double *arg, double *dest ) const
738  {
739  typedef AdaptiveDiscreteFunction< DomainSpaceType > DomainFunctionType;
740  typedef AdaptiveDiscreteFunction< RangeSpaceType > RangeFunctionType;
741 
742  DomainFunctionType farg( "multOEM arg", domainSpace_, arg );
743  RangeFunctionType fdest( "multOEM dest", rangeSpace_, dest );
744  apply( farg, fdest );
745  }
746 
748  void resort()
749  {
750  matrix_.resort();
751  }
752 
754  {
755  }
756 
758  template < class DiscreteFunctionType >
759  void extractDiagonal( DiscreteFunctionType& diag ) const
760  {
761  // this only works for matrices with same number of rows,cols
762  assert( matrix_.rows() == matrix_.cols() );
763  typedef typename DiscreteFunctionType :: DofIteratorType DofIteratorType ;
764  const DofIteratorType dofEnd = diag.dend();
765  unsigned int row = 0;
766  for( DofIteratorType dofIt = diag.dbegin();
767  dofIt != dofEnd; ++ dofIt, ++row )
768  {
769  assert( row < ( unsigned int )matrix_.rows() );
770  (*dofIt) = matrix_( row, row );
771  }
772  }
773 
775  template <class HangingNodesType>
776  void changeHangingNodes(const HangingNodesType& hangingNodes)
777  {
778  {
779  typedef typename HangingNodesType :: IteratorType IteratorType;
780  const IteratorType end = hangingNodes.end();
781  for( IteratorType it = hangingNodes.begin(); it != end; ++it)
782  {
783  insertHangingRow( hangingNodes, (*it).first , (*it).second );
784  }
785  }
786 
787  /*
788  {
789  typedef typename HangingNodesType :: SlaveNodesType SlaveNodesType;
790  typedef typename SlaveNodesType :: const_iterator iterator;
791  iterator end = hangingNodes.slaveNodes().end();
792  for( iterator it = hangingNodes.slaveNodes().begin();
793  it != end; ++it )
794  {
795  matrix().unitRow( *it );
796  matrix().set( *it, *it, 0.0 );
797  }
798  }
799  */
800  }
801  protected:
803  template <class HangingNodesType, class ColumnVectorType>
804  void insertHangingRow( const HangingNodesType& hangingNodes,
805  const int row, const ColumnVectorType& colVec)
806  {
807  const size_t cols = colVec.size();
808 
809  // distribute row to associated rows
810  const int nonZeros = matrix().numNonZeros( row );
811  for( int c = 0; c < nonZeros; ++c)
812  {
813  std::pair< double, int > val = matrix().realValue( row, c );
814  for( size_t j = 0; j < cols; ++ j)
815  {
816  const double value = colVec[j].second * val.first;
817  assert( ! hangingNodes.isHangingNode( colVec[j].first ) );
818  matrix().add( colVec[j].first , val.second, value );
819  }
820  }
821 
822  // replace hanging row
823  matrix().unitRow( row );
824  for( size_t j = 0; j < cols; ++ j)
825  {
826  matrix().add( row, colVec[j].first , -colVec[j].second );
827  }
828  }
829  };
830 
831 
832 
833  template< class DomainSpace, class RangeSpace, class Matrix >
834  template< class MatrixObject >
835  struct SparseRowMatrixObject< DomainSpace, RangeSpace, Matrix >::LocalMatrixTraits
836  {
837  typedef DomainSpace DomainSpaceType;
838  typedef RangeSpace RangeSpaceType;
839 
842 
843  typedef typename SparseRowMatrixObjectType :: template LocalMatrix< MatrixObject > LocalMatrixType;
844 
845  typedef typename RangeSpaceType :: RangeFieldType RangeFieldType;
846  typedef RangeFieldType LittleBlockType;
847 
850  };
851 
852 
853 
855  template< class DomainSpace, class RangeSpace, class Matrix >
856  template< class MatrixObject >
857  class SparseRowMatrixObject< DomainSpace, RangeSpace, Matrix > :: LocalMatrix
858  : public LocalMatrixDefault< LocalMatrixTraits< MatrixObject > >
859  {
860  public:
862  typedef MatrixObject MatrixObjectType;
863 
866 
867  private:
869 
870  public:
872  typedef typename MatrixObjectType :: MatrixType MatrixType;
873 
876 
878  typedef RangeFieldType DofType;
879 
882 
887 
888  protected:
889  MatrixType &matrix_;
890  const DomainMapperType& domainMapper_;
891  const RangeMapperType& rangeMapper_;
892 
893  typedef std :: vector< typename RangeMapperType :: SizeType > RowIndicesType ;
895  RowIndicesType rowIndices_;
896 
897  typedef std :: vector< typename DomainMapperType :: SizeType > ColumnIndicesType ;
899  ColumnIndicesType columnIndices_;
900 
901  using BaseType :: domainSpace_;
902  using BaseType :: rangeSpace_;
903 
904  public:
907  inline LocalMatrix( const MatrixObjectType &matrixObject,
908  const DomainSpaceType &domainSpace,
909  const RangeSpaceType &rangeSpace,
910  const DomainMapperType& domainMapper,
911  const RangeMapperType& rangeMapper )
912  : BaseType( domainSpace, rangeSpace),
913  matrix_( matrixObject.matrix() ),
914  domainMapper_( domainMapper ),
915  rangeMapper_( rangeMapper )
916  {
917  }
918 
919  LocalMatrix( const LocalMatrix & ) = delete;
920 
921  void init( const DomainEntityType &domainEntity, const RangeEntityType &rangeEntity )
922  {
923  /*******************************************************************
924  * Rows belong to the DomainSpace and Columns to the RangeSpace *
925  *******************************************************************/
926 
927  // initialize base functions sets
928  BaseType::init( domainEntity, rangeEntity );
929 
930  // rows are determined by the range space
931  rowIndices_.resize( rangeMapper_.numDofs( rangeEntity ) );
932  rangeMapper_.mapEach( rangeEntity, AssignFunctor< RowIndicesType >( rowIndices_ ) );
933 
934  // columns are determind by the domain space
935  columnIndices_.resize( domainMapper_.numDofs( domainEntity ) );
936  domainMapper_.mapEach( domainEntity, AssignFunctor< ColumnIndicesType >( columnIndices_ ) );
937  }
938 
940  int rows () const
941  {
942  return rowIndices_.size();
943  }
944 
946  int columns () const
947  {
948  return columnIndices_.size();
949  }
950 
952  void add( int localRow, int localCol, const DofType value )
953  {
954  assert( value == value );
955  assert( (localRow >= 0) && (localRow < rows()) );
956  assert( (localCol >= 0) && (localCol < columns()) );
957 
958  matrix_.add( rowIndices_[ localRow ], columnIndices_[ localCol ], value );
959  }
960 
962  DofType get( int localRow, int localCol ) const
963  {
964  assert( (localRow >= 0) && (localRow < rows()) );
965  assert( (localCol >= 0) && (localCol < columns()) );
966 
967  return matrix_( rowIndices_[ localRow ], columnIndices_[ localCol ] );
968  }
969 
971  void set( int localRow, int localCol, const DofType value )
972  {
973  assert( (localRow >= 0) && (localRow < rows()) );
974  assert( (localCol >= 0) && (localCol < columns()) );
975 
976  matrix_.set( rowIndices_[ localRow ], columnIndices_[ localCol ], value );
977  }
978 
980  void unitRow( const int localRow )
981  {
982  assert( (localRow >= 0) && (localRow < rows()) );
983  matrix_.unitRow( rowIndices_[ localRow ] );
984  }
985 
987  void clearRow( const int localRow )
988  {
989  assert( (localRow >= 0) && (localRow < rows()) );
990  matrix_.clearRow( rowIndices_[localRow]);
991  }
992 
994  void clearCol ( const int localCol )
995  {
996  assert( (localCol >= 0) && (localCol < columns()) );
997  matrix_.clearCol( columnIndices_[localCol] );
998  }
999 
1001  void clear ()
1002  {
1003  const int row = rows();
1004  for( int i = 0; i < row; ++i )
1005  matrix_.clearRow( rowIndices_[ i ] );
1006  }
1007 
1009  void scale ( const DofType& value )
1010  {
1011  const int row = rows();
1012  for( int i = 0; i < row; ++i )
1013  matrix_.scaleRow( rowIndices_[ i ] , value );
1014  }
1015 
1017  void resort ()
1018  {
1019  const int row = rows();
1020  for( int i = 0; i < row; ++i )
1021  matrix_.resortRow( rowIndices_[ i ] );
1022  }
1023  };
1024 
1025  } // namespace Fem
1026 
1027 } // namespace Dune
1028 
1029 #include "spmatrix.cc"
1030 
1031 #endif // #ifndef DUNE_FEM_SPMATRIX_HH
Traits::RangeMapperType RangeMapperType
type of nonblocked domain mapper
Definition: spmatrix.hh:886
int memSize_
number of nonzeros per row
Definition: spmatrix.hh:123
Traits::DomainSpaceType DomainSpaceType
type of domain discrete function space
Definition: localmatrix.hh:51
void apply(const DomainFunction &arg, RangeFunction &dest) const
apply matrix to discrete function
Definition: spmatrix.hh:684
T popValue(int i)
return value and clear matrix entry
Definition: spmatrix.hh:182
void insertHangingRow(const HangingNodesType &hangingNodes, const int row, const ColumnVectorType &colVec)
insert row to be a row for a hanging node
Definition: spmatrix.hh:804
void multOEM(const double *arg, double *dest) const
mult method of matrix object used by oem solver
Definition: spmatrix.hh:737
MatrixType & matrix() const
return reference to stability matrix
Definition: spmatrix.hh:541
std::vector< int > newIndices_
Definition: spmatrix.hh:127
const DomainSpaceType & domainSpace_
Definition: spmatrix.hh:482
const double omega_
Definition: spmatrix.hh:131
bool sorted_
Definition: spmatrix.hh:124
void getLocalMatrix(const DomainEntityType &domainEntity, const RangeEntityType &rangeEntity, LocalMatrix &localMat) const
Definition: spmatrix.hh:607
int cols() const
return number of columns
Definition: spmatrix.hh:246
LocalMatrixStackType localMatrixStack_
Definition: spmatrix.hh:493
Traits::RangeFieldType RangeFieldType
type of entries of little blocks
Definition: spmatrix.hh:875
void clearRow(const int localRow)
set matrix row to zero
Definition: spmatrix.hh:987
RowIndicesType rowIndices_
global index in the DomainSpace
Definition: spmatrix.hh:895
Definition: localmatrixwrapper.hh:16
const T & val(int i) const
return reference to value on given entry
Definition: spmatrix.hh:175
int * col_
data values (nz_ * dim_[0] elements)
Definition: spmatrix.hh:118
RangeSpaceType::EntityType RangeEntityType
Definition: localmatrix.hh:65
MatrixType matrix_
Definition: spmatrix.hh:490
DomainSpaceType::EntityType ColumnEntityType
Definition: spmatrix.hh:445
void apply(const AdaptiveDiscreteFunction< DomainSpaceType > &arg, AdaptiveDiscreteFunction< RangeSpaceType > &dest) const
apply matrix to discrete function
Definition: spmatrix.hh:694
void unitRow(const int localRow)
set matrix row to zero except diagonla entry
Definition: spmatrix.hh:980
T * values_
Definition: spmatrix.hh:117
Matrix MatrixType
Definition: spmatrix.hh:455
RangeSpaceType::RangeFieldType RangeFieldType
Definition: spmatrix.hh:845
std::pair< const T, int > realValue(int index) const
return pair< value, column >, used by BlockMatrix
Definition: spmatrix.hh:226
void add(int localRow, int localCol, const DofType value)
add value to matrix entry
Definition: spmatrix.hh:952
int sequence_
Definition: spmatrix.hh:488
default implementation for a general operator stencil
Definition: stencil.hh:31
Fem::ObjectStack< LocalMatrixFactoryType > LocalMatrixStackType
Definition: spmatrix.hh:475
static constexpr T max(T a)
Definition: utility.hh:65
MatrixType PreconditionMatrixType
Definition: spmatrix.hh:467
const RangeSpaceType & rangeSpace_
Definition: spmatrix.hh:483
virtual int numIterations() const
Definition: spmatrix.hh:57
Definition: adaptivefunction/adaptivefunction.hh:34
void multiply(const FieldMatrix< Field1, m, n > &A, const FieldVector< Field2, n > &x, FieldVector< Field3, m > &y)
Definition: fieldmatrixhelper.hh:18
NonBlockMapper< RangeBlockMapperType, RangeSpaceType::localBlockSize > RangeMapperType
Definition: spmatrix.hh:453
void precondition(const T *u, T *x) const
apply preconditioning, calls ssorPreconditioning at the moment
Definition: spmatrix.hh:402
SparseRowMatrixObject< RowSpaceType, ColumnSpaceType > MatrixObjectType
Definition: spmatrix.hh:425
T & val(int i)
matrix value taken from real array
Definition: spmatrix.hh:167
const DomainMapperType & domainMapper_
Definition: spmatrix.hh:890
double ddotOEM(const double *v, const double *w) const
mult method of matrix object used by oem solver
Definition: spmatrix.hh:728
LocalColumnObjectType localColumn(const DomainEntityType &domainEntity) const
Definition: spmatrix.hh:562
const RangeMapperType & rangeMapper_
Definition: spmatrix.hh:891
Definition: spmatrix.hh:417
static bool readParameter(std::istream &file, const std::string keyword, T &data, bool verbose=true, bool warn=true)
Definition: asciiparser.hh:18
Default implementation for local matrix classes.
Definition: localmatrix.hh:266
void addScaledLocalMatrix(const DomainEntityType &domainEntity, const RangeEntityType &rangeEntity, const LocalMatrix &localMat, const Scalar &s)
Definition: spmatrix.hh:581
const DomainSpaceType & domainSpace() const
return domain space (i.e. space that builds the rows)
Definition: spmatrix.hh:535
RangeSpaceType::BlockMapperType RangeBlockMapperType
Definition: spmatrix.hh:451
LocalMatrix.
Definition: spmatrix.hh:467
bool rightPrecondition() const
returns true if preconditioing is called before matrix multiply
Definition: spmatrix.hh:399
std::pair< T, int > realValue(int row, int fakeCol)
return pair< value, column >, used by BlockMatrix
Definition: spmatrix.hh:202
int maxNonZerosEstimate() const
Return an upper bound for the maximum number of non-zero entries in all row.
Definition: stencil.hh:99
SparseRowMatrixObject(const DomainSpaceType &domainSpace, const RangeSpaceType &rangeSpace, const MatrixParameter &param=SparseRowMatrixParameter())
construct matrix object
Definition: spmatrix.hh:521
std::vector< T > newValues_
Definition: spmatrix.hh:128
int numNonZeros() const
return max number of non zeros
Definition: spmatrix.hh:249
const RangeSpaceType & rangeSpace() const
return range space (i.e. space that builds the columns)
Definition: spmatrix.hh:538
void changeHangingNodes(const HangingNodesType &hangingNodes)
delete all row belonging to a hanging node and rebuild them
Definition: spmatrix.hh:776
SparseRowMatrixTraits< RowSpaceType, ColumnSpaceType > ThisType
Definition: spmatrix.hh:424
ColumnIndicesType columnIndices_
global index in the RangeSpace
Definition: spmatrix.hh:899
void apply_t(const RangeFunction &arg, DomainFunction &dest) const
apply transposed matrix to discrete function
Definition: spmatrix.hh:706
ThisType MatrixBaseType
type of the base matrix (for consistency with ISTLMatrixObject)
Definition: spmatrix.hh:111
SparseRowMatrixObjectType::DomainMapperType DomainMapperType
Definition: spmatrix.hh:848
LocalMatrixTraits< MatrixObjectType > Traits
type of the traits
Definition: spmatrix.hh:865
Definition: misc/functor.hh:30
ThisType LocalMatrixFactoryType
Definition: spmatrix.hh:474
Definition: nonblockmapper.hh:19
bool hasPreconditionMatrix() const
return true if precoditioning matrix is provided
Definition: spmatrix.hh:625
ColSpaceImp ColumnSpaceType
Definition: spmatrix.hh:423
SparseRowMatrixParameter(const std::string keyPrefix="spmatrix.")
Definition: spmatrix.hh:84
int realCol(int row, int fakeCol) const
return real column number for (row,localCol)
Definition: spmatrix.hh:193
RangeSpaceType::EntityType RowEntityType
Definition: spmatrix.hh:446
void resort()
resort all global rows of matrix to have ascending numbering
Definition: spmatrix.hh:1017
int size(int i) const
return number of rows = 0, cols = 1
Definition: spmatrix.hh:240
MatrixObject MatrixObjectType
type of matrix object
Definition: spmatrix.hh:862
RangeSpaceType::EntityType RangeEntityType
Definition: spmatrix.hh:444
ColumnObject< ThisType > LocalColumnObjectType
Definition: spmatrix.hh:479
std::vector< typename DomainMapperType::SizeType > ColumnIndicesType
Definition: spmatrix.hh:897
DomainSpace DomainSpaceType
Definition: spmatrix.hh:437
int numNonZeros(int i) const
return number of non zeros in row
Definition: spmatrix.hh:252
SparseRowMatrixObjectType::template LocalMatrix< MatrixObject > LocalMatrixType
Definition: spmatrix.hh:843
Definition: io/parameter.hh:549
void clear()
resize all matrices and clear them
Definition: spmatrix.hh:619
MatrixType & matrix_
Definition: spmatrix.hh:889
Definition: coordinate.hh:4
SparseRowMatrix.
Definition: spmatrix.hh:104
LocalMatrix< ThisType > ObjectType
type of local matrix
Definition: spmatrix.hh:473
virtual std::string preconditionName() const
Definition: spmatrix.hh:72
std::vector< typename RangeMapperType::SizeType > RowIndicesType
Definition: spmatrix.hh:893
void communicate()
do default communication of space for this discrete function
Definition: discretefunction.hh:831
RangeMapperType rangeMapper_
Definition: spmatrix.hh:486
void resort()
resort row numbering in matrix to have ascending numbering
Definition: spmatrix.hh:748
Traits::LittleBlockType LittleBlockType
type of little blocks
Definition: spmatrix.hh:881
int rows() const
return number of rows
Definition: spmatrix.hh:243
const PreconditionMatrixType & preconditionMatrix() const
return reference to preconditioner
Definition: spmatrix.hh:631
void clearCol(const int localCol)
set matrix column to zero
Definition: spmatrix.hh:994
ColSpaceImp DomainSpaceType
Definition: spmatrix.hh:428
Traits::DomainMapperType DomainMapperType
type of nonblocked domain mapper
Definition: spmatrix.hh:884
void apply_t(const AdaptiveDiscreteFunction< RangeSpaceType > &arg, AdaptiveDiscreteFunction< DomainSpaceType > &dest) const
apply transposed matrix to discrete function
Definition: spmatrix.hh:716
DomainSpace DomainSpaceType
Definition: spmatrix.hh:837
MatrixParameter.
Definition: spmatrix.hh:45
DomainSpaceType::GridType GridType
Definition: spmatrix.hh:460
DofType * leakPointer()
Definition: adaptivefunction/adaptivefunction.hh:95
void scale(const DofType &value)
scale local matrix with a certain value
Definition: spmatrix.hh:1009
Traits::RangeSpaceType RangeSpaceType
type of range discrete function space
Definition: localmatrix.hh:54
MatrixParameter(const std::string keyPrefix="")
Definition: spmatrix.hh:49
NonBlockMapper< DomainBlockMapperType, DomainSpaceType::localBlockSize > DomainMapperType
Definition: spmatrix.hh:450
ObjectType * newObject() const
interface method from LocalMatrixFactory
Definition: spmatrix.hh:547
void setLocalMatrix(const DomainEntityType &domainEntity, const RangeEntityType &rangeEntity, const LocalMatrix &localMat)
Definition: spmatrix.hh:594
int rows() const
return number of rows
Definition: spmatrix.hh:940
DomainSpaceType::BlockMapperType DomainBlockMapperType
Definition: spmatrix.hh:448
PairFunctor< Mapper, Entity, Functor > makePairFunctor(const Mapper &mapper, const Entity &entity, Functor functor)
Definition: operator/matrix/functor.hh:65
std::pair< T, int > realValue(int index)
return pair< value, column >, used by BlockMatrix
Definition: spmatrix.hh:220
void extractDiagonal(DiscreteFunctionType &diag) const
extract diagonal entries from matrix into discrete function
Definition: spmatrix.hh:759
virtual double relaxation() const
Definition: spmatrix.hh:62
interface for matrices to be used with OEM sovlers
Definition: oemsolver/oemsolver.hh:31
Definition: spmatrix.hh:420
SparseRowMatrix< T > ThisType
remember the value type
Definition: spmatrix.hh:109
int nz_
dim_[0] x dim_[1] Matrix
Definition: spmatrix.hh:121
T Ttype
Definition: spmatrix.hh:107
RowSpaceImp RangeSpaceType
Definition: spmatrix.hh:427
DomainSpaceType::EntityType DomainEntityType
Definition: spmatrix.hh:443
RangeSpace RangeSpaceType
Definition: spmatrix.hh:438
int * nonZeros_
row_ptr (dim_[0]+1 elements)
Definition: spmatrix.hh:119
bool preconditioning_
Definition: spmatrix.hh:491
void clear()
clear all entries belonging to local matrix
Definition: spmatrix.hh:1001
DomainMapperType domainMapper_
Definition: spmatrix.hh:485
SparseRowMatrixObject(const DomainSpaceType &domainSpace, const RangeSpaceType &rangeSpace, const std::string &paramfile)
construct matrix object
Definition: spmatrix.hh:498
int numberOfValues() const
length of array used to store matrix
Definition: spmatrix.hh:165
RangeFieldType DofType
type of the DoFs
Definition: spmatrix.hh:878
static T getValue(const std::string &key)
get a mandatory parameter from the container
Definition: io/parameter.hh:332
MatrixParameter BaseType
Definition: spmatrix.hh:82
LocalMatrixType localMatrix(const DomainEntityType &domainEntity, const RangeEntityType &rangeEntity) const
return local matrix
Definition: spmatrix.hh:553
void mult(const MatrixImp &m, const VectorType *x, VectorType *ret)
Definition: oemsolver/oemsolver.hh:165
Definition: columnobject.hh:11
int columns() const
return number of columns
Definition: spmatrix.hh:946
RangeFieldType LittleBlockType
Definition: spmatrix.hh:846
MatrixObjectType::MatrixType MatrixType
type of matrix
Definition: spmatrix.hh:872
virtual double overflowFraction() const
Definition: spmatrix.hh:52
void reserve(const Stencil &stencil, bool verbose=false)
reserve memory for assemble based on the provided stencil
Definition: spmatrix.hh:638
int dim(int i) const
return number of rows = 0, cols = 1
Definition: spmatrix.hh:238
SparseRowMatrixObject< DomainSpaceType, RangeSpaceType, Matrix > SparseRowMatrixObjectType
Definition: spmatrix.hh:841
virtual int method() const
Definition: spmatrix.hh:67
std::pair< const T, int > realValue(int row, int fakeCol) const
return pair< value, column >, used by BlockMatrix
Definition: spmatrix.hh:211
LocalMatrixWrapper< LocalMatrixStackType > LocalMatrixType
type of local matrix
Definition: spmatrix.hh:477
RangeSpace RangeSpaceType
Definition: spmatrix.hh:838
void createPreconditionMatrix()
Definition: spmatrix.hh:753
LocalMatrix(const MatrixObjectType &matrixObject, const DomainSpaceType &domainSpace, const RangeSpaceType &rangeSpace, const DomainMapperType &domainMapper, const RangeMapperType &rangeMapper)
Definition: spmatrix.hh:907
Definition: spmatrix.hh:79
void addLocalMatrix(const DomainEntityType &domainEntity, const RangeEntityType &rangeEntity, const LocalMatrix &localMat)
Definition: spmatrix.hh:568
void init(const DomainEntityType &domainEntity, const RangeEntityType &rangeEntity)
Definition: spmatrix.hh:921
DomainSpaceType::EntityType DomainEntityType
Definition: localmatrix.hh:64
RowSpaceImp RowSpaceType
Definition: spmatrix.hh:422
SparseRowMatrixObjectType::RangeMapperType RangeMapperType
Definition: spmatrix.hh:849