1 #ifndef DUNE_FEM_LOCALMASSMATRIX_HH 2 #define DUNE_FEM_LOCALMASSMATRIX_HH 5 #include <dune/common/dynmatrix.hh> 6 #include <dune/common/fmatrix.hh> 7 #include <dune/common/nullptr.hh> 10 #include <dune/geometry/typeindex.hh> 30 template<
class DiscreteFunctionSpace,
class VolumeQuadrature >
37 typedef typename DiscreteFunctionSpaceType :: RangeFieldType
ctype;
39 typedef typename DiscreteFunctionSpaceType :: RangeType
RangeType;
44 typedef Dune::FieldMatrix< ctype, dgNumDofs, dgNumDofs >
DGMatrixType;
47 typedef typename DiscreteFunctionSpaceType :: GridPartType
GridPartType;
49 typedef typename DiscreteFunctionSpaceType :: IndexSetType
IndexSetType;
50 typedef typename IndexSetType :: IndexType
IndexType;
54 typedef typename GridPartType :: GridType
GridType;
55 typedef typename DiscreteFunctionSpaceType :: EntityType
EntityType;
56 typedef typename EntityType :: Geometry
Geometry;
69 typedef Dune::DynamicMatrix< RangeFieldType >
MatrixType;
72 const DiscreteFunctionSpaceType&
spc_;
83 mutable Dune::DynamicVector< RangeFieldType >
rhs_,
row_;
86 mutable std::vector< RangeType >
phi_;
102 static const int dimRange = DiscreteFunctionSpaceType::dimRange;
109 void mass (
const EntityType &,
const VolumeQuadratureType &,
int, MassFactorType & )
const 113 template<
class BasisFunctionSet >
117 const GeometryType geomType = geo.type();
118 typedef typename MassMatrixStorageType::iterator iterator;
119 MassMatrixStorageType &massMap = localInverseMassMatrix_[ GlobalGeometryTypeIndex::index( geomType ) ];
121 std::pair< iterator, bool > insertPair = massMap.insert( std::make_pair( numBasisFct,
nullptr ) );
122 MatrixType *&matrix = insertPair.first->second;
124 if( insertPair.second )
126 matrix =
new MatrixType( numBasisFct, numBasisFct, 0.0 );
135 template<
class MassCaller,
class BasisFunctionSet >
139 const int numDofs = basisSet.
size();
144 if( numDofs !=
int( matrix_.rows() ) )
145 matrix_.resize( numDofs, numDofs );
147 buildMatrix( caller, entity, geo, basisSet, numDofs, matrix_ );
152 assert( numDofs ==
int( matrix_.rows() ) );
159 return (volumeQuadOrd_ < 0 ? 2 * spc_.order( entity ) :
volumeQuadOrd_);
178 , indexSet_( spc.indexSet() )
179 , geoInfo_( indexSet_ )
180 , volumeQuadOrd_ ( volQuadOrd )
181 , affine_ (
setup() )
182 , rhs_(), row_(), matrix_()
185 , localInverseMassMatrix_( GlobalGeometryTypeIndex :: size( GridType::dimension ) )
186 , lastEntityIndex_( -1 )
187 , lastTopologyId_( ~0u )
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_ ),
200 phiMass_( other.phiMass_ ),
201 localInverseMassMatrix_( GlobalGeometryTypeIndex :: size( GridType::dimension ) ),
202 lastEntityIndex_( other.lastEntityIndex_ ),
203 lastTopologyId_( other.lastTopologyId_ ),
204 sequence_( other.sequence_ )
209 typedef typename MassMatrixStorageType::iterator iterator;
210 for(
unsigned int i = 0; i < localInverseMassMatrix_.size(); ++i )
212 const iterator end = localInverseMassMatrix_[ i ].end();
213 for( iterator it = localInverseMassMatrix_[ i ].begin(); it != end; ++it )
230 template<
class MassCaller,
class LocalFunction >
233 Geometry geo = entity.geometry();
234 if(
affine() || geo.affine() )
241 template<
class LocalFunction >
244 NoMassDummyCaller caller;
249 template<
class LocalFunction >
255 template<
class LocalMatrix >
258 const EntityType &entity = localMatrix.domainEntity();
259 Geometry geo = entity.geometry();
260 if(
affine() || geo.affine() )
266 template<
class LocalMatrix >
269 const EntityType &entity = localMatrix.domainEntity();
270 Geometry geo = entity.geometry();
271 if(
affine() || geo.affine() )
285 template<
class MassCaller,
class LocalFunction >
288 Geometry geo = entity.geometry();
292 const bool isAffine =
affine() || geo.affine();
294 assert(
affine() ? geo.affine() : true );
297 if( isAffine && !caller.hasMass() )
303 lf[ l ] *= massVolInv;
309 dgRhs_[ l ] = lf[ l ];
312 dgMatrix_.solve( dgX_, dgRhs_ );
320 template<
class LocalMatrix >
323 const EntityType &entity = localMatrix.domainEntity();
324 Geometry geo = entity.geometry();
325 assert(
dgNumDofs == localMatrix.columns() );
328 if(
affine() || geo.affine() )
332 NoMassDummyCaller caller;
336 const int rows = localMatrix.rows();
337 for(
int i = 0; i < rows; ++i )
340 dgRhs_[ j ] = localMatrix.get( i, j );
341 dgMatrix_.mtv( dgRhs_, dgX_ );
343 localMatrix.set( i, j, dgX_[ j ] );
348 template<
class LocalMatrix >
351 const EntityType &entity = localMatrix.domainEntity();
352 Geometry geo = entity.geometry();
353 assert(
dgNumDofs == localMatrix.columns() );
356 if(
affine() || geo.affine() )
360 NoMassDummyCaller caller;
364 leftMultiplyScaled( dgMatrix_, localMatrix );
372 const int currentSequence = spc_.sequence();
373 const unsigned int topologyId = entity.type().id();
374 const IndexType entityIndex = indexSet_.index( entity ) ;
377 if( sequence_ != currentSequence ||
378 lastEntityIndex_ != entityIndex ||
379 lastTopologyId_ != topologyId )
382 lastEntityIndex_ = entityIndex ;
383 sequence_ = currentSequence;
384 lastTopologyId_ = topologyId ;
398 template<
class MassCaller,
class LocalFunction >
403 MatrixType &invMassMatrix
407 const int numDofs = lf.
numDofs();
408 rhs_.resize( numDofs );
409 for(
int l = 0; l < numDofs; ++l )
413 multiply( numDofs, invMassMatrix, rhs_, lf );
416 template<
class LocalMatrix >
419 NoMassDummyCaller caller;
420 MatrixType &invMassMatrix
423 const int cols = localMatrix.columns();
427 const int rows = localMatrix.rows();
428 for(
int i = 0; i < rows; ++i )
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 ] );
438 template<
class LocalMatrix >
441 NoMassDummyCaller caller;
442 MatrixType &invMassMatrix
445 const int cols = localMatrix.columns();
449 const int rows = localMatrix.rows();
450 for(
int i = 0; i < rows; ++i )
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 ] );
466 template<
class MassCaller,
class LocalFunction >
470 const int numDofs = lf.
numDofs();
473 MatrixType &invMassMatrix =
480 rhs_.resize( numDofs );
481 for(
int l = 0; l < numDofs; ++l )
482 rhs_[ l ] = lf[ l ] * massVolInv;
485 multiply( numDofs, invMassMatrix, rhs_, lf );
488 template<
class LocalMatrix >
491 const int cols = localMatrix.colums();
492 MatrixType &invMassMatrix =
500 const int rows = localMatrix.rows();
501 for(
int i = 0; i < rows; ++i )
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 ] );
511 template<
class LocalMatrix >
514 const int cols = localMatrix.columns();
515 MatrixType &invMassMatrix =
523 const int rows = localMatrix.rows();
524 for(
int i = 0; i < rows; ++i )
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 ] );
543 const std::vector<Dune::GeometryType>& geomTypes = geoInfo_.
geomTypes(0);
546 if( (geomTypes.size() == 1) && geomTypes[0].isSimplex() )
556 template<
class MassCaller,
class Matrix >
558 const Geometry &geo,
const BasisFunctionSetType &
set,
559 std::size_t numDofs, Matrix &matrix )
const 561 assert( numDofs ==
set.size() );
569 if( caller.hasMass() )
576 template <
class Matrix>
578 const EntityType& en,
580 const BasisFunctionSetType&
set,
581 const VolumeQuadratureType& volQuad,
584 const bool applyIntegrationElement =
true )
const 586 const int volNop = volQuad.nop();
587 for(
int qp=0; qp<volNop; ++qp)
590 const double intel = ( applyIntegrationElement ) ?
591 ( volQuad.weight(qp) * geo.integrationElement(volQuad.point(qp)) ) : volQuad.weight(qp) ;
594 set.evaluateAll(volQuad[qp], phi_);
596 for(
int m=0; m<numDofs; ++m)
598 const RangeType& phi_m = phi_[m];
599 const ctype val = intel * (phi_m * phi_m);
602 for(
int k=m+1; k<numDofs; ++k)
604 const ctype val = intel * (phi_m * phi_[k]);
613 template <
class MassCallerType,
class Matrix>
615 MassCallerType& caller,
616 const EntityType& en,
618 const BasisFunctionSetType&
set,
619 const VolumeQuadratureType& volQuad,
621 Matrix& matrix)
const 626 const int volNop = volQuad.nop();
627 for(
int qp=0; qp<volNop; ++qp)
630 const double intel = volQuad.weight(qp)
631 * geo.integrationElement(volQuad.point(qp));
634 set.evaluateAll( volQuad[qp], phi_);
637 caller.mass( en, volQuad, qp, mass);
640 for(
int m=0; m<numDofs; ++m)
642 mass.mv( phi_[m], phiMass_[m] );
646 for(
int m=0; m<numDofs; ++m)
648 for(
int k=0; k<numDofs; ++k)
650 matrix[m][k] += intel * (phiMass_[m] * phi_[k]);
657 template <
class Matrix,
class Rhs,
class X>
659 const Matrix& matrix,
663 assert( (
int) matrix.rows() == size );
664 assert( (
int) matrix.cols() == size );
665 assert( (
int) rhs.size() == size );
667 for(
int row = 0; row < size; ++ row )
669 RangeFieldType
sum = 0;
671 typedef typename Matrix :: const_row_reference MatRow;
672 MatRow matRow = matrix[ row ];
675 for(
int col = 0; col < size; ++ col )
677 sum += matRow[ col ] * rhs[ col ];
692 template<
class DiscreteFunctionSpace,
class VolumeQuadrature >
700 : BaseType( spc, volQuadOrd )
713 template<
class DiscreteFunctionSpace,
class VolumeQuadrature >
723 : BaseType( spc, volQuadOrd )
728 template <
class MassCallerType,
class LocalFunctionType>
730 const EntityType& en,
731 LocalFunctionType& lf)
const 733 BaseType :: applyInverseDgOrthoNormalBasis( caller, en, lf );
737 template <
class LocalFunctionType>
739 LocalFunctionType& lf)
const 741 typename BaseType :: NoMassDummyCaller caller;
746 template<
class LocalFunction >
752 template<
class LocalMatrix >
755 BaseType::rightMultiplyInverseDgOrthoNormalBasis( localMatrix );
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
Definition: localmassmatrix.hh:100
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
Definition: localmassmatrix.hh:42
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's geometries.
Definition: checkgeomaffinity.hh:22
VolumeQuadrature VolumeQuadratureType
Definition: localmassmatrix.hh:58
Definition: localmassmatrix.hh:41
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
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: localmassmatrix.hh:63
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