1 #ifndef DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI1Q2DLOCALBASIS_HH
2 #define DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI1Q2DLOCALBASIS_HH
6 #include <dune/common/fmatrix.hh>
8 #include "../../common/localbasis.hh"
21 template<
class D,
class R>
32 sign0 = sign1 = sign2 = sign3 = 1.0;
42 sign0 = sign1 = sign2 = sign3 = 1.0;
74 std::vector<typename Traits::RangeType>& out)
const
78 out[0][0] = sign0*(in[0] - 1.0);
80 out[1][0] = 6.0*in[0]*in[1] - 3.0*in[0]-6*in[1] + 3.0;
81 out[1][1] = -3.0*in[1]*in[1] + 3.0*in[1];
82 out[2][0] = sign1*(in[0]);
84 out[3][0] = -6.0*in[0]*in[1] + 3.0*in[0];
85 out[3][1] = 3.0*in[1]*in[1] - 3.0*in[1];
87 out[4][1] = sign2*(in[1] - 1.0);
88 out[5][0] = 3.0*in[0]*in[0] - 3.0*in[0];
89 out[5][1] = -6.0*in[0]*in[1] + 6.0*in[0] + 3.0*in[1] - 3.0;
91 out[6][1] = sign3*(in[1]);
92 out[7][0] = -3.0*in[0]*in[0] + 3.0*in[0];
93 out[7][1] = 6.0*in[0]*in[1] - 3.0*in[1];
103 std::vector<typename Traits::JacobianType>& out)
const
107 out[0][0][0] = sign0;
112 out[1][0][0] = 6.0*in[1] - 3.0;
113 out[1][0][1] = 6.0*in[0] - 6.0;
115 out[1][1][1] = -6.0*in[1] + 3.0;
117 out[2][0][0] = sign1;
122 out[3][0][0] = -6.0*in[1] + 3.0;
123 out[3][0][1] = -6.0*in[0];
125 out[3][1][1] = 6.0*in[1] - 3.0;
130 out[4][1][1] = sign2;
132 out[5][0][0] = 6.0*in[0] - 3.0;
134 out[5][1][0] = -6.0*in[1] + 6.0;
135 out[5][1][1] = -6.0*in[0] + 3.0;
140 out[6][1][1] = sign3;
142 out[7][0][0] = -6.0*in[0] + 3.0;
144 out[7][1][0] = 6.0*in[1];
145 out[7][1][1] = 6.0*in[0] - 3.0;
155 R sign0, sign1, sign2, sign3;
158 #endif // DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI1Q2DLOCALBASIS_HH