dune-localfunctions  2.6-git
pk1d.hh
Go to the documentation of this file.
1 // -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2 // vi: set et ts=4 sw=2 sts=2:
3 #ifndef DUNE_PK1DLOCALFINITEELEMENT_HH
4 #define DUNE_PK1DLOCALFINITEELEMENT_HH
5 
6 #include <cstddef>
7 
8 #include <dune/geometry/type.hh>
9 
12 #include "pk1d/pk1dlocalbasis.hh"
15 
16 namespace Dune
17 {
18 
21  template<class D, class R, unsigned int k>
23  {
24  public:
30 
34  {}
35 
38  Pk1DLocalFiniteElement (int variant) :
39  coefficients(variant)
40  {}
41 
46  Pk1DLocalFiniteElement (const unsigned int vertexmap[3]) :
47  coefficients(vertexmap)
48  {}
49 
52  const typename Traits::LocalBasisType& localBasis () const
53  {
54  return basis;
55  }
56 
60  {
61  return coefficients;
62  }
63 
67  {
68  return interpolation;
69  }
70 
72  unsigned int size () const
73  {
74  return basis.size();
75  }
76 
79  static constexpr GeometryType type ()
80  {
81  return GeometryTypes::line;
82  }
83 
84  private:
86  Pk1DLocalCoefficients<k> coefficients;
88  };
89 
91 
98  template<class Geometry, class RF, std::size_t k>
100  typedef typename Geometry::ctype DF;
103 
104  public:
108  struct Traits {
112  typename Basis::Traits
115  };
116 
117  private:
118  static const GeometryType gt;
119  static const LocalBasis localBasis;
120  static const LocalInterpolation localInterpolation;
121 
122  typename Traits::Basis basis_;
123  typename Traits::Interpolation interpolation_;
124  typename Traits::Coefficients coefficients_;
125 
126  public:
128 
141  template<class VertexOrder>
142  Pk1DFiniteElement(const Geometry &geometry,
143  const VertexOrder& vertexOrder) :
144  basis_(localBasis, geometry), interpolation_(localInterpolation),
145  coefficients_(vertexOrder.begin(0, 0))
146  { }
147 
148  const typename Traits::Basis& basis() const { return basis_; }
149  const typename Traits::Interpolation& interpolation() const
150  { return interpolation_; }
151  const typename Traits::Coefficients& coefficients() const
152  { return coefficients_; }
153  const GeometryType &type() const { return gt; }
154  };
155 
156  template<class Geometry, class RF, std::size_t k>
157  const GeometryType
158  Pk1DFiniteElement<Geometry, RF, k>::gt(GeometryTypes::simplex(2));
159 
160  template<class Geometry, class RF, std::size_t k>
161  const typename Pk1DFiniteElement<Geometry, RF, k>::LocalBasis
162  Pk1DFiniteElement<Geometry, RF, k>::localBasis = LocalBasis();
163 
164  template<class Geometry, class RF, std::size_t k>
165  const typename Pk1DFiniteElement<Geometry, RF, k>::LocalInterpolation
166  Pk1DFiniteElement<Geometry, RF, k>::localInterpolation =
167  LocalInterpolation();
168 
170 
180  template<class Geometry, class RF, std::size_t k>
183 
185 
199  template<class VertexOrder>
200  const FiniteElement make(const Geometry& geometry,
201  const VertexOrder& vertexOrder)
202  { return FiniteElement(geometry, vertexOrder); }
203  };
204 }
205 
206 #endif
pk1dlocalbasis.hh
Dune::Pk1DLocalFiniteElement::localInterpolation
const Traits::LocalInterpolationType & localInterpolation() const
Definition: pk1d.hh:66
Dune::Pk1DFiniteElement::Traits
Definition: pk1d.hh:108
Dune::Pk1DLocalFiniteElement::Traits
LocalFiniteElementTraits< Pk1DLocalBasis< D, R, k >, Pk1DLocalCoefficients< k >, Pk1DLocalInterpolation< Pk1DLocalBasis< D, R, k > > > Traits
Definition: pk1d.hh:29
Dune::Pk1DLocalFiniteElement
Definition: pk1d.hh:22
Dune::ScalarLocalToGlobalBasisAdaptor
Convert a simple scalar local basis into a global basis.
Definition: localtoglobaladaptors.hh:63
Dune::Pk1DLocalBasis
Lagrange shape functions of arbitrary order on the 1D reference triangle.
Definition: pk1dlocalbasis.hh:25
Dune::LocalFiniteElementTraits
traits helper struct
Definition: localfiniteelementtraits.hh:10
Dune::Pk1DFiniteElement::coefficients
const Traits::Coefficients & coefficients() const
Definition: pk1d.hh:151
Dune::Pk1DFiniteElement::Traits::Basis
ScalarLocalToGlobalBasisAdaptor< LocalBasis, Geometry > Basis
Definition: pk1d.hh:109
Dune::LocalFiniteElementTraits::LocalInterpolationType
LI LocalInterpolationType
Definition: localfiniteelementtraits.hh:22
Dune::Pk1DLocalInterpolation
Definition: pk1dlocalinterpolation.hh:11
pk1dlocalinterpolation.hh
Dune::Pk1DLocalFiniteElement::Pk1DLocalFiniteElement
Pk1DLocalFiniteElement(const unsigned int vertexmap[3])
Definition: pk1d.hh:46
localtoglobaladaptors.hh
Dune::Pk1DFiniteElement::Traits::Coefficients
Pk1DLocalCoefficients< k > Coefficients
Definition: pk1d.hh:114
Dune::LocalFiniteElementTraits::LocalCoefficientsType
LC LocalCoefficientsType
Definition: localfiniteelementtraits.hh:18
Dune::Pk1DLocalFiniteElement::localCoefficients
const Traits::LocalCoefficientsType & localCoefficients() const
Definition: pk1d.hh:59
Dune::Pk1DLocalFiniteElement::size
unsigned int size() const
Number of shape functions in this finite element.
Definition: pk1d.hh:72
Dune::Pk1DLocalFiniteElement::localBasis
const Traits::LocalBasisType & localBasis() const
Definition: pk1d.hh:52
Dune::Pk1DFiniteElement::type
const GeometryType & type() const
Definition: pk1d.hh:153
Dune::Pk1DFiniteElement
Langrange finite element of arbitrary order on triangles.
Definition: pk1d.hh:99
Dune::LocalToGlobalInterpolationAdaptor
Convert a local interpolation into a global interpolation.
Definition: localtoglobaladaptors.hh:147
pk1dlocalcoefficients.hh
Dune::Pk1DLocalFiniteElement::Pk1DLocalFiniteElement
Pk1DLocalFiniteElement()
Definition: pk1d.hh:33
localfiniteelementtraits.hh
Dune::Pk1DFiniteElement::Pk1DFiniteElement
Pk1DFiniteElement(const Geometry &geometry, const VertexOrder &vertexOrder)
construct a Pk1DFiniteElement
Definition: pk1d.hh:142
Dune::Pk1DFiniteElementFactory
Factory for Pk1DFiniteElement objects.
Definition: pk1d.hh:181
Dune::FiniteElementFactoryInterface::FiniteElement
ImplementationDefined FiniteElement
Type of the finite element.
Definition: interface.hh:115
Dune::Pk1DFiniteElement::Traits::Interpolation
LocalToGlobalInterpolationAdaptor< LocalInterpolation, typename Basis::Traits > Interpolation
Definition: pk1d.hh:113
Dune::Pk1DFiniteElement::basis
const Traits::Basis & basis() const
Definition: pk1d.hh:148
Dune::Pk1DFiniteElementFactory::FiniteElement
Pk1DFiniteElement< Geometry, RF, k > FiniteElement
Definition: pk1d.hh:182
Dune::Pk1DLocalFiniteElement::Pk1DLocalFiniteElement
Pk1DLocalFiniteElement(int variant)
Definition: pk1d.hh:38
Dune::Pk1DLocalFiniteElement::type
static constexpr GeometryType type()
Definition: pk1d.hh:79
Dune::Pk1DLocalCoefficients
Layout map for Pk elements.
Definition: pk1dlocalcoefficients.hh:22
Dune::LocalFiniteElementTraits::LocalBasisType
LB LocalBasisType
Definition: localfiniteelementtraits.hh:14
Dune
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:15
Dune::Pk1DFiniteElement::interpolation
const Traits::Interpolation & interpolation() const
Definition: pk1d.hh:149
Dune::Pk1DFiniteElementFactory::make
const FiniteElement make(const Geometry &geometry, const VertexOrder &vertexOrder)
construct Pk1DFiniteElementFactory
Definition: pk1d.hh:200
Dune::ScalarLocalToGlobalBasisAdaptor::Traits
LocalToGlobalBasisAdaptorTraits< typename LocalBasis::Traits, Geometry::coorddimension > Traits
Definition: localtoglobaladaptors.hh:82