DUNE-FEM (unstable)

localinterpolation.hh
Go to the documentation of this file.
1#ifndef DUNE_FEM_SPACE_DISCONTINUOUSGALERKIN_LOCALINTERPOLATION_HH
2#define DUNE_FEM_SPACE_DISCONTINUOUSGALERKIN_LOCALINTERPOLATION_HH
3
4// dune-fem includes
6#include <dune/fem/operator/1order/localmassmatrix.hh>
7#include <dune/fem/quadrature/cachingquadrature.hh>
8#include <dune/fem/quadrature/agglomerationquadrature.hh>
9
16namespace Dune
17{
18
19 namespace Fem
20 {
21
22 // DiscontinuousGalerkinLocalInterpolation
23 // ---------------------------------------
24
28 template< class DiscreteFunctionSpace >
30 {
32
33 public:
35 typedef typename DiscreteFunctionSpaceType::GridType GridType;
36 typedef typename DiscreteFunctionSpaceType::EntityType EntityType;
37
38 static const bool isAlwaysAffine = Dune::Capabilities::isCartesian< GridType >::v ||
40 // always true for orthonormal spaces
41 //static const bool isAlwaysAffine = true;
42
43 private:
44 typedef typename DiscreteFunctionSpaceType::RangeType RangeType;
45 typedef typename RangeType::value_type RangeFieldType;
46
47 typedef typename DiscreteFunctionSpaceType::GridPartType GridPartType;
48
49 typedef typename DiscreteFunctionSpaceType :: LocalMassMatrixType
51
52 typedef typename LocalMassMatrixType::VolumeQuadratureType QuadratureType;
53
54 public:
56 : space_( space )
57 {}
58
59 DiscontinuousGalerkinLocalInterpolation ( const ThisType &other ) = default;
61
62 void bind( const EntityType& ) {}
63 void unbind() {}
64
65 ThisType &operator= ( const ThisType &other ) = delete;
66
67 template< class LocalFunction, class LocalDofVector >
68 void operator () ( const LocalFunction &localFunction, LocalDofVector &dofs ) const
69 {
70 // set all dofs to zero
71 std::fill( dofs.begin(), dofs.end(), typename LocalDofVector::value_type(0) );
72
73 // get entity and geometry
74 const EntityType &entity = localFunction.entity();
75
76 if( entity.type().isNone() )
77 {
79 ElementQuadratureType quadrature( entity, localFunction.order() + space_.order( entity ) );
80 bool isAffine = computeInterpolation( entity, quadrature, localFunction, dofs );
81 if( ! isAffine )
82 {
84 AggloMassMatrix massMat( space_, massMatrix().volumeQuadratureOrder( entity ) );
85 // apply inverse of mass matrix
86 auto basisFunctionSet = space_.basisFunctionSet(entity);
87 massMat.applyInverse( entity, basisFunctionSet, dofs );
88 }
89 }
90 else
91 {
92 QuadratureType quadrature( entity, localFunction.order() + space_.order( entity ) );
93 bool isAffine = computeInterpolation( entity, quadrature, localFunction, dofs );
94 if( ! isAffine )
95 {
96 // apply inverse of mass matrix
97 auto basisFunctionSet = space_.basisFunctionSet(entity);
98 massMatrix().applyInverse( entity, basisFunctionSet, dofs );
99 }
100 }
101
102 }
103
104 private:
105 template<class QuadImpl, class LocalFunction, class LocalDofVector >
106 bool computeInterpolation( const EntityType& entity,
107 const QuadImpl& quadrature,
108 const LocalFunction &localFunction,
109 LocalDofVector &dofs ) const
110 {
111 const int nop = quadrature.nop();
112
113 auto& values = space_.localMassMatrixStorage().second;
114
115 // adjust size of values
116 values.resize( nop );
117
118 // evaluate local function for all quadrature points
119 localFunction.evaluateQuadrature( quadrature, values );
120
121 bool isAffine = isAlwaysAffine ;
122 if( ! isAlwaysAffine )
123 {
124 const auto geometry = entity.geometry();
125 isAffine = geometry.affine();
126
127 if( ! isAffine )
128 {
129 // apply weight
130 for(auto qp : quadrature )
131 values[ qp.index() ] *= qp.weight() * geometry.integrationElement( qp.position() );
132 }
133 }
134
135 if( isAffine )
136 {
137 // apply weight only
138 for(auto qp : quadrature )
139 values[ qp.index() ] *= qp.weight();
140 }
141
142 // add values to local function
143 space_.basisFunctionSet(entity).axpy( quadrature, values, dofs );
144
145 return isAffine;
146 }
147
148 const LocalMassMatrixType &massMatrix () const
149 {
150 return space_.localMassMatrixStorage().first;
151 }
152
153 const DiscreteFunctionSpaceType& space_;
154 };
155
156 } // namespace Fem
157
158} // namespace Dune
159
160#endif // #ifndef DUNE_FEM_SPACE_DISCONTINUOUSGALERKIN_LOCALINTERPOLATION_HH
discrete function space
Definition: localinterpolation.hh:30
GridPartType::GridType GridType
type of underlying dune grid
Definition: discretefunctionspace.hh:214
quadrature on the codim-0 reference element
Definition: elementquadrature.hh:58
interface for local functions
Definition: localfunction.hh:77
void evaluateQuadrature(const QuadratureType &quad, Vectors &... vec) const
evaluate all basisfunctions for all quadrature points and store the results in the result vector
Definition: localfunction.hh:384
const EntityType & entity() const
obtain the entity, this local function lives on
Definition: localfunction.hh:302
int order() const
obtain the order of this local function
Definition: localfunction.hh:290
void applyInverse(MassCaller &caller, const EntityType &entity, const BasisFunctionSet &basisFunctionSet, LocalFunction &lf) const
Definition: localmassmatrix.hh:295
Local Mass Matrix for arbitrary spaces.
Definition: localmassmatrix.hh:919
A set of traits classes to store static information about grid implementation.
Dune namespace.
Definition: alignedallocator.hh:13
Specialize with 'true' for if the codimension 0 entity of the grid has only one possible geometry typ...
Definition: capabilities.hh:27
Specialize with 'true' if the grid is a Cartesian grid. Cartesian grids satisfy the following propert...
Definition: capabilities.hh:48
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 27, 22:29, 2024)