DUNE-FEM (unstable)

Dune::Fem::LocalMassMatrix< DiscreteFunctionSpace, VolumeQuadrature > Class Template Reference

Local Mass Matrix for arbitrary spaces. More...

#include <dune/fem/operator/1order/localmassmatrix.hh>

Public Types

enum  
 is true if grid is structured grid
 

Public Member Functions

int volumeQuadratureOrder (const EntityType &entity) const
 return appropriate quadrature order, default is 2 * order(entity)
 
bool affine () const
 returns true if geometry mapping is affine
 
double getAffineMassFactor (const Geometry &geo) const
 return mass factor for diagonal mass matrix
 
template<class MassCaller , class BasisFunctionSet , class LocalFunction >
void applyInverse (MassCaller &caller, const EntityType &entity, const BasisFunctionSet &basisFunctionSet, LocalFunction &lf) const
 
template<class MassCaller , class LocalFunction >
void applyInverse (MassCaller &caller, const EntityType &entity, LocalFunction &lf) const
 
template<class LocalFunction >
void applyInverse (const EntityType &entity, LocalFunction &lf) const
 apply local dg mass matrix to local function lf without mass factor
 
template<class LocalFunction >
void applyInverse (LocalFunction &lf) const
 apply local dg mass matrix to local function lf without mass factor
 
template<class LocalMatrix >
void rightMultiplyInverse (LocalMatrix &localMatrix) const
 compute localMatrix * M^-1
 
template<class LocalMatrix >
void leftMultiplyInverse (LocalMatrix &localMatrix) const
 compute M^-1 * localMatrix
 

Protected Member Functions

int maxVolumeQuadratureOrder () const
 return appropriate quadrature order, default is 2 * order()
 
template<class LocalMatrix >
void rightMultiplyInverseDgOrthoNormalBasis (LocalMatrix &localMatrix) const
 compute localMatrix * M^-1
 
template<class LocalMatrix >
void leftMultiplyInverseDgOrthoNormalBasis (LocalMatrix &localMatrix) const
 compute M^-1 * localMatrix
 
bool entityHasChanged (const EntityType &entity) const
 returns true if the entity has been changed
 
template<class MassCaller , class BasisFunctionSet , class LocalFunction >
void applyInverseDefault (MassCaller &caller, const EntityType &entity, const Geometry &geo, const BasisFunctionSet &basisFunctionSet, LocalFunction &lf) const
 
template<class LocalMatrix >
void rightMultiplyInverseDefault (const EntityType &entity, const Geometry &geo, LocalMatrix &localMatrix) const
 compute localMatrix * M^-1
 
template<class LocalMatrix >
void leftMultiplyInverseDefault (const EntityType &entity, const Geometry &geo, LocalMatrix &localMatrix) const
 compute M^-1 * localMatrix
 
template<class BasisFunctionSet , class LocalFunction >
void applyInverseLocally (const EntityType &entity, const Geometry &geo, const BasisFunctionSet &basisFunctionSet, LocalFunction &lf) const
 apply local mass matrix to local function lf
 
template<class LocalMatrix >
void leftMultiplyInverseLocally (const EntityType &entity, const Geometry &geo, LocalMatrix &localMatrix) const
 compute M^-1 * localMatrix
 
bool setup () const
 setup and return affinity
 
template<class MassCaller , class Matrix >
void buildMatrix (MassCaller &caller, const EntityType &entity, const Geometry &geo, const BasisFunctionSetType &set, std::size_t numDofs, Matrix &matrix) const
 build local mass matrix
 
template<class Matrix >
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
 
template<class MassCallerType , class Matrix >
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
 

Detailed Description

template<class DiscreteFunctionSpace, class VolumeQuadrature>
class Dune::Fem::LocalMassMatrix< DiscreteFunctionSpace, VolumeQuadrature >

Local Mass Matrix for arbitrary spaces.

Member Function Documentation

◆ applyInverse() [1/2]

template<class DiscreteFunctionSpace , class VolumeQuadrature , bool refElemScaling = true>
template<class MassCaller , class BasisFunctionSet , class LocalFunction >
void Dune::Fem::LocalMassMatrixImplementation< DiscreteFunctionSpace, VolumeQuadrature, refElemScaling >::applyInverse ( MassCaller &  caller,
const EntityType &  entity,
const BasisFunctionSet basisFunctionSet,
LocalFunction lf 
) const
inlineinherited

◆ applyInverse() [2/2]

template<class DiscreteFunctionSpace , class VolumeQuadrature , bool refElemScaling = true>
template<class MassCaller , class LocalFunction >
void Dune::Fem::LocalMassMatrixImplementation< DiscreteFunctionSpace, VolumeQuadrature, refElemScaling >::applyInverse ( MassCaller &  caller,
const EntityType &  entity,
LocalFunction lf 
) const
inlineinherited

◆ applyInverseDefault()

template<class DiscreteFunctionSpace , class VolumeQuadrature , bool refElemScaling = true>
template<class MassCaller , class BasisFunctionSet , class LocalFunction >
void Dune::Fem::LocalMassMatrixImplementation< DiscreteFunctionSpace, VolumeQuadrature, refElemScaling >::applyInverseDefault ( MassCaller &  caller,
const EntityType &  entity,
const Geometry &  geo,
const BasisFunctionSet basisFunctionSet,
LocalFunction lf 
) const
inlineprotectedinherited

The documentation for this class was generated from the following file:
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 13, 23:29, 2024)