3 #ifndef DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI2_SIMPLEX2D_LOCALBASIS_HH
4 #define DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI2_SIMPLEX2D_LOCALBASIS_HH
9 #include <dune/common/fmatrix.hh>
11 #include "../../common/localbasis.hh"
24 template<
class D,
class R>
30 R,2,Dune::FieldVector<R,2>,
36 for (
size_t i=0; i<3; i++)
47 for (
size_t i=0; i<3; i++)
48 sign_[i] = s[i] ? -1.0 : 1.0;
64 std::vector<typename Traits::RangeType>& out)
const
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]);
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];
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]);
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]);
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];
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]);
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]);
93 out[7][0] = -3*in[0] + 6*in[0]*in[0];
94 out[7][1] = 3*in[1] - 6*in[1]*in[1];
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]);
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];
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];
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];
118 std::vector<typename Traits::JacobianType>& out)
const
122 out[0][0][0] = sign_[0]*(-2*in[1] + 2*in[0]);
123 out[0][0][1] = sign_[0]*(-2*in[0]);
125 out[0][1][0] = sign_[0]*(-2*in[1]);
126 out[0][1][1] = sign_[0]*(6 -2*in[0] - 10*in[1]);
129 out[1][0][0] = 1.5 + 3*in[1] - 9*in[0];
130 out[1][0][1] = 3*in[0];
132 out[1][1][0] = 6 - 15*in[1];
133 out[1][1][1] = 10.5 - 15*in[0] - 15*in[1];
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]);
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]);
144 out[3][0][0] = sign_[1]*(6 - 2*in[1] - 10*in[0]);
145 out[3][0][1] = sign_[1]*(-2*in[0]);
147 out[3][1][0] = sign_[1]*(-2*in[1]);
148 out[3][1][1] = sign_[1]*(-2*in[0] + 2*in[1]);
151 out[4][0][0] = -10.5 + 15*in[1] + 15*in[0];
152 out[4][0][1] = -6 + 15*in[0];
154 out[4][1][0] = -3*in[1];
155 out[4][1][1] = -1.5 - 3*in[0] + 9*in[1];
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]);
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]);
166 out[6][0][0] = sign_[2]*(-3 + 4*in[1] + 8*in[0]);
167 out[6][0][1] = sign_[2]*(4*in[0]);
169 out[6][1][0] = sign_[2]*(4*in[1]);
170 out[6][1][1] = sign_[2]*(-3 + 4*in[0] + 8*in[1]);
173 out[7][0][0] = -3 + 12*in[0];
177 out[7][1][1] = 3 - 12*in[1];
180 out[8][0][0] = sign_[2]*(-10*in[1] + 10*in[0]);
181 out[8][0][1] = sign_[2]*(-10*in[0]);
183 out[8][1][0] = sign_[2]*(-10*in[1]);
184 out[8][1][1] = sign_[2]*(-10*in[0] + 10*in[1]);
187 out[9][0][0] = 18 - 12*in[1] - 36*in[0];
188 out[9][0][1] = -12*in[0];
190 out[9][1][0] = -12*in[1];
191 out[9][1][1] = 6 - 12*in[0] - 12*in[1];
193 out[10][0][0] = 6 - 12*in[1] - 12*in[0];
194 out[10][0][1] = -12*in[0];
196 out[10][1][0] = -12*in[1];
197 out[10][1][1] = 18 - 12*in[0] - 36*in[1];
199 out[11][0][0] = 90 - 180*in[1] - 180*in[0];
200 out[11][0][1] = -180*in[0];
202 out[11][1][0] = 180*in[1];
203 out[11][1][1] = -90 + 180*in[0] + 180*in[1];
216 #endif // DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI2_SIMPLEX2D_LOCALBASIS_HH