3 #ifndef DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI1_CUBE2D_LOCALBASIS_HH
4 #define DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI1_CUBE2D_LOCALBASIS_HH
9 #include <dune/common/fmatrix.hh>
11 #include "../../common/localbasis.hh"
24 template<
class D,
class R>
35 for (
size_t i=0; i<4; i++)
46 for (
size_t i=0; i<4; i++)
47 sign_[i] = s[i] ? -1.0 : 1.0;
63 std::vector<typename Traits::RangeType>& out)
const
67 out[0][0] = sign_[0]*(in[0] - 1.0);
69 out[1][0] = 6.0*in[0]*in[1] - 3.0*in[0]-6*in[1] + 3.0;
70 out[1][1] = -3.0*in[1]*in[1] + 3.0*in[1];
71 out[2][0] = sign_[1]*(in[0]);
73 out[3][0] = -6.0*in[0]*in[1] + 3.0*in[0];
74 out[3][1] = 3.0*in[1]*in[1] - 3.0*in[1];
76 out[4][1] = sign_[2]*(in[1] - 1.0);
77 out[5][0] = 3.0*in[0]*in[0] - 3.0*in[0];
78 out[5][1] = -6.0*in[0]*in[1] + 6.0*in[0] + 3.0*in[1] - 3.0;
80 out[6][1] = sign_[3]*(in[1]);
81 out[7][0] = -3.0*in[0]*in[0] + 3.0*in[0];
82 out[7][1] = 6.0*in[0]*in[1] - 3.0*in[1];
92 std::vector<typename Traits::JacobianType>& out)
const
96 out[0][0][0] = sign_[0];
101 out[1][0][0] = 6.0*in[1] - 3.0;
102 out[1][0][1] = 6.0*in[0] - 6.0;
104 out[1][1][1] = -6.0*in[1] + 3.0;
106 out[2][0][0] = sign_[1];
111 out[3][0][0] = -6.0*in[1] + 3.0;
112 out[3][0][1] = -6.0*in[0];
114 out[3][1][1] = 6.0*in[1] - 3.0;
119 out[4][1][1] = sign_[2];
121 out[5][0][0] = 6.0*in[0] - 3.0;
123 out[5][1][0] = -6.0*in[1] + 6.0;
124 out[5][1][1] = -6.0*in[0] + 3.0;
129 out[6][1][1] = sign_[3];
131 out[7][0][0] = -6.0*in[0] + 3.0;
133 out[7][1][0] = 6.0*in[1];
134 out[7][1][1] = 6.0*in[0] - 3.0;
147 #endif // DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI1_CUBE2D_LOCALBASIS_HH