dune-localfunctions  2.4.1-rc2
brezzidouglasmarini2simplex2dlocalbasis.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_BREZZIDOUGLASMARINI2_SIMPLEX2D_LOCALBASIS_HH
4 #define DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI2_SIMPLEX2D_LOCALBASIS_HH
5 
6 #include <array>
7 #include <bitset>
8 #include <vector>
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:
31  R,2,Dune::FieldVector<R,2>,
32  Dune::FieldMatrix<R,2,2> > Traits;
33 
36  {
37  for (size_t i=0; i<3; i++)
38  sign_[i] = 1.0;
39  }
40 
46  BDM2Simplex2DLocalBasis(std::bitset<3> s)
47  {
48  for (size_t i=0; i<3; i++)
49  sign_[i] = s[i] ? -1.0 : 1.0;
50  }
51 
53  unsigned int size() const
54  {
55  return 12;
56  }
57 
64  inline void evaluateFunction(const typename Traits::DomainType& in,
65  std::vector<typename Traits::RangeType>& out) const
66  {
67  out.resize(size());
68 
69  out[0][0] = sign_[0]*(-2*in[0]*in[1] + in[0]*in[0]);
70  out[0][1] = sign_[0]*(-1 + 6*in[1] -2*in[0]*in[1] - 5*in[1]*in[1]);
71 
72  out[1][0] = 1.5*in[0] + 3*in[0]*in[1] - 4.5*in[0]*in[0];
73  out[1][1] = -3 + 6*in[0] + 10.5*in[1] - 15*in[0]*in[1] - 7.5*in[1]*in[1];
74 
75  out[2][0] = sign_[0]*(-7.5*in[0] + 5*in[0]*in[1] + 12.5*in[0]*in[0]);
76  out[2][1] = sign_[0]*(-5 + 30*in[0] + 7.5*in[1] - 25*in[0]*in[1] - 30*in[0]*in[0] - 2.5*in[1]*in[1]);
77 
78 
79 
80  out[3][0] = sign_[1]*(-1 + 6*in[0] - 2*in[0]*in[1] - 5*in[0]*in[0]);
81  out[3][1] = sign_[1]*(-2*in[0]*in[1] + in[1]*in[1]);
82 
83  out[4][0] = 3 - 10.5*in[0] - 6*in[1] + 15*in[0]*in[1] + 7.5*in[0]*in[0];
84  out[4][1] = -1.5*in[1] - 3*in[0]*in[1] + 4.5*in[1]*in[1];
85 
86  out[5][0] = sign_[1]*(-5 + 7.5*in[0] + 30*in[1] - 25*in[0]*in[1] - 2.5*in[0]*in[0] - 30*in[1]*in[1]);
87  out[5][1] = sign_[1]*(-7.5*in[1] + 5*in[0]*in[1] + 12.5*in[1]*in[1]);
88 
89 
90 
91  out[6][0] = sign_[2]*(-3*in[0] + 4*in[0]*in[1] + 4*in[0]*in[0]);
92  out[6][1] = sign_[2]*(-3*in[1] + 4*in[0]*in[1] + 4*in[1]*in[1]);
93 
94  out[7][0] = -3*in[0] + 6*in[0]*in[0];
95  out[7][1] = 3*in[1] - 6*in[1]*in[1];
96 
97  out[8][0] = sign_[2]*(-10*in[0]*in[1] + 5*in[0]*in[0]);
98  out[8][1] = sign_[2]*(-10*in[0]*in[1] + 5*in[1]*in[1]);
99 
100 
101 
102  out[9][0] = 18*in[0] - 12*in[0]*in[1] - 18*in[0]*in[0];
103  out[9][1] = 6*in[1] - 12*in[0]*in[1] - 6*in[1]*in[1];
104 
105  out[10][0] = 6*in[0] - 12*in[0]*in[1] - 6*in[0]*in[0];
106  out[10][1] = 18*in[1] - 12*in[0]*in[1] - 18*in[1]*in[1];
107 
108  out[11][0] = 90*in[0] - 180*in[0]*in[1] - 90*in[0]*in[0];
109  out[11][1] = -90*in[1] + 180*in[0]*in[1] + 90*in[1]*in[1];
110  }
111 
118  inline void evaluateJacobian(const typename Traits::DomainType& in,
119  std::vector<typename Traits::JacobianType>& out) const
120  {
121  out.resize(size());
122 
123  out[0][0][0] = sign_[0]*(-2*in[1] + 2*in[0]);
124  out[0][0][1] = sign_[0]*(-2*in[0]);
125 
126  out[0][1][0] = sign_[0]*(-2*in[1]);
127  out[0][1][1] = sign_[0]*(6 -2*in[0] - 10*in[1]);
128 
129 
130  out[1][0][0] = 1.5 + 3*in[1] - 9*in[0];
131  out[1][0][1] = 3*in[0];
132 
133  out[1][1][0] = 6 - 15*in[1];
134  out[1][1][1] = 10.5 - 15*in[0] - 15*in[1];
135 
136 
137  out[2][0][0] = sign_[0]*(-7.5 + 5*in[1] + 25*in[0]);
138  out[2][0][1] = sign_[0]*(5*in[0]);
139 
140  out[2][1][0] = sign_[0]*(30 - 25*in[1] - 60*in[0]);
141  out[2][1][1] = sign_[0]*(7.5 - 25*in[0] - 5*in[1]);
142 
143 
144 
145  out[3][0][0] = sign_[1]*(6 - 2*in[1] - 10*in[0]);
146  out[3][0][1] = sign_[1]*(-2*in[0]);
147 
148  out[3][1][0] = sign_[1]*(-2*in[1]);
149  out[3][1][1] = sign_[1]*(-2*in[0] + 2*in[1]);
150 
151 
152  out[4][0][0] = -10.5 + 15*in[1] + 15*in[0];
153  out[4][0][1] = -6 + 15*in[0];
154 
155  out[4][1][0] = -3*in[1];
156  out[4][1][1] = -1.5 - 3*in[0] + 9*in[1];
157 
158 
159  out[5][0][0] = sign_[1]*(7.5 - 25*in[1] - 5*in[0]);
160  out[5][0][1] = sign_[1]*(30 - 25*in[0] - 60*in[1]);
161 
162  out[5][1][0] = sign_[1]*(5*in[1]);
163  out[5][1][1] = sign_[1]*(-7.5 + 5*in[0] + 25*in[1]);
164 
165 
166 
167  out[6][0][0] = sign_[2]*(-3 + 4*in[1] + 8*in[0]);
168  out[6][0][1] = sign_[2]*(4*in[0]);
169 
170  out[6][1][0] = sign_[2]*(4*in[1]);
171  out[6][1][1] = sign_[2]*(-3 + 4*in[0] + 8*in[1]);
172 
173 
174  out[7][0][0] = -3 + 12*in[0];
175  out[7][0][1] = 0;
176 
177  out[7][1][0] = 0;
178  out[7][1][1] = 3 - 12*in[1];
179 
180 
181  out[8][0][0] = sign_[2]*(-10*in[1] + 10*in[0]);
182  out[8][0][1] = sign_[2]*(-10*in[0]);
183 
184  out[8][1][0] = sign_[2]*(-10*in[1]);
185  out[8][1][1] = sign_[2]*(-10*in[0] + 10*in[1]);
186 
187 
188  out[9][0][0] = 18 - 12*in[1] - 36*in[0];
189  out[9][0][1] = -12*in[0];
190 
191  out[9][1][0] = -12*in[1];
192  out[9][1][1] = 6 - 12*in[0] - 12*in[1];
193 
194  out[10][0][0] = 6 - 12*in[1] - 12*in[0];
195  out[10][0][1] = -12*in[0];
196 
197  out[10][1][0] = -12*in[1];
198  out[10][1][1] = 18 - 12*in[0] - 36*in[1];
199 
200  out[11][0][0] = 90 - 180*in[1] - 180*in[0];
201  out[11][0][1] = -180*in[0];
202 
203  out[11][1][0] = 180*in[1];
204  out[11][1][1] = -90 + 180*in[0] + 180*in[1];
205  }
206 
208  unsigned int order() const
209  {
210  return 2; // TODO: check whether this is not order 3
211  }
212 
213  private:
214  std::array<R,3> sign_;
215  };
216 } // end namespace Dune
217 #endif // DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI2_SIMPLEX2D_LOCALBASIS_HH
Type traits for LocalBasisVirtualInterface.
Definition: localbasis.hh:37
First order Brezzi-Douglas-Marini shape functions on quadrilaterals.
Definition: brezzidouglasmarini2simplex2dlocalbasis.hh:26
unsigned int order() const
Polynomial order of the shape functions.
Definition: brezzidouglasmarini2simplex2dlocalbasis.hh:208
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:14
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: brezzidouglasmarini2simplex2dlocalbasis.hh:118
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: brezzidouglasmarini2simplex2dlocalbasis.hh:64
BDM2Simplex2DLocalBasis(std::bitset< 3 > s)
Make set number s, where 0 <= s < 8.
Definition: brezzidouglasmarini2simplex2dlocalbasis.hh:46
BDM2Simplex2DLocalBasis()
Standard constructor.
Definition: brezzidouglasmarini2simplex2dlocalbasis.hh:35
unsigned int size() const
number of shape functions
Definition: brezzidouglasmarini2simplex2dlocalbasis.hh:53
LocalBasisTraits< D, 2, Dune::FieldVector< D, 2 >, R, 2, Dune::FieldVector< R, 2 >, Dune::FieldMatrix< R, 2, 2 > > Traits
Definition: brezzidouglasmarini2simplex2dlocalbasis.hh:32
D DomainType
domain type
Definition: localbasis.hh:49