3 #ifndef DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI2_CUBE2D_LOCALBASIS_HH
4 #define DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI2_CUBE2D_LOCALBASIS_HH
10 #include <dune/common/fmatrix.hh>
12 #include "../../common/localbasis.hh"
25 template<
class D,
class R>
31 R,2,Dune::FieldVector<R,2>,
37 for (
size_t i=0; i<4; i++)
48 for (
size_t i=0; i<4; i++)
49 sign_[i] = s[i] ? -1.0 : 1.0;
65 std::vector<typename Traits::RangeType>& out)
const
69 out[0][0] = sign_[0]*(-2.25 + 5.25*in[0] + 7.5*in[1] - 7.5*in[0]*in[1] - 3.0*in[0]*in[0] - 7.5*in[1]*in[1] + 7.5*in[0]*in[1]*in[1]);
70 out[0][1] = sign_[0]*(-1.25*in[1] + 3.75*in[1]*in[1] - 2.5*in[1]*in[1]*in[1]);
71 out[1][0] = 3.0 - 3.0*in[0]-6.0*in[1] + 6.0*in[0]*in[1];
73 out[2][0] = sign_[0]*(-3.75 + 3.75*in[0] + 22.5*in[1] - 22.5*in[0]*in[1] - 22.5*in[1]*in[1] + 22.5*in[0]*in[1]*in[1]);
74 out[2][1] = sign_[0]*(-3.75*in[1] + 11.25*in[1]*in[1] - 7.5*in[1]*in[1]*in[1]);
75 out[3][0] = sign_[1]*(-0.75*in[0] - 7.5*in[0]*in[1] + 3.0*in[0]*in[0] + 7.5*in[0]*in[1]*in[1]);
76 out[3][1] = sign_[1]*(-1.25*in[1] + 3.75*in[1]*in[1] - 2.5*in[1]*in[1]*in[1]);
77 out[4][0] = 3.0*in[0] - 6.0*in[0]*in[1];
79 out[5][0] = sign_[1]*(+3.75*in[0] - 22.5*in[0]*in[1] + 22.5*in[0]*in[1]*in[1]);
80 out[5][1] = sign_[1]*(-3.75*in[1] + 11.25*in[1]*in[1] - 7.5*in[1]*in[1]*in[1]);
81 out[6][0] = sign_[2]*(-1.25*in[0] + 3.75*in[0]*in[0] - 2.5*in[0]*in[0]*in[0]);
82 out[6][1] = sign_[2]*(-2.25 + 7.5*in[0] + 5.25*in[1] - 7.5*in[0]*in[1] - 7.5*in[0]*in[0] - 3.0*in[1]*in[1] + 7.5*in[0]*in[0]*in[1]);
84 out[7][1] = -3.0 + 6.0*in[0] + 3.0*in[1] - 6.0*in[0]*in[1];
85 out[8][0] = sign_[2]*(-3.75*in[0] + 11.25*in[0]*in[0] - 7.5*in[0]*in[0]*in[0]);
86 out[8][1] = sign_[2]*(-3.75 + 22.5*in[0] + 3.75*in[1] - 22.5*in[0]*in[1] - 22.5*in[0]*in[0] + 22.5*in[0]*in[0]*in[1]);
87 out[9][0] = sign_[3]*(-1.25*in[0] + 3.75*in[0]*in[0] - 2.5*in[0]*in[0]*in[0]);
88 out[9][1] = sign_[3]*(-0.75*in[1] - 7.5*in[0]*in[1] + 3.0*in[1]*in[1] + 7.5*in[0]*in[0]*in[1]);
90 out[10][1] = -3.0*in[1] + 6.0*in[0]*in[1];
91 out[11][0] = sign_[3]*(-3.75*in[0] + 11.25*in[0]*in[0] - 7.5*in[0]*in[0]*in[0]);
92 out[11][1] = sign_[3]*(3.75*in[1] - 22.5*in[0]*in[1] + 22.5*in[0]*in[0]*in[1]);
93 out[12][0] = 6.0*in[0] - 6.0*in[0]*in[0];
96 out[13][1] = 6.0*in[1] - 6.0*in[1]*in[1];
106 std::vector<typename Traits::JacobianType>& out)
const
110 out[0][0][0] = sign_[0]*(5.25 - 7.5*in[1] - 6.0*in[0] + 7.5*in[1]*in[1]);
111 out[0][0][1] = sign_[0]*(7.5 - 7.5*in[0] - 15.0*in[1] + 15.0*in[0]*in[1]);
113 out[0][1][1] = sign_[0]*(-1.25 + 7.5*in[1] - 7.5*in[1]*in[1]);
115 out[1][0][0] = -3.0 + 6.0*in[1];
116 out[1][0][1] = -6.0 + 6.0*in[0];
120 out[2][0][0] = sign_[0]*(3.75 - 22.5*in[1] + 22.5*in[1]*in[1]);
121 out[2][0][1] = sign_[0]*(22.5 - 22.5*in[0] - 45.0*in[1] + 45.0*in[0]*in[1]);
123 out[2][1][1] = sign_[0]*(-3.75 + 22.5*in[1] - 22.5*in[1]*in[1]);
125 out[3][0][0] = sign_[1]*(-0.75 - 7.5*in[1] + 6.0*in[0] + 7.5*in[1]*in[1]);
126 out[3][0][1] = sign_[1]*(-7.5*in[0] + 15.0*in[0]*in[1]);
128 out[3][1][1] = sign_[1]*(-1.25 + 7.5*in[1] - 7.5*in[1]*in[1]);
130 out[4][0][0] = 3.0 - 6.0*in[1];
131 out[4][0][1] = -6.0*in[0];
135 out[5][0][0] = sign_[1]*(3.75 - 22.5*in[1] + 22.5*in[1]*in[1]);
136 out[5][0][1] = sign_[1]*(-22.5*in[0] + 45.0*in[0]*in[1]);
138 out[5][1][1] = sign_[1]*(-3.75 + 22.5*in[1] - 22.5*in[1]*in[1]);
140 out[6][0][0] = sign_[2]*(-1.25 + 7.5*in[0] - 7.5*in[0]*in[0]);
142 out[6][1][0] = sign_[2]*(7.5 - 7.5*in[1] - 15.0*in[0] + 15.0*in[0]*in[1]);
143 out[6][1][1] = sign_[2]*(5.25 - 7.5*in[0]- 6.0*in[1] + 7.5*in[0]*in[0]);
147 out[7][1][0] = 6.0 - 6.0*in[1];
148 out[7][1][1] = 3.0 - 6.0*in[0];
150 out[8][0][0] = sign_[2]*(-3.75 + 22.5*in[0] - 22.5*in[0]*in[0]);
152 out[8][1][0] = sign_[2]*(22.5 - 22.5*in[1] - 45.0*in[0] + 45.0*in[0]*in[1]);
153 out[8][1][1] = sign_[2]*(3.75 - 22.5*in[0] + 22.5*in[0]*in[0]);
155 out[9][0][0] = sign_[3]*(-1.25 + 7.5*in[0] - 7.5*in[0]*in[0]);
157 out[9][1][0] = sign_[3]*(-7.5*in[1] + 15.0*in[0]*in[1]);
158 out[9][1][1] = sign_[3]*(-0.75 - 7.5*in[0] + 6.0*in[1] + 7.5*in[0]*in[0]);
162 out[10][1][0] = 6.0*in[1];
163 out[10][1][1] = -3.0 + 6.0*in[0];
165 out[11][0][0] = sign_[3]*(-3.75 + 22.5*in[0] - 22.5*in[0]*in[0]);
167 out[11][1][0] = sign_[3]*(-22.5*in[1] + 45*in[0]*in[1]);
168 out[11][1][1] = sign_[3]*(3.75 - 22.5*in[0] + 22.5*in[0]*in[0]);
170 out[12][0][0] = 6.0 - 12.0*in[0];
178 out[13][1][1] = 6.0 - 12.0*in[1];
188 std::array<R,4> sign_;
191 #endif // DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI2_CUBE2D_LOCALBASIS_HH
First order Brezzi-Douglas-Marini shape functions on quadrilaterals.
Definition: brezzidouglasmarini2cube2dlocalbasis.hh:26
unsigned int order() const
Polynomial order of the shape functions.
Definition: brezzidouglasmarini2cube2dlocalbasis.hh:182
LocalBasisTraits< D, 2, Dune::FieldVector< D, 2 >, R, 2, Dune::FieldVector< R, 2 >, Dune::FieldMatrix< R, 2, 2 > > Traits
Definition: brezzidouglasmarini2cube2dlocalbasis.hh:32
unsigned int size() const
number of shape functions
Definition: brezzidouglasmarini2cube2dlocalbasis.hh:53
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: brezzidouglasmarini2cube2dlocalbasis.hh:64
Type traits for LocalBasisVirtualInterface.
Definition: localbasis.hh:37
BDM2Cube2DLocalBasis()
Standard constructor.
Definition: brezzidouglasmarini2cube2dlocalbasis.hh:35
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: brezzidouglasmarini2cube2dlocalbasis.hh:105
D DomainType
domain type
Definition: localbasis.hh:49
BDM2Cube2DLocalBasis(std::bitset< 4 > s)
Make set number s, where 0 <= s < 16.
Definition: brezzidouglasmarini2cube2dlocalbasis.hh:46