dune-fem  2.4.1-rc
localmassmatrix.hh
Go to the documentation of this file.
1 #ifndef DUNE_FEM_LOCALMASSMATRIX_HH
2 #define DUNE_FEM_LOCALMASSMATRIX_HH
3 
4 //- dune-common includes
5 #include <dune/common/dynmatrix.hh>
6 #include <dune/common/fmatrix.hh>
7 #include <dune/common/nullptr.hh>
8 
9 //- dune-geometry includes
10 #include <dune/geometry/typeindex.hh>
11 
12 //- dune-fem includes
16 #include <dune/fem/version.hh>
17 
18 namespace Dune
19 {
20 
21  namespace Fem
22  {
23 
30  template< class DiscreteFunctionSpace, class VolumeQuadrature >
32  {
34 
35  public:
37  typedef typename DiscreteFunctionSpaceType :: RangeFieldType ctype;
38  typedef typename DiscreteFunctionSpaceType :: RangeFieldType RangeFieldType;
39  typedef typename DiscreteFunctionSpaceType :: RangeType RangeType;
40 
41  enum { localBlockSize = DiscreteFunctionSpaceType :: localBlockSize };
43 
44  typedef Dune::FieldMatrix< ctype, dgNumDofs, dgNumDofs > DGMatrixType;
45  typedef Dune::FieldVector< ctype, dgNumDofs > DGVectorType;
46 
47  typedef typename DiscreteFunctionSpaceType :: GridPartType GridPartType;
48 
49  typedef typename DiscreteFunctionSpaceType :: IndexSetType IndexSetType;
50  typedef typename IndexSetType :: IndexType IndexType;
51 
52  typedef typename DiscreteFunctionSpaceType :: BasisFunctionSetType BasisFunctionSetType;
53 
54  typedef typename GridPartType :: GridType GridType;
55  typedef typename DiscreteFunctionSpaceType :: EntityType EntityType;
56  typedef typename EntityType :: Geometry Geometry;
57 
58  typedef VolumeQuadrature VolumeQuadratureType;
59 
61 
63  enum { StructuredGrid = Dune::Capabilities::isCartesian< GridType >::v };
64 
67 
68  // use dynamic matrix from dune-common
69  typedef Dune::DynamicMatrix< RangeFieldType > MatrixType;
70 
71  protected:
72  const DiscreteFunctionSpaceType& spc_;
73  const IndexSetType& indexSet_;
74 
75  GeometryInformationType geoInfo_;
76  const int volumeQuadOrd_;
77  const bool affine_;
78 
79  mutable DGMatrixType dgMatrix_;
80  mutable DGVectorType dgX_, dgRhs_;
81 
82  // use dynamic vector from dune-common
83  mutable Dune::DynamicVector< RangeFieldType > rhs_, row_;
84  mutable MatrixType matrix_;
85 
86  mutable std::vector< RangeType > phi_;
87  mutable std::vector< RangeType > phiMass_;
88 
89  typedef std::map< const int, MatrixType* > MassMatrixStorageType;
90  typedef std::vector< MassMatrixStorageType > LocalInverseMassMatrixStorageType;
91 
92  mutable LocalInverseMassMatrixStorageType localInverseMassMatrix_;
93 
94  // index of entity from index set, don't setup mass matrix for the same entity twice
95  mutable IndexType lastEntityIndex_;
96  mutable unsigned int lastTopologyId_ ;
97  // sequence number (obtained from DofManager via the space)
98  mutable int sequence_;
99 
101  {
102  static const int dimRange = DiscreteFunctionSpaceType::dimRange;
103 
104  typedef Dune::FieldMatrix< ctype, dimRange, dimRange > MassFactorType;
105 
106  // return false since we don;t have a mass term
107  bool hasMass () const { return false; }
108 
109  void mass ( const EntityType &, const VolumeQuadratureType &, int, MassFactorType & ) const
110  {}
111  };
112 
113  template< class BasisFunctionSet >
114  MatrixType &getLocalInverseMassMatrix ( const EntityType &entity, const Geometry &geo,
115  const BasisFunctionSet &basisSet, int numBasisFct ) const
116  {
117  const GeometryType geomType = geo.type();
118  typedef typename MassMatrixStorageType::iterator iterator;
119  MassMatrixStorageType &massMap = localInverseMassMatrix_[ GlobalGeometryTypeIndex::index( geomType ) ];
120 
121  std::pair< iterator, bool > insertPair = massMap.insert( std::make_pair( numBasisFct, nullptr ) );
122  MatrixType *&matrix = insertPair.first->second;
123 
124  if( insertPair.second )
125  {
126  matrix = new MatrixType( numBasisFct, numBasisFct, 0.0 );
127  VolumeQuadratureType volQuad( entity, volumeQuadratureOrder( entity ) );
128  buildMatrixNoMassFactor( entity, geo, basisSet, volQuad, numBasisFct, *matrix, false );
129  matrix->invert();
130  }
131 
132  return *matrix;
133  }
134 
135  template< class MassCaller, class BasisFunctionSet >
136  MatrixType &getLocalInverseMassMatrixDefault ( MassCaller &caller, const EntityType &entity,
137  const Geometry &geo, const BasisFunctionSet &basisSet ) const
138  {
139  const int numDofs = basisSet.size();
140  // if sequence changed or entity index changed: recompute mass matrix
141  if( entityHasChanged( entity ) || (numDofs != int( matrix_.rows())) )
142  {
143  // resize temporary memory if necessary
144  if( numDofs != int( matrix_.rows() ) )
145  matrix_.resize( numDofs, numDofs );
146 
147  buildMatrix( caller, entity, geo, basisSet, numDofs, matrix_ );
148  matrix_.invert();
149  }
150 
151  // make sure that rhs_ has the correct size
152  assert( numDofs == int( matrix_.rows() ) );
153  return matrix_;
154  }
155 
157  int volumeQuadratureOrder ( const EntityType &entity ) const
158  {
159  return (volumeQuadOrd_ < 0 ? 2 * spc_.order( entity ) : volumeQuadOrd_);
160  }
161 
164  {
165  return (volumeQuadOrd_ < 0 ? 2 * spc_.order() : volumeQuadOrd_);
166  }
167 
168  // return number of max non blocked dofs
169  int maxNumDofs () const
170  {
171  return spc_.blockMapper().maxNumDofs() * localBlockSize;
172  }
173 
174  public:
176  explicit LocalMassMatrixImplementation ( const DiscreteFunctionSpaceType &spc, int volQuadOrd = -1 )
177  : spc_(spc)
178  , indexSet_( spc.indexSet() )
179  , geoInfo_( indexSet_ )
180  , volumeQuadOrd_ ( volQuadOrd )
181  , affine_ ( setup() )
182  , rhs_(), row_(), matrix_()
183  , phi_( maxNumDofs() )
184  , phiMass_( maxNumDofs() )
185  , localInverseMassMatrix_( GlobalGeometryTypeIndex :: size( GridType::dimension ) )
186  , lastEntityIndex_( -1 )
187  , lastTopologyId_( ~0u )
188  , sequence_( -1 )
189  {}
190 
192  LocalMassMatrixImplementation ( const ThisType &other )
193  : spc_(other.spc_),
194  indexSet_( spc_.indexSet() ),
195  geoInfo_( indexSet_ ),
196  volumeQuadOrd_( other.volumeQuadOrd_ ),
197  affine_( other.affine_ ),
198  rhs_( other.rhs_ ), row_( other.row_ ), matrix_( other.matrix_ ),
199  phi_( other.phi_ ),
200  phiMass_( other.phiMass_ ),
201  localInverseMassMatrix_( GlobalGeometryTypeIndex :: size( GridType::dimension ) ),
202  lastEntityIndex_( other.lastEntityIndex_ ),
203  lastTopologyId_( other.lastTopologyId_ ),
204  sequence_( other.sequence_ )
205  {}
206 
208  {
209  typedef typename MassMatrixStorageType::iterator iterator;
210  for( unsigned int i = 0; i < localInverseMassMatrix_.size(); ++i )
211  {
212  const iterator end = localInverseMassMatrix_[ i ].end();
213  for( iterator it = localInverseMassMatrix_[ i ].begin(); it != end; ++it )
214  delete it->second;
215  }
216  }
217 
218  public:
220  bool affine () const { return affine_; }
221 
223  double getAffineMassFactor(const Geometry& geo) const
224  {
225  return geoInfo_.referenceVolume( geo.type() ) / geo.volume();
226  }
227 
230  template< class MassCaller, class LocalFunction >
231  void applyInverse ( MassCaller &caller, const EntityType &entity, LocalFunction &lf ) const
232  {
233  Geometry geo = entity.geometry();
234  if( affine() || geo.affine() )
235  applyInverseLocally( caller, entity, geo, lf );
236  else
237  applyInverseDefault( caller, entity, geo, lf );
238  }
239 
241  template< class LocalFunction >
242  void applyInverse ( const EntityType &entity, LocalFunction &lf ) const
243  {
244  NoMassDummyCaller caller;
245  applyInverse( caller, entity, lf );
246  }
247 
249  template< class LocalFunction >
250  void applyInverse ( LocalFunction &lf ) const
251  {
252  applyInverse( lf.entity(), lf );
253  }
254 
255  template< class LocalMatrix >
256  void rightMultiplyInverse ( LocalMatrix &localMatrix ) const
257  {
258  const EntityType &entity = localMatrix.domainEntity();
259  Geometry geo = entity.geometry();
260  if( affine() || geo.affine() )
261  rightMultiplyInverseLocally( entity, geo, localMatrix );
262  else
263  rightMultiplyInverseDefault( entity, geo, localMatrix );
264  }
265 
266  template< class LocalMatrix >
267  void leftMultiplyInverse ( LocalMatrix &localMatrix ) const
268  {
269  const EntityType &entity = localMatrix.domainEntity();
270  Geometry geo = entity.geometry();
271  if( affine() || geo.affine() )
272  leftMultiplyInverseLocally( entity, geo, localMatrix );
273  else
274  leftMultiplyInverseDefault( entity, geo, localMatrix );
275  }
276 
278  // end of public methods
280 
281  protected:
283  // applyInverse for DG spaces
285  template< class MassCaller, class LocalFunction >
286  void applyInverseDgOrthoNormalBasis ( MassCaller &caller, const EntityType &entity, LocalFunction &lf ) const
287  {
288  Geometry geo = entity.geometry();
289  assert( dgNumDofs == lf.numDofs() );
290 
291  // affine_ can be a static information
292  const bool isAffine = affine() || geo.affine();
293  // make sure that for affine grids the geometry info is also correct
294  assert( affine() ? geo.affine() : true );
295 
296  // in case of affine mappings we only have to multiply with a factor
297  if( isAffine && !caller.hasMass() )
298  {
299  const double massVolInv = getAffineMassFactor( geo );
300 
301  // apply inverse mass matrix
302  for( int l = 0; l < dgNumDofs; ++l )
303  lf[ l ] *= massVolInv;
304  }
305  else
306  {
307  // copy local function to right hand side
308  for( int l = 0; l < dgNumDofs; ++l )
309  dgRhs_[ l ] = lf[ l ];
310 
311  buildMatrix( caller, entity, geo, lf.basisFunctionSet(), dgNumDofs, dgMatrix_ );
312  dgMatrix_.solve( dgX_, dgRhs_ );
313 
314  // copy back to local function
315  for( int l = 0; l < dgNumDofs; ++l )
316  lf[ l ] = dgX_[ l ];
317  }
318  }
319 
320  template< class LocalMatrix >
321  void rightMultiplyInverseDgOrthoNormalBasis ( LocalMatrix &localMatrix ) const
322  {
323  const EntityType &entity = localMatrix.domainEntity();
324  Geometry geo = entity.geometry();
325  assert( dgNumDofs == localMatrix.columns() );
326 
327  // in case of affine mappings we only have to multiply with a factor
328  if( affine() || geo.affine() )
329  localMatrix.scale( getAffineMassFactor( geo ) );
330  else
331  {
332  NoMassDummyCaller caller;
333  buildMatrix( caller, entity, geo, localMatrix.domainBasisFunctionSet(), dgNumDofs, dgMatrix_ );
334  dgMatrix_.invert();
335 
336  const int rows = localMatrix.rows();
337  for( int i = 0; i < rows; ++i )
338  {
339  for( int j = 0; j < dgNumDofs; ++j )
340  dgRhs_[ j ] = localMatrix.get( i, j );
341  dgMatrix_.mtv( dgRhs_, dgX_ );
342  for( int j = 0; j < dgNumDofs; ++j )
343  localMatrix.set( i, j, dgX_[ j ] );
344  }
345  }
346  }
347 
348  template< class LocalMatrix >
349  void leftMultiplyInverseDgOrthoNormalBasis ( LocalMatrix &localMatrix ) const
350  {
351  const EntityType &entity = localMatrix.domainEntity();
352  Geometry geo = entity.geometry();
353  assert( dgNumDofs == localMatrix.columns() );
354 
355  // in case of affine mappings we only have to multiply with a factor
356  if( affine() || geo.affine() )
357  localMatrix.scale( getAffineMassFactor( geo ) );
358  else
359  {
360  NoMassDummyCaller caller;
361  buildMatrix( caller, entity, geo, localMatrix.domainBasisFunctionSet(), dgNumDofs, dgMatrix_ );
362  dgMatrix_.invert();
363 
364  leftMultiplyScaled( dgMatrix_, localMatrix );
365  }
366  }
367 
369  bool entityHasChanged( const EntityType& entity ) const
370  {
371  // don't compute matrix new for the same entity
372  const int currentSequence = spc_.sequence();
373  const unsigned int topologyId = entity.type().id();
374  const IndexType entityIndex = indexSet_.index( entity ) ;
375 
376  // check whether sequence has been updated
377  if( sequence_ != currentSequence ||
378  lastEntityIndex_ != entityIndex ||
379  lastTopologyId_ != topologyId )
380  {
381  // update identifiers
382  lastEntityIndex_ = entityIndex ;
383  sequence_ = currentSequence;
384  lastTopologyId_ = topologyId ;
385 
386  return true ;
387  }
388  else
389  // the entity did not change
390  return false ;
391  }
392 
394  // standard applyInverse method
398  template< class MassCaller, class LocalFunction >
399  void applyInverseDefault ( MassCaller &caller, const EntityType &entity,
400  const Geometry &geo, LocalFunction &lf ) const
401  {
402  // get local inverted mass matrix
403  MatrixType &invMassMatrix
404  = getLocalInverseMassMatrixDefault ( caller, entity, geo, lf.basisFunctionSet() );
405 
406  // copy local function to right hand side
407  const int numDofs = lf.numDofs();
408  rhs_.resize( numDofs );
409  for( int l = 0; l < numDofs; ++l )
410  rhs_[ l ] = lf[ l ];
411 
412  // apply inverse to right hand side and store in lf
413  multiply( numDofs, invMassMatrix, rhs_, lf );
414  }
415 
416  template< class LocalMatrix >
417  void rightMultiplyInverseDefault ( const EntityType &entity, const Geometry &geo, LocalMatrix &localMatrix ) const
418  {
419  NoMassDummyCaller caller;
420  MatrixType &invMassMatrix
421  = getLocalInverseMassMatrixDefault ( caller, entity, geo, localMatrix.domainBasisFunctionSet() );
422 
423  const int cols = localMatrix.columns();
424  rhs_.resize( cols );
425  row_.resize( cols );
426 
427  const int rows = localMatrix.rows();
428  for( int i = 0; i < rows; ++i )
429  {
430  for( int j = 0; j < cols; ++j )
431  rhs_[ j ] = localMatrix.get( i, j );
432  invMassMatrix.mtv( rhs_, row_ );
433  for( int j = 0; j < cols; ++j )
434  localMatrix.set( i, j, row_[ j ] );
435  }
436  }
437 
438  template< class LocalMatrix >
439  void leftMultiplyInverseDefault ( const EntityType &entity, const Geometry &geo, LocalMatrix &localMatrix ) const
440  {
441  NoMassDummyCaller caller;
442  MatrixType &invMassMatrix
443  = getLocalInverseMassMatrixDefault ( caller, entity, geo, localMatrix.rangeBasisFunctionSet() );
444 
445  const int cols = localMatrix.columns();
446  rhs_.resize( cols );
447  row_.resize( cols );
448 
449  const int rows = localMatrix.rows();
450  for( int i = 0; i < rows; ++i )
451  {
452  for( int j = 0; j < cols; ++j )
453  rhs_[ j ] = localMatrix.get( i, j );
454  invMassMatrix.mv( rhs_, row_ );
455  for( int j = 0; j < cols; ++j )
456  localMatrix.set( i, j, row_[ j ] );
457  }
458  }
459 
460 
462  // local applyInverse method for affine geometries
466  template< class MassCaller, class LocalFunction >
467  void applyInverseLocally ( MassCaller &caller, const EntityType &entity,
468  const Geometry &geo, LocalFunction &lf ) const
469  {
470  const int numDofs = lf.numDofs();
471 
472  // get local inverted mass matrix
473  MatrixType &invMassMatrix =
474  getLocalInverseMassMatrix( entity, geo, lf.basisFunctionSet(), numDofs );
475 
476  const double massVolInv = getAffineMassFactor( geo );
477 
478  // copy local function to right hand side
479  // and apply inverse mass volume fraction
480  rhs_.resize( numDofs );
481  for( int l = 0; l < numDofs; ++l )
482  rhs_[ l ] = lf[ l ] * massVolInv;
483 
484  // apply inverse local mass matrix and store in lf
485  multiply( numDofs, invMassMatrix, rhs_, lf );
486  }
487 
488  template< class LocalMatrix >
489  void rightMultiplyInverseLocally ( const EntityType &entity, const Geometry &geo, LocalMatrix &localMatrix ) const
490  {
491  const int cols = localMatrix.colums();
492  MatrixType &invMassMatrix =
493  getLocalInverseMassMatrix( entity, geo, localMatrix.domainBasisFunctionSet(), cols );
494 
495  const double massVolInv = getAffineMassFactor( geo );
496 
497  rhs_.resize( cols );
498  row_.resize( cols );
499 
500  const int rows = localMatrix.rows();
501  for( int i = 0; i < rows; ++i )
502  {
503  for( int j = 0; j < cols; ++j )
504  rhs_[ j ] = localMatrix.get( i, j ) * massVolInv;
505  invMassMatrix.mtv( rhs_, row_ );
506  for( int j = 0; j < cols; ++j )
507  localMatrix.set( i, j, row_[ j ] );
508  }
509  }
510 
511  template< class LocalMatrix >
512  void leftMultiplyInverseLocally ( const EntityType &entity, const Geometry &geo, LocalMatrix &localMatrix ) const
513  {
514  const int cols = localMatrix.columns();
515  MatrixType &invMassMatrix =
516  getLocalInverseMassMatrix( entity, geo, localMatrix.rangeBasisFunctionSet(), cols );
517 
518  const double massVolInv = getAffineMassFactor( geo );
519 
520  rhs_.resize( cols );
521  row_.resize( cols );
522 
523  const int rows = localMatrix.rows();
524  for( int i = 0; i < rows; ++i )
525  {
526  for( int j = 0; j < cols; ++j )
527  rhs_[ j ] = localMatrix.get( i, j ) * massVolInv;
528  invMassMatrix.mv( rhs_, row_ );
529  for( int j = 0; j < cols; ++j )
530  localMatrix.set( i, j, row_[ j ] );
531  }
532  }
533 
534 
536  bool setup () const
537  {
538  // for structured grids this is always true
539  if( StructuredGrid )
540  return true;
541 
542  // get types for codim 0
543  const std::vector<Dune::GeometryType>& geomTypes = geoInfo_.geomTypes(0);
544 
545  // for simplices we also have affine mappings
546  if( (geomTypes.size() == 1) && geomTypes[0].isSimplex() )
547  {
548  return true;
549  }
550 
551  // otherwise use geometry affinity
552  return false ;
553  }
554 
556  template< class MassCaller, class Matrix >
557  void buildMatrix ( MassCaller &caller, const EntityType &entity,
558  const Geometry &geo, const BasisFunctionSetType &set,
559  std::size_t numDofs, Matrix &matrix ) const
560  {
561  assert( numDofs == set.size() );
562 
563  // clear matrix
564  matrix = 0;
565 
566  // create quadrature
567  VolumeQuadratureType volQuad( entity, volumeQuadratureOrder( entity ) );
568 
569  if( caller.hasMass() )
570  buildMatrixWithMassFactor( caller, entity, geo, set, volQuad, numDofs, matrix );
571  else
572  buildMatrixNoMassFactor( entity, geo, set, volQuad, numDofs, matrix );
573  }
574 
576  template <class Matrix>
578  const EntityType& en,
579  const Geometry& geo,
580  const BasisFunctionSetType& set,
581  const VolumeQuadratureType& volQuad,
582  const int numDofs,
583  Matrix& matrix,
584  const bool applyIntegrationElement = true ) const
585  {
586  const int volNop = volQuad.nop();
587  for(int qp=0; qp<volNop; ++qp)
588  {
589  // calculate integration weight
590  const double intel = ( applyIntegrationElement ) ?
591  ( volQuad.weight(qp) * geo.integrationElement(volQuad.point(qp)) ) : volQuad.weight(qp) ;
592 
593  // eval basis functions
594  set.evaluateAll(volQuad[qp], phi_);
595 
596  for(int m=0; m<numDofs; ++m)
597  {
598  const RangeType& phi_m = phi_[m];
599  const ctype val = intel * (phi_m * phi_m);
600  matrix[m][m] += val;
601 
602  for(int k=m+1; k<numDofs; ++k)
603  {
604  const ctype val = intel * (phi_m * phi_[k]);
605  matrix[m][k] += val;
606  matrix[k][m] += val;
607  }
608  }
609  }
610  }
611 
613  template <class MassCallerType, class Matrix>
615  MassCallerType& caller,
616  const EntityType& en,
617  const Geometry& geo,
618  const BasisFunctionSetType& set,
619  const VolumeQuadratureType& volQuad,
620  const int numDofs,
621  Matrix& matrix) const
622  {
623  typedef typename MassCallerType :: MassFactorType MassFactorType;
624  MassFactorType mass;
625 
626  const int volNop = volQuad.nop();
627  for(int qp=0; qp<volNop; ++qp)
628  {
629  // calculate integration weight
630  const double intel = volQuad.weight(qp)
631  * geo.integrationElement(volQuad.point(qp));
632 
633  // eval basis functions
634  set.evaluateAll( volQuad[qp], phi_);
635 
636  // call mass factor
637  caller.mass( en, volQuad, qp, mass);
638 
639  // apply mass matrix to all basis functions
640  for(int m=0; m<numDofs; ++m)
641  {
642  mass.mv( phi_[m], phiMass_[m] );
643  }
644 
645  // add values to matrix
646  for(int m=0; m<numDofs; ++m)
647  {
648  for(int k=0; k<numDofs; ++k)
649  {
650  matrix[m][k] += intel * (phiMass_[m] * phi_[k]);
651  }
652  }
653  }
654  }
655 
656  // implement matvec with matrix (mv of densematrix is too stupid)
657  template <class Matrix, class Rhs, class X>
658  void multiply( const int size,
659  const Matrix& matrix,
660  const Rhs& rhs,
661  X& x ) const
662  {
663  assert( (int) matrix.rows() == size );
664  assert( (int) matrix.cols() == size );
665  assert( (int) rhs.size() == size );
666 
667  for( int row = 0; row < size; ++ row )
668  {
669  RangeFieldType sum = 0;
670  // get matrix row
671  typedef typename Matrix :: const_row_reference MatRow;
672  MatRow matRow = matrix[ row ];
673 
674  // multiply row with right hand side
675  for( int col = 0; col < size; ++ col )
676  {
677  sum += matRow[ col ] * rhs[ col ];
678  }
679 
680  // set to result to result vector
681  x[ row ] = sum;
682  }
683  }
684  };
685 
686 
687 
688  // LocalMassMatrix
689  // ---------------
690 
692  template< class DiscreteFunctionSpace, class VolumeQuadrature >
694  : public LocalMassMatrixImplementation< DiscreteFunctionSpace, VolumeQuadrature >
695  {
697 
698  public:
699  explicit LocalMassMatrix ( const DiscreteFunctionSpace &spc, int volQuadOrd = -1 )
700  : BaseType( spc, volQuadOrd )
701  {}
702  };
703 
704 
705 
707  //
708  // DG LocalMassMatrix Implementation
709  //
711 
713  template< class DiscreteFunctionSpace, class VolumeQuadrature >
715  : public LocalMassMatrixImplementation< DiscreteFunctionSpace, VolumeQuadrature >
716  {
718 
719  public:
721 
722  explicit LocalMassMatrixImplementationDgOrthoNormal ( const DiscreteFunctionSpace &spc, int volQuadOrd = -1 )
723  : BaseType( spc, volQuadOrd )
724  {}
725 
728  template <class MassCallerType, class LocalFunctionType>
729  void applyInverse(MassCallerType& caller,
730  const EntityType& en,
731  LocalFunctionType& lf) const
732  {
733  BaseType :: applyInverseDgOrthoNormalBasis( caller, en, lf );
734  }
735 
737  template <class LocalFunctionType>
738  void applyInverse(const EntityType& en,
739  LocalFunctionType& lf) const
740  {
741  typename BaseType :: NoMassDummyCaller caller;
742  applyInverse(caller, en, lf );
743  }
744 
746  template< class LocalFunction >
747  void applyInverse ( LocalFunction &lf ) const
748  {
749  applyInverse( lf.entity(), lf );
750  }
751 
752  template< class LocalMatrix >
753  void rightMultiplyInverse ( LocalMatrix &localMatrix ) const
754  {
755  BaseType::rightMultiplyInverseDgOrthoNormalBasis( localMatrix );
756  }
757  };
758 
760 
761  } // namespace Fem
762 
763 } // namespace Dune
764 
765 #endif // #ifndef DUNE_FEM_LOCALMASSMATRIX_HH
static const int dimRange
Definition: localmassmatrix.hh:102
AllGeomTypes< typename GridPartType::IndexSetType, GridType > GeometryInformationType
Definition: localmassmatrix.hh:65
DGMatrixType dgMatrix_
Definition: localmassmatrix.hh:79
DiscreteFunctionSpaceType::EntityType EntityType
Definition: localmassmatrix.hh:55
void rightMultiplyInverseLocally(const EntityType &entity, const Geometry &geo, LocalMatrix &localMatrix) const
Definition: localmassmatrix.hh:489
int numDofs() const
obtain the number of local DoFs
Definition: localfunction.hh:340
void rightMultiplyInverseDefault(const EntityType &entity, const Geometry &geo, LocalMatrix &localMatrix) const
Definition: localmassmatrix.hh:417
void mass(const EntityType &, const VolumeQuadratureType &, int, MassFactorType &) const
Definition: localmassmatrix.hh:109
Dune::DynamicVector< RangeFieldType > rhs_
Definition: localmassmatrix.hh:83
void rightMultiplyInverse(LocalMatrix &localMatrix) const
Definition: localmassmatrix.hh:753
void applyInverseDgOrthoNormalBasis(MassCaller &caller, const EntityType &entity, LocalFunction &lf) const
Definition: localmassmatrix.hh:286
GridPartType::GridType GridType
Definition: localmassmatrix.hh:54
void buildMatrixNoMassFactor(const EntityType &en, const Geometry &geo, const BasisFunctionSetType &set, const VolumeQuadratureType &volQuad, const int numDofs, Matrix &matrix, const bool applyIntegrationElement=true) const
build local mass matrix with mass factor
Definition: localmassmatrix.hh:577
void multiply(const int size, const Matrix &matrix, const Rhs &rhs, X &x) const
Definition: localmassmatrix.hh:658
DiscreteFunctionSpaceType::RangeType RangeType
Definition: localmassmatrix.hh:39
DiscreteFunctionSpaceType::BasisFunctionSetType BasisFunctionSetType
Definition: localmassmatrix.hh:52
static constexpr T sum(T a)
Definition: utility.hh:33
const IndexSetType & indexSet_
Definition: localmassmatrix.hh:73
DiscreteFunctionSpaceType::IndexSetType IndexSetType
Definition: localmassmatrix.hh:49
void leftMultiplyInverseDefault(const EntityType &entity, const Geometry &geo, LocalMatrix &localMatrix) const
Definition: localmassmatrix.hh:439
Local Mass Matrix inversion implementation, select the correct method in your implementation.
Definition: localmassmatrix.hh:31
LocalMassMatrix(const DiscreteFunctionSpace &spc, int volQuadOrd=-1)
Definition: localmassmatrix.hh:699
LocalInverseMassMatrixStorageType localInverseMassMatrix_
Definition: localmassmatrix.hh:92
GeometryInformationType geoInfo_
Definition: localmassmatrix.hh:75
Dune::DynamicMatrix< RangeFieldType > MatrixType
Definition: localmassmatrix.hh:69
int maxNumDofs() const
Definition: localmassmatrix.hh:169
void leftMultiplyInverseLocally(const EntityType &entity, const Geometry &geo, LocalMatrix &localMatrix) const
Definition: localmassmatrix.hh:512
void buildMatrix(MassCaller &caller, const EntityType &entity, const Geometry &geo, const BasisFunctionSetType &set, std::size_t numDofs, Matrix &matrix) const
build local mass matrix
Definition: localmassmatrix.hh:557
void applyInverse(MassCallerType &caller, const EntityType &en, LocalFunctionType &lf) const
Definition: localmassmatrix.hh:729
const BasisFunctionSetType & basisFunctionSet() const
obtain the basis function set for this local function
Definition: localfunction.hh:279
void applyInverse(MassCaller &caller, const EntityType &entity, LocalFunction &lf) const
Definition: localmassmatrix.hh:231
DiscreteFunctionSpace DiscreteFunctionSpaceType
Definition: localmassmatrix.hh:36
bool entityHasChanged(const EntityType &entity) const
returns true if the entity has been changed
Definition: localmassmatrix.hh:369
interface for local functions
Definition: localfunction.hh:41
Local Mass Matrix for arbitrary spaces.
Definition: localmassmatrix.hh:693
int volumeQuadratureOrder(const EntityType &entity) const
return appropriate quadrature order, default is 2 * order(entity)
Definition: localmassmatrix.hh:157
Dune::DynamicVector< RangeFieldType > row_
Definition: localmassmatrix.hh:83
void applyInverse(LocalFunction &lf) const
apply local dg mass matrix to local function lf without mass factor
Definition: localmassmatrix.hh:250
MatrixType matrix_
Definition: localmassmatrix.hh:84
const DiscreteFunctionSpaceType & spc_
Definition: localmassmatrix.hh:72
EntityType::Geometry Geometry
Definition: localmassmatrix.hh:56
IndexSetType::IndexType IndexType
Definition: localmassmatrix.hh:50
double getAffineMassFactor(const Geometry &geo) const
return mass factor for diagonal mass matrix
Definition: localmassmatrix.hh:223
std::map< const int, MatrixType * > MassMatrixStorageType
Definition: localmassmatrix.hh:89
MatrixType & getLocalInverseMassMatrixDefault(MassCaller &caller, const EntityType &entity, const Geometry &geo, const BasisFunctionSet &basisSet) const
Definition: localmassmatrix.hh:136
DGVectorType dgX_
Definition: localmassmatrix.hh:80
Dune::FieldVector< ctype, dgNumDofs > DGVectorType
Definition: localmassmatrix.hh:45
void buildMatrixWithMassFactor(MassCallerType &caller, const EntityType &en, const Geometry &geo, const BasisFunctionSetType &set, const VolumeQuadratureType &volQuad, const int numDofs, Matrix &matrix) const
build local mass matrix with mass factor
Definition: localmassmatrix.hh:614
double referenceVolume(const GeometryType &type) const
return volume of reference element for geometry of type type
Definition: allgeomtypes.hh:67
Definition: coordinate.hh:4
std::vector< RangeType > phi_
Definition: localmassmatrix.hh:86
Helper class to check affinity of the grid&#39;s geometries.
Definition: checkgeomaffinity.hh:22
VolumeQuadrature VolumeQuadratureType
Definition: localmassmatrix.hh:58
void leftMultiplyInverse(LocalMatrix &localMatrix) const
Definition: localmassmatrix.hh:267
LocalMassMatrixImplementationDgOrthoNormal(const DiscreteFunctionSpace &spc, int volQuadOrd=-1)
Definition: localmassmatrix.hh:722
const EntityType & entity() const
obtain the entity, this local function lives on
Definition: localfunction.hh:285
void applyInverseDefault(MassCaller &caller, const EntityType &entity, const Geometry &geo, LocalFunction &lf) const
Definition: localmassmatrix.hh:399
void applyInverse(const EntityType &en, LocalFunctionType &lf) const
apply local dg mass matrix to local function lf without mass factor
Definition: localmassmatrix.hh:738
DiscreteFunctionSpaceType::GridPartType GridPartType
Definition: localmassmatrix.hh:47
LocalMassMatrixImplementation(const ThisType &other)
copy constructor
Definition: localmassmatrix.hh:192
int maxVolumeQuadratureOrder() const
return appropriate quadrature order, default is 2 * order()
Definition: localmassmatrix.hh:163
const bool affine_
Definition: localmassmatrix.hh:77
void applyInverse(LocalFunction &lf) const
apply local dg mass matrix to local function lf without mass factor
Definition: localmassmatrix.hh:747
std::vector< MassMatrixStorageType > LocalInverseMassMatrixStorageType
Definition: localmassmatrix.hh:90
void rightMultiplyInverse(LocalMatrix &localMatrix) const
Definition: localmassmatrix.hh:256
void rightMultiplyInverseDgOrthoNormalBasis(LocalMatrix &localMatrix) const
Definition: localmassmatrix.hh:321
GeometryInformationType::DomainType DomainType
Definition: localmassmatrix.hh:66
const int volumeQuadOrd_
Definition: localmassmatrix.hh:76
Dune::FieldMatrix< ctype, dimRange, dimRange > MassFactorType
Definition: localmassmatrix.hh:104
BaseType::EntityType EntityType
Definition: localmassmatrix.hh:720
MatrixType & getLocalInverseMassMatrix(const EntityType &entity, const Geometry &geo, const BasisFunctionSet &basisSet, int numBasisFct) const
Definition: localmassmatrix.hh:114
discrete function space
bool setup() const
setup and return affinity
Definition: localmassmatrix.hh:536
~LocalMassMatrixImplementation()
Definition: localmassmatrix.hh:207
Dune::FieldMatrix< ctype, dgNumDofs, dgNumDofs > DGMatrixType
Definition: localmassmatrix.hh:44
void applyInverseLocally(MassCaller &caller, const EntityType &entity, const Geometry &geo, LocalFunction &lf) const
Definition: localmassmatrix.hh:467
const std::vector< GeometryType > & geomTypes(unsigned int codim) const
returns vector with geometry tpyes this index set has indices for
Definition: allgeomtypes.hh:145
std::vector< RangeType > phiMass_
Definition: localmassmatrix.hh:87
void leftMultiplyInverseDgOrthoNormalBasis(LocalMatrix &localMatrix) const
Definition: localmassmatrix.hh:349
FieldVector< ctype, dim > DomainType
type of domain vector
Definition: allgeomtypes.hh:46
DiscreteFunctionSpaceType::RangeFieldType RangeFieldType
Definition: localmassmatrix.hh:38
bool affine() const
returns true if geometry mapping is affine
Definition: localmassmatrix.hh:220
Fem::GeometryAffinityCheck< VolumeQuadratureType > GeometryAffinityCheckType
Definition: localmassmatrix.hh:60
unsigned int lastTopologyId_
Definition: localmassmatrix.hh:96
void applyInverse(const EntityType &entity, LocalFunction &lf) const
apply local dg mass matrix to local function lf without mass factor
Definition: localmassmatrix.hh:242
DG Local Mass Matrix for arbitrary spaces.
Definition: localmassmatrix.hh:714
IndexType lastEntityIndex_
Definition: localmassmatrix.hh:95
DiscreteFunctionSpaceType::RangeFieldType ctype
Definition: localmassmatrix.hh:37
int sequence_
Definition: localmassmatrix.hh:98
Definition: basisfunctionset/basisfunctionset.hh:31
LocalMassMatrixImplementation(const DiscreteFunctionSpaceType &spc, int volQuadOrd=-1)
constructor taking space and volume quadrature order
Definition: localmassmatrix.hh:176
DGVectorType dgRhs_
Definition: localmassmatrix.hh:80
bool hasMass() const
Definition: localmassmatrix.hh:107
std::size_t size() const
return size of basis function set