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