dune-localfunctions  2.5-git
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  gt.makeTriangle();
38  }
39 
45  BDM1Simplex2DLocalFiniteElement (int s) : basis(s), interpolation(s)
46  {
47  gt.makeTriangle();
48  }
49 
50  const typename Traits::LocalBasisType& localBasis () const
51  {
52  return basis;
53  }
54 
56  {
57  return coefficients;
58  }
59 
61  {
62  return interpolation;
63  }
64 
66  unsigned int size () const
67  {
68  return basis.size();
69  }
70 
71  GeometryType type () const
72  {
73  return gt;
74  }
75 
76  private:
77  BDM1Simplex2DLocalBasis<D,R> basis;
78  BDM1Simplex2DLocalCoefficients coefficients;
80  GeometryType gt;
81  };
82 }
83 #endif // DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI1_SIMPLEX2D_LOCALFINITEELEMENT_HH
First order Brezzi-Douglas-Marini shape functions on triangles.
Definition: brezzidouglasmarini1simplex2d.hh:25
LI LocalInterpolationType
Definition: localfiniteelementtraits.hh:22
const Traits::LocalBasisType & localBasis() const
Definition: brezzidouglasmarini1simplex2d.hh:50
GeometryType type() const
Definition: brezzidouglasmarini1simplex2d.hh:71
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:15
BDM1Simplex2DLocalFiniteElement()
Standard constructor.
Definition: brezzidouglasmarini1simplex2d.hh:35
traits helper struct
Definition: localfiniteelementtraits.hh:10
First order Brezzi-Douglas-Marini shape functions on the reference triangle.
Definition: brezzidouglasmarini1simplex2dlocalbasis.hh:27
const Traits::LocalInterpolationType & localInterpolation() const
Definition: brezzidouglasmarini1simplex2d.hh:60
LocalFiniteElementTraits< BDM1Simplex2DLocalBasis< D, R >, BDM1Simplex2DLocalCoefficients, BDM1Simplex2DLocalInterpolation< BDM1Simplex2DLocalBasis< D, R > > > Traits
Definition: brezzidouglasmarini1simplex2d.hh:32
const Traits::LocalCoefficientsType & localCoefficients() const
Definition: brezzidouglasmarini1simplex2d.hh:55
unsigned int size() const
Number of shape functions in this finite element.
Definition: brezzidouglasmarini1simplex2d.hh:66
LC LocalCoefficientsType
Definition: localfiniteelementtraits.hh:18
First order Brezzi-Douglas-Marini shape functions on the reference triangle.
Definition: brezzidouglasmarini1simplex2dlocalinterpolation.hh:21
LB LocalBasisType
Definition: localfiniteelementtraits.hh:14
BDM1Simplex2DLocalFiniteElement(int s)
Make set number s, where 0 <= s < 8.
Definition: brezzidouglasmarini1simplex2d.hh:45
Layout map for Brezzi-Douglas-Marini-1 elements on triangles.
Definition: brezzidouglasmarini1simplex2dlocalcoefficients.hh:21