3 #ifndef DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI1_CUBE3D_LOCALBASIS_HH
4 #define DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI1_CUBE3D_LOCALBASIS_HH
10 #include <dune/common/fmatrix.hh>
12 #include "../../common/localbasis.hh"
26 template<
class D,
class R>
32 R,3,Dune::FieldVector<R,3>,
38 for (
size_t i=0; i<6; i++)
49 for (
size_t i=0; i<6; i++)
50 sign_[i] = s[i] ? -1.0 : 1.0;
66 std::vector<typename Traits::RangeType>& out)
const
70 out[0][0] = sign_[0] * (in[0] - 1.0);
73 out[1][0] = sign_[1] * in[0];
77 out[2][1] = sign_[2] * (in[1] - 1.0);
80 out[3][1] = sign_[3] * in[1];
84 out[4][2] = sign_[4] * (in[2] - 1.0);
87 out[5][2] = sign_[5] * in[2];
88 out[6][0] = 6.0 * in[0] * in[1] - 3 * in[0]-6 * in[1] + 3.0;
89 out[6][1] = -3.0 * in[1] * in[1] + 3 * in[1];
91 out[7][0] = -6.0 * in[0] * in[1] + 3 * in[0];
92 out[7][1] = 3.0 * in[1] * in[1] - 3 * in[1];
94 out[8][0] = 3.0 * in[0] * in[0] - 3 * in[0];
95 out[8][1] = -6.0 * in[0] * in[1] + 3 * in[1]+6 * in[0]-3.0;
97 out[9][0] = -3.0 * in[0] * in[0] + 3 * in[0];
98 out[9][1] = 6.0 * in[0] * in[1] - 3 * in[1];
100 out[10][0] = -3.0 * in[0] * in[0] + 3 * in[0];
102 out[10][2] = 6.0 * in[0] * in[2]-6 * in[0]-3 * in[2] + 3.0;
103 out[11][0] = 3.0 * in[0] * in[0]-3 * in[0];
105 out[11][2] = -6.0 * in[0] * in[2] + 3 * in[2];
106 out[12][0] = -6.0 * in[0] * in[2]+6 * in[2] + 3 * in[0]-3.0;
108 out[12][2] = 3.0 * in[2] * in[2]-3 * in[2];
109 out[13][0] = -3 * in[0]+6 * in[0] * in[2];
111 out[13][2] = -3.0 * in[2] * in[2] + 3 * in[2];
113 out[14][1] = 6.0 * in[1] * in[2]-3 * in[1]-6 * in[2] + 3.0;
114 out[14][2] = -3 * in[2] * in[2] + 3 * in[2];
116 out[15][1] = -6.0 * in[1] * in[2] + 3 * in[1];
117 out[15][2] = 3.0 * in[2] * in[2]-3 * in[2];
119 out[16][1] = 3.0 * in[1] * in[1]-3 * in[1];
120 out[16][2] = -6.0 * in[1] * in[2] + 3 * in[2]+6 * in[1]-3.0;
122 out[17][1] = -3.0 * in[1] * in[1] + 3 * in[1];
123 out[17][2] = 6.0 * in[1] * in[2] - 3.0 * in[2];
133 std::vector<typename Traits::JacobianType>& out)
const
137 out[0][0][0] = sign_[0]; out[0][0][1] = 0; out[0][0][2] = 0;
138 out[0][1][0] = 0; out[0][1][1] = 0; out[0][1][2] = 0;
139 out[0][2][0] = 0; out[0][2][1] = 0; out[0][2][2] = 0;
141 out[1][0][0] = sign_[1]; out[1][0][1] = 0; out[1][0][2] = 0;
142 out[1][1][0] = 0; out[1][1][1] = 0; out[1][1][2] = 0;
143 out[1][2][0] = 0; out[1][2][1] = 0; out[1][2][2] = 0;
145 out[2][0][0] = 0; out[2][0][1] = 0; out[2][0][2] = 0;
146 out[2][1][0] = 0; out[2][1][1] = sign_[2]; out[2][1][2] = 0;
147 out[2][2][0] = 0; out[2][2][1] = 0; out[2][2][2] = 0;
149 out[3][0][0] = 0; out[3][0][1] = 0; out[3][0][2] = 0;
150 out[3][1][0] = 0; out[3][1][1] = sign_[3]; out[3][1][2] = 0;
151 out[3][2][0] = 0; out[3][2][1] = 0; out[3][2][2] = 0;
153 out[4][0][0] = 0; out[4][0][1] = 0; out[4][0][2] = 0;
154 out[4][1][0] = 0; out[4][1][1] = 0; out[4][1][2] = 0;
155 out[4][2][0] = 0; out[4][2][1] = 0; out[4][2][2] = sign_[4];
157 out[5][0][0] = 0; out[5][0][1] = 0; out[5][0][2] = 0;
158 out[5][1][0] = 0; out[5][1][1] = 0; out[5][1][2] = 0;
159 out[5][2][0] = 0; out[5][2][1] = 0; out[5][2][2] = sign_[5];
163 out[6][0][0] = 6*in[1]-3; out[6][0][1] = 6*in[0]-6; out[6][0][2] = 0.0;
164 out[6][1][0] = 0.0; out[6][1][1] =-6*in[1]+3; out[6][1][2] = 0.0;
165 out[6][2][0] = 0.0; out[6][2][1] = 0.0; out[6][2][2] = 0.0;
167 out[7][0][0] = -6*in[1]+3; out[7][0][1] =-6*in[0]; out[7][0][2] = 0.0;
168 out[7][1][0] = 0.0; out[7][1][1] = 6*in[1]-3; out[7][1][2] = 0.0;
169 out[7][2][0] = 0.0; out[7][2][1] = 0.0; out[7][2][2] = 0.0;
172 out[8][0][0] = 6*in[0]-3; out[8][0][1] = 0.0; out[8][0][2] = 0.0;
173 out[8][1][0] =-6*in[1]+6; out[8][1][1] =-6*in[0]+3; out[8][1][2] = 0.0;
174 out[8][2][0] = 0.0; out[8][2][1] = 0.0; out[8][2][2] = 0.0;
176 out[9][0][0] = -6*in[0]+3; out[9][0][1] = 0.0; out[9][0][2] = 0.0;
177 out[9][1][0] = 6*in[1]; out[9][1][1] = 6*in[0]-3; out[9][1][2] = 0.0;
178 out[9][2][0] = 0.0; out[9][2][1] = 0.0; out[9][2][2] = 0.0;
180 out[10][0][0] = -6*in[0]+3; out[10][0][1] = 0.0; out[10][0][2] = 0.0;
181 out[10][1][0] = 0.0; out[10][1][1] = 0.0; out[10][1][2] = 0.0;
182 out[10][2][0] =6*in[2]-6; out[10][2][1] = 0.0; out[10][2][2] =6*in[0]-3;
184 out[11][0][0] =6*in[0]-3; out[11][0][1] = 0.0; out[11][0][2] = 0.0;
185 out[11][1][0] = 0.0; out[11][1][1] = 0.0; out[11][1][2] = 0.0;
186 out[11][2][0] = -6*in[2]; out[11][2][1] = 0.0; out[11][2][2] = -6*in[0]+3;
189 out[12][0][0] = -6*in[2]+3; out[12][0][1] =0.0; out[12][0][2] = -6*in[0]+6;
190 out[12][1][0] = 0.0; out[12][1][1] = 0.0; out[12][1][2] = 0.0;
191 out[12][2][0] = 0.0; out[12][2][1] = 0.0; out[12][2][2] =6*in[2]-3;
194 out[13][0][0] =6*in[2]-3; out[13][0][1] = 0.0; out[13][0][2] = 6*in[0];
195 out[13][1][0] = 0.0; out[13][1][1] = 0.0; out[13][1][2] = 0.0;
196 out[13][2][0] = 0.0; out[13][2][1] = 0.0; out[13][2][2] =-6*in[2]+3;
198 out[14][0][0] = 0.0; out[14][0][1] = 0.0; out[14][0][2] = 0.0;
199 out[14][1][0] = 0.0; out[14][1][1] = 6*in[2]-3; out[14][1][2] = 6*in[1]-6;
200 out[14][2][0] = 0.0; out[14][2][1] = 0.0; out[14][2][2] =-6*in[2]+3;
202 out[15][0][0] = 0.0; out[15][0][1] = 0.0; out[15][0][2] = 0.0;
203 out[15][1][0] = 0.0; out[15][1][1] =-6*in[2]+3; out[15][1][2] =-6*in[1];
204 out[15][2][0] = 0.0; out[15][2][1] = 0.0; out[15][2][2] = 6*in[2]-3;
206 out[16][0][0] = 0.0; out[16][0][1] = 0.0; out[16][0][2] = 0.0;
207 out[16][1][0] = 0.0; out[16][1][1] = 6*in[1]-3; out[16][1][2] = 0.0;
208 out[16][2][0] = 0.0; out[16][2][1] =-6*in[2]+6; out[16][2][2] =-6*in[1]+3;
210 out[17][0][0] = 0.0; out[17][0][1] = 0.0; out[17][0][2] = 0.0;
211 out[17][1][0] = 0.0; out[17][1][1] =-6*in[1]+3; out[17][1][2] = 0.0;
212 out[17][2][0] = 0.0; out[17][2][1] = 6*in[2]; out[17][2][2] = 6*in[1]-3;
222 std::array<R,6> sign_;
225 #endif // DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI1_CUBE3D_LOCALBASIS_HH
unsigned int size() const
number of shape functions
Definition: brezzidouglasmarini1cube3dlocalbasis.hh:54
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: brezzidouglasmarini1cube3dlocalbasis.hh:65
BDM1Cube3DLocalBasis(std::bitset< 6 > s)
Make set number s, where 0 <= s < 64.
Definition: brezzidouglasmarini1cube3dlocalbasis.hh:47
BDM1Cube3DLocalBasis()
Standard constructor.
Definition: brezzidouglasmarini1cube3dlocalbasis.hh:36
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: brezzidouglasmarini1cube3dlocalbasis.hh:132
Type traits for LocalBasisVirtualInterface.
Definition: localbasis.hh:37
LocalBasisTraits< D, 3, Dune::FieldVector< D, 3 >, R, 3, Dune::FieldVector< R, 3 >, Dune::FieldMatrix< R, 3, 3 > > Traits
Definition: brezzidouglasmarini1cube3dlocalbasis.hh:33
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:14
unsigned int order() const
Polynomial order of the shape functions.
Definition: brezzidouglasmarini1cube3dlocalbasis.hh:216
First order Brezzi-Douglas-Marini shape functions on the reference hexahedron.
Definition: brezzidouglasmarini1cube3dlocalbasis.hh:27
D DomainType
domain type
Definition: localbasis.hh:49