dune-fem  2.4.1-rc
discontinuousgalerkin/interpolation.hh
Go to the documentation of this file.
1 #ifndef DUNE_FEM_SPACE_DISCONTINUOUSGALERKIN_INTERPOLATION_HH
2 #define DUNE_FEM_SPACE_DISCONTINUOUSGALERKIN_INTERPOLATION_HH
3 
4 #include <type_traits>
5 #include <utility>
6 
7 #include <dune/grid/common/capabilities.hh>
8 
12 
13 namespace Dune
14 {
15 
16  namespace Fem
17  {
18 
19  // Internal forward declaration
20  // ----------------------------
21 
22  template< class GridPart, class BasisFunctionSet,
23  class Quadrature = CachingQuadrature< GridPart, BasisFunctionSet::EntityType::codimension > >
25 
26 
27 
28  // DiscontinuousGalerkinLocalL2Projection
29  // --------------------------------------
30 
31  template< class GridPart, class BasisFunctionSet, class Quadrature >
33  : public LocalL2Projection< BasisFunctionSet, DiscontinuousGalerkinLocalL2Projection< GridPart, BasisFunctionSet, Quadrature > >
34  {
37 
38  public:
41 
42  private:
43  typedef GridPart GridPartType;
44  typedef typename GridPartType::GridType GridType;
45  static const bool cartesian = Dune::Capabilities::isCartesian< GridType >::v;
46 
47  typedef typename std::conditional< cartesian,
50  >::type LocalRieszProjectionType;
51 
53 
54  public:
59  explicit DiscontinuousGalerkinLocalL2Projection ( const BasisFunctionSetType &basisFunctionSet )
60  : impl_( LocalRieszProjectionType( basisFunctionSet ) )
61  {}
62 
63  explicit DiscontinuousGalerkinLocalL2Projection ( BasisFunctionSetType &&basisFunctionSet )
64  : impl_( LocalRieszProjectionType( std::forward< BasisFunctionSetType >( basisFunctionSet ) ) )
65  {}
66 
73  DiscontinuousGalerkinLocalL2Projection ( const ThisType & ) = default;
74 
76  : impl_( std::move( other.impl_ ) )
77  {}
78 
79  ThisType &operator= ( const ThisType & ) = default;
80 
81  ThisType &operator= ( ThisType &&other )
82  {
83  impl_ = std::move( other.impl_ );
84  return *this;
85  }
86 
94  BasisFunctionSet basisFunctionSet () const
95  {
96  return impl_.basisFunctionSet();
97  }
98 
100  template< class LocalFunction, class LocalDofVector >
101  void apply ( const LocalFunction &localFunction, LocalDofVector &localDofVector ) const
102  {
103  impl_( localFunction, localDofVector );
104  }
105 
108  private:
109  Implementation impl_;
110  };
111 
112  } // namespace Fem
113 
114 } // namespace Dune
115 
116 #endif // #ifndef DUNE_FEM_SPACE_DISCONTINUOUSGALERKIN_INTERPOLATION_HH
BaseType::BasisFunctionSetType BasisFunctionSetType
basis function set type
Definition: discontinuousgalerkin/interpolation.hh:40
DiscontinuousGalerkinLocalL2Projection(BasisFunctionSetType &&basisFunctionSet)
Definition: discontinuousgalerkin/interpolation.hh:63
void apply(const LocalFunction &localFunction, LocalDofVector &localDofVector) const
please doc me
Definition: discontinuousgalerkin/interpolation.hh:101
interface for local functions
Definition: localfunction.hh:41
Definition: coordinate.hh:4
DiscontinuousGalerkinLocalL2Projection(const BasisFunctionSetType &basisFunctionSet)
Definition: discontinuousgalerkin/interpolation.hh:59
BasisFunctionSet basisFunctionSet() const
return basis function set
Definition: discontinuousgalerkin/interpolation.hh:94
Definition: operator/projection/local/riesz/orthonormal.hh:23
STL namespace.
Definition: discontinuousgalerkin/interpolation.hh:24
void move(ArrayInterface< T > &array, const unsigned int oldOffset, const unsigned int newOffset, const unsigned int length)
Definition: array_inline.hh:38
please doc me
Definition: local/l2projection.hh:27
ThisType & operator=(const ThisType &)=default
BasisFunctionSetType basisFunctionSet() const
return basis function set
Definition: local/l2projection.hh:175
DiscontinuousGalerkinLocalL2Projection(ThisType &&other)
Definition: discontinuousgalerkin/interpolation.hh:75
Definition: basisfunctionset/basisfunctionset.hh:31