dune-localfunctions  2.3.1-rc1
brezzidouglasmarini1simplex2dlocalbasis.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_LOCALBASIS_HH
4 #define DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI1_SIMPLEX2D_LOCALBASIS_HH
5 
6 #include <vector>
7 #include <bitset>
8 
9 #include <dune/common/fmatrix.hh>
10 
11 #include "../../common/localbasis.hh"
12 
13 namespace Dune
14 {
24  template<class D, class R>
26  {
27 
28  public:
29  typedef LocalBasisTraits<D,2,Dune::FieldVector<D,2>,R,2,Dune::FieldVector<R,2>,
30  Dune::FieldMatrix<R,2,2> > Traits;
31 
34  {
35  for (size_t i=0; i<3; i++)
36  sign_[i] = 1.0;
37  }
38 
44  BDM1Simplex2DLocalBasis (std::bitset<3> s)
45  {
46  for (size_t i=0; i<3; i++)
47  sign_[i] = s[i] ? -1.0 : 1.0;
48  }
49 
51  unsigned int size () const
52  {
53  return 6;
54  }
55 
62  inline void evaluateFunction (const typename Traits::DomainType& in,
63  std::vector<typename Traits::RangeType>& out) const
64  {
65  out.resize(6);
66 
67  out[0][0] = sign_[0]*in[0];
68  out[0][1] = sign_[0]*(in[1] - 1.0);
69  out[1][0] = sign_[1]*(in[0] - 1.0);
70  out[1][1] = sign_[1]*in[1];
71  out[2][0] = sign_[2]*in[0];
72  out[2][1] = sign_[2]*in[1];
73  out[3][0] = 3.0*in[0];
74  out[3][1] = 3.0 - 6.0*in[0] - 3.0*in[1];
75  out[4][0] = -3.0 + 3.0*in[0] + 6.0*in[1];
76  out[4][1] = -3.0*in[1];
77  out[5][0] = -3.0*in[0];
78  out[5][1] = 3.0*in[1];
79  }
80 
87  inline void evaluateJacobian (const typename Traits::DomainType& in,
88  std::vector<typename Traits::JacobianType>& out) const
89  {
90  out.resize(6);
91 
92  out[0][0][0] = sign_[0];
93  out[0][0][1] = 0.0;
94  out[0][1][0] = 0.0;
95  out[0][1][1] = sign_[0];
96 
97  out[1][0][0] = sign_[1];
98  out[1][0][1] = 0.0;
99  out[1][1][0] = 0.0;
100  out[1][1][1] = sign_[1];
101 
102  out[2][0][0] = sign_[2];
103  out[2][0][1] = 0.0;
104  out[2][1][0] = 0.0;
105  out[2][1][1] = sign_[2];
106 
107  out[3][0][0] = 3.0;
108  out[3][0][1] = 0.0;
109  out[3][1][0] = -6.0;
110  out[3][1][1] = -3.0;
111 
112  out[4][0][0] = 3.0;
113  out[4][0][1] = 6.0;
114  out[4][1][0] = 0.0;
115  out[4][1][1] = -3.0;
116 
117  out[5][0][0] = -3.0;
118  out[5][0][1] = 0.0;
119  out[5][1][0] = 0.0;
120  out[5][1][1] = 3.0;
121  }
122 
124  unsigned int order () const
125  {
126  return 1;
127  }
128 
129  private:
130  array<R,3> sign_;
131  };
132 }
133 #endif // DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI1_SIMPLEX2D_LOCALBASIS_HH
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: brezzidouglasmarini1simplex2dlocalbasis.hh:87
First order Brezzi-Douglas-Marini shape functions on the reference triangle.
Definition: brezzidouglasmarini1simplex2dlocalbasis.hh:25
Type traits for LocalBasisVirtualInterface.
Definition: localbasis.hh:39
BDM1Simplex2DLocalBasis()
Standard constructor.
Definition: brezzidouglasmarini1simplex2dlocalbasis.hh:33
unsigned int order() const
Polynomial order of the shape functions.
Definition: brezzidouglasmarini1simplex2dlocalbasis.hh:124
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: brezzidouglasmarini1simplex2dlocalbasis.hh:62
unsigned int size() const
number of shape functions
Definition: brezzidouglasmarini1simplex2dlocalbasis.hh:51
LocalBasisTraits< D, 2, Dune::FieldVector< D, 2 >, R, 2, Dune::FieldVector< R, 2 >, Dune::FieldMatrix< R, 2, 2 > > Traits
Definition: brezzidouglasmarini1simplex2dlocalbasis.hh:30
BDM1Simplex2DLocalBasis(std::bitset< 3 > s)
Make set number s, where 0 <= s < 8.
Definition: brezzidouglasmarini1simplex2dlocalbasis.hh:44
D DomainType
domain type
Definition: localbasis.hh:51