3 #ifndef DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI1_CUBE2D_LOCALBASIS_HH 4 #define DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI1_CUBE2D_LOCALBASIS_HH 10 #include <dune/common/fmatrix.hh> 12 #include "../../common/localbasis.hh" 25 template<
class D,
class R>
36 for (
size_t i=0; i<4; i++)
47 for (
size_t i=0; i<4; i++)
48 sign_[i] = s[i] ? -1.0 : 1.0;
64 std::vector<typename Traits::RangeType>& out)
const 68 out[0][0] = sign_[0]*(in[0] - 1.0);
70 out[1][0] = 6.0*in[0]*in[1] - 3.0*in[0]-6*in[1] + 3.0;
71 out[1][1] = -3.0*in[1]*in[1] + 3.0*in[1];
72 out[2][0] = sign_[1]*(in[0]);
74 out[3][0] = -6.0*in[0]*in[1] + 3.0*in[0];
75 out[3][1] = 3.0*in[1]*in[1] - 3.0*in[1];
77 out[4][1] = sign_[2]*(in[1] - 1.0);
78 out[5][0] = 3.0*in[0]*in[0] - 3.0*in[0];
79 out[5][1] = -6.0*in[0]*in[1] + 6.0*in[0] + 3.0*in[1] - 3.0;
81 out[6][1] = sign_[3]*(in[1]);
82 out[7][0] = -3.0*in[0]*in[0] + 3.0*in[0];
83 out[7][1] = 6.0*in[0]*in[1] - 3.0*in[1];
93 std::vector<typename Traits::JacobianType>& out)
const 97 out[0][0][0] = sign_[0];
102 out[1][0][0] = 6.0*in[1] - 3.0;
103 out[1][0][1] = 6.0*in[0] - 6.0;
105 out[1][1][1] = -6.0*in[1] + 3.0;
107 out[2][0][0] = sign_[1];
112 out[3][0][0] = -6.0*in[1] + 3.0;
113 out[3][0][1] = -6.0*in[0];
115 out[3][1][1] = 6.0*in[1] - 3.0;
120 out[4][1][1] = sign_[2];
122 out[5][0][0] = 6.0*in[0] - 3.0;
124 out[5][1][0] = -6.0*in[1] + 6.0;
125 out[5][1][1] = -6.0*in[0] + 3.0;
130 out[6][1][1] = sign_[3];
132 out[7][0][0] = -6.0*in[0] + 3.0;
134 out[7][1][0] = 6.0*in[1];
135 out[7][1][1] = 6.0*in[0] - 3.0;
145 std::array<R,4> sign_;
148 #endif // DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI1_CUBE2D_LOCALBASIS_HH unsigned int size() const
number of shape functions
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:52
Type traits for LocalBasisVirtualInterface.
Definition: localbasis.hh:37
BDM1Cube2DLocalBasis(std::bitset< 4 > s)
Make set number s, where 0 <= s < 16.
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:45
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:92
First order Brezzi-Douglas-Marini shape functions on the reference quadrilateral. ...
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:26
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:14
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:63
unsigned int order() const
Polynomial order of the shape functions.
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:139
BDM1Cube2DLocalBasis()
Standard constructor.
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:34
D DomainType
domain type
Definition: localbasis.hh:49
LocalBasisTraits< D, 2, Dune::FieldVector< D, 2 >, R, 2, Dune::FieldVector< R, 2 >, Dune::FieldMatrix< R, 2, 2 > > Traits
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:31