dune-localfunctions  2.3.1
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 <bitset>
7 #include <vector>
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:
30  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  BDM2Simplex2DLocalBasis(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 12;
55  }
56 
63  inline void evaluateFunction(const typename Traits::DomainType& in,
64  std::vector<typename Traits::RangeType>& out) const
65  {
66  out.resize(size());
67 
68  out[0][0] = sign_[0]*(-2*in[0]*in[1] + in[0]*in[0]);
69  out[0][1] = sign_[0]*(-1 + 6*in[1] -2*in[0]*in[1] - 5*in[1]*in[1]);
70 
71  out[1][0] = 1.5*in[0] + 3*in[0]*in[1] - 4.5*in[0]*in[0];
72  out[1][1] = -3 + 6*in[0] + 10.5*in[1] - 15*in[0]*in[1] - 7.5*in[1]*in[1];
73 
74  out[2][0] = sign_[0]*(-7.5*in[0] + 5*in[0]*in[1] + 12.5*in[0]*in[0]);
75  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]);
76 
77 
78 
79  out[3][0] = sign_[1]*(-1 + 6*in[0] - 2*in[0]*in[1] - 5*in[0]*in[0]);
80  out[3][1] = sign_[1]*(-2*in[0]*in[1] + in[1]*in[1]);
81 
82  out[4][0] = 3 - 10.5*in[0] - 6*in[1] + 15*in[0]*in[1] + 7.5*in[0]*in[0];
83  out[4][1] = -1.5*in[1] - 3*in[0]*in[1] + 4.5*in[1]*in[1];
84 
85  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]);
86  out[5][1] = sign_[1]*(-7.5*in[1] + 5*in[0]*in[1] + 12.5*in[1]*in[1]);
87 
88 
89 
90  out[6][0] = sign_[2]*(-3*in[0] + 4*in[0]*in[1] + 4*in[0]*in[0]);
91  out[6][1] = sign_[2]*(-3*in[1] + 4*in[0]*in[1] + 4*in[1]*in[1]);
92 
93  out[7][0] = -3*in[0] + 6*in[0]*in[0];
94  out[7][1] = 3*in[1] - 6*in[1]*in[1];
95 
96  out[8][0] = sign_[2]*(-10*in[0]*in[1] + 5*in[0]*in[0]);
97  out[8][1] = sign_[2]*(-10*in[0]*in[1] + 5*in[1]*in[1]);
98 
99 
100 
101  out[9][0] = 18*in[0] - 12*in[0]*in[1] - 18*in[0]*in[0];
102  out[9][1] = 6*in[1] - 12*in[0]*in[1] - 6*in[1]*in[1];
103 
104  out[10][0] = 6*in[0] - 12*in[0]*in[1] - 6*in[0]*in[0];
105  out[10][1] = 18*in[1] - 12*in[0]*in[1] - 18*in[1]*in[1];
106 
107  out[11][0] = 90*in[0] - 180*in[0]*in[1] - 90*in[0]*in[0];
108  out[11][1] = -90*in[1] + 180*in[0]*in[1] + 90*in[1]*in[1];
109  }
110 
117  inline void evaluateJacobian(const typename Traits::DomainType& in,
118  std::vector<typename Traits::JacobianType>& out) const
119  {
120  out.resize(size());
121 
122  out[0][0][0] = sign_[0]*(-2*in[1] + 2*in[0]);
123  out[0][0][1] = sign_[0]*(-2*in[0]);
124 
125  out[0][1][0] = sign_[0]*(-2*in[1]);
126  out[0][1][1] = sign_[0]*(6 -2*in[0] - 10*in[1]);
127 
128 
129  out[1][0][0] = 1.5 + 3*in[1] - 9*in[0];
130  out[1][0][1] = 3*in[0];
131 
132  out[1][1][0] = 6 - 15*in[1];
133  out[1][1][1] = 10.5 - 15*in[0] - 15*in[1];
134 
135 
136  out[2][0][0] = sign_[0]*(-7.5 + 5*in[1] + 25*in[0]);
137  out[2][0][1] = sign_[0]*(5*in[0]);
138 
139  out[2][1][0] = sign_[0]*(30 - 25*in[1] - 60*in[0]);
140  out[2][1][1] = sign_[0]*(7.5 - 25*in[0] - 5*in[1]);
141 
142 
143 
144  out[3][0][0] = sign_[1]*(6 - 2*in[1] - 10*in[0]);
145  out[3][0][1] = sign_[1]*(-2*in[0]);
146 
147  out[3][1][0] = sign_[1]*(-2*in[1]);
148  out[3][1][1] = sign_[1]*(-2*in[0] + 2*in[1]);
149 
150 
151  out[4][0][0] = -10.5 + 15*in[1] + 15*in[0];
152  out[4][0][1] = -6 + 15*in[0];
153 
154  out[4][1][0] = -3*in[1];
155  out[4][1][1] = -1.5 - 3*in[0] + 9*in[1];
156 
157 
158  out[5][0][0] = sign_[1]*(7.5 - 25*in[1] - 5*in[0]);
159  out[5][0][1] = sign_[1]*(30 - 25*in[0] - 60*in[1]);
160 
161  out[5][1][0] = sign_[1]*(5*in[1]);
162  out[5][1][1] = sign_[1]*(-7.5 + 5*in[0] + 25*in[1]);
163 
164 
165 
166  out[6][0][0] = sign_[2]*(-3 + 4*in[1] + 8*in[0]);
167  out[6][0][1] = sign_[2]*(4*in[0]);
168 
169  out[6][1][0] = sign_[2]*(4*in[1]);
170  out[6][1][1] = sign_[2]*(-3 + 4*in[0] + 8*in[1]);
171 
172 
173  out[7][0][0] = -3 + 12*in[0];
174  out[7][0][1] = 0;
175 
176  out[7][1][0] = 0;
177  out[7][1][1] = 3 - 12*in[1];
178 
179 
180  out[8][0][0] = sign_[2]*(-10*in[1] + 10*in[0]);
181  out[8][0][1] = sign_[2]*(-10*in[0]);
182 
183  out[8][1][0] = sign_[2]*(-10*in[1]);
184  out[8][1][1] = sign_[2]*(-10*in[0] + 10*in[1]);
185 
186 
187  out[9][0][0] = 18 - 12*in[1] - 36*in[0];
188  out[9][0][1] = -12*in[0];
189 
190  out[9][1][0] = -12*in[1];
191  out[9][1][1] = 6 - 12*in[0] - 12*in[1];
192 
193  out[10][0][0] = 6 - 12*in[1] - 12*in[0];
194  out[10][0][1] = -12*in[0];
195 
196  out[10][1][0] = -12*in[1];
197  out[10][1][1] = 18 - 12*in[0] - 36*in[1];
198 
199  out[11][0][0] = 90 - 180*in[1] - 180*in[0];
200  out[11][0][1] = -180*in[0];
201 
202  out[11][1][0] = 180*in[1];
203  out[11][1][1] = -90 + 180*in[0] + 180*in[1];
204  }
205 
207  unsigned int order() const
208  {
209  return 2; // TODO: check whether this is not order 3
210  }
211 
212  private:
213  array<R,3> sign_;
214  };
215 } // end namespace Dune
216 #endif // DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI2_SIMPLEX2D_LOCALBASIS_HH