dune-localfunctions  2.7.0
brezzidouglasmarini1simplex2d.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_LOCALFUNCTIONS_BREZZIDOUGLASMARINI1_SIMPLEX2D_LOCALFINITEELEMENT_HH
4 #define DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI1_SIMPLEX2D_LOCALFINITEELEMENT_HH
5 
6 #include <dune/geometry/type.hh>
7 
8 #include "../common/localfiniteelementtraits.hh"
12 
13 namespace Dune
14 {
15 
24  template<class D, class R>
26  {
27 
28  public:
33 
36  {}
37 
44  basis(s),
45  interpolation(s)
46  {}
47 
48  const typename Traits::LocalBasisType& localBasis () const
49  {
50  return basis;
51  }
52 
54  {
55  return coefficients;
56  }
57 
59  {
60  return interpolation;
61  }
62 
64  unsigned int size () const
65  {
66  return basis.size();
67  }
68 
69  static constexpr GeometryType type ()
70  {
71  return GeometryTypes::triangle;
72  }
73 
74  private:
76  BDM1Simplex2DLocalCoefficients coefficients;
78  };
79 }
80 #endif // DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI1_SIMPLEX2D_LOCALFINITEELEMENT_HH
brezzidouglasmarini1simplex2dlocalbasis.hh
Dune::BDM1Simplex2DLocalFiniteElement
First order Brezzi-Douglas-Marini shape functions on triangles.
Definition: brezzidouglasmarini1simplex2d.hh:25
Dune::BDM1Simplex2DLocalFiniteElement::BDM1Simplex2DLocalFiniteElement
BDM1Simplex2DLocalFiniteElement(int s)
Make set number s, where 0 <= s < 8.
Definition: brezzidouglasmarini1simplex2d.hh:43
Dune::BDM1Simplex2DLocalFiniteElement::type
static constexpr GeometryType type()
Definition: brezzidouglasmarini1simplex2d.hh:69
Dune
Definition: bdfmcube.hh:15
Dune::BDM1Simplex2DLocalInterpolation
First order Brezzi-Douglas-Marini shape functions on the reference triangle.
Definition: brezzidouglasmarini1simplex2dlocalinterpolation.hh:22
Dune::BDM1Simplex2DLocalBasis
First order Brezzi-Douglas-Marini shape functions on the reference triangle.
Definition: brezzidouglasmarini1simplex2dlocalbasis.hh:27
Dune::BDM1Simplex2DLocalFiniteElement::localCoefficients
const Traits::LocalCoefficientsType & localCoefficients() const
Definition: brezzidouglasmarini1simplex2d.hh:53
Dune::LocalFiniteElementTraits
traits helper struct
Definition: localfiniteelementtraits.hh:10
Dune::BDM1Simplex2DLocalCoefficients
Layout map for Brezzi-Douglas-Marini-1 elements on triangles.
Definition: brezzidouglasmarini1simplex2dlocalcoefficients.hh:21
Dune::LocalFiniteElementTraits::LocalInterpolationType
LI LocalInterpolationType
Definition: localfiniteelementtraits.hh:22
Dune::LocalFiniteElementTraits::LocalCoefficientsType
LC LocalCoefficientsType
Definition: localfiniteelementtraits.hh:18
Dune::BDM1Simplex2DLocalFiniteElement::localBasis
const Traits::LocalBasisType & localBasis() const
Definition: brezzidouglasmarini1simplex2d.hh:48
Dune::BDM1Simplex2DLocalFiniteElement::size
unsigned int size() const
Number of shape functions in this finite element.
Definition: brezzidouglasmarini1simplex2d.hh:64
brezzidouglasmarini1simplex2dlocalcoefficients.hh
brezzidouglasmarini1simplex2dlocalinterpolation.hh
Dune::BDM1Simplex2DLocalFiniteElement::localInterpolation
const Traits::LocalInterpolationType & localInterpolation() const
Definition: brezzidouglasmarini1simplex2d.hh:58
Dune::BDM1Simplex2DLocalFiniteElement::BDM1Simplex2DLocalFiniteElement
BDM1Simplex2DLocalFiniteElement()
Standard constructor.
Definition: brezzidouglasmarini1simplex2d.hh:35
Dune::BDM1Simplex2DLocalFiniteElement::Traits
LocalFiniteElementTraits< BDM1Simplex2DLocalBasis< D, R >, BDM1Simplex2DLocalCoefficients, BDM1Simplex2DLocalInterpolation< BDM1Simplex2DLocalBasis< D, R > > > Traits
Definition: brezzidouglasmarini1simplex2d.hh:32
Dune::LocalFiniteElementTraits::LocalBasisType
LB LocalBasisType
Definition: localfiniteelementtraits.hh:14