3 #ifndef DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI1_SIMPLEX2D_LOCALBASIS_HH
4 #define DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI1_SIMPLEX2D_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<3; i++)
47 for (
size_t i=0; i<3; 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];
69 out[0][1] = sign_[0]*(in[1] - 1.0);
70 out[1][0] = sign_[1]*(in[0] - 1.0);
71 out[1][1] = sign_[1]*in[1];
72 out[2][0] = sign_[2]*in[0];
73 out[2][1] = sign_[2]*in[1];
74 out[3][0] = 3.0*in[0];
75 out[3][1] = 3.0 - 6.0*in[0] - 3.0*in[1];
76 out[4][0] = -3.0 + 3.0*in[0] + 6.0*in[1];
77 out[4][1] = -3.0*in[1];
78 out[5][0] = -3.0*in[0];
79 out[5][1] = 3.0*in[1];
89 std::vector<typename Traits::JacobianType>& out)
const
93 out[0][0][0] = sign_[0];
96 out[0][1][1] = sign_[0];
98 out[1][0][0] = sign_[1];
101 out[1][1][1] = sign_[1];
103 out[2][0][0] = sign_[2];
106 out[2][1][1] = sign_[2];
131 std::array<R,3> sign_;
134 #endif // DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI1_SIMPLEX2D_LOCALBASIS_HH
LocalBasisTraits< D, 2, Dune::FieldVector< D, 2 >, R, 2, Dune::FieldVector< R, 2 >, Dune::FieldMatrix< R, 2, 2 > > Traits
Definition: brezzidouglasmarini1simplex2dlocalbasis.hh:31
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: brezzidouglasmarini1simplex2dlocalbasis.hh:88
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:14
BDM1Simplex2DLocalBasis()
Standard constructor.
Definition: brezzidouglasmarini1simplex2dlocalbasis.hh:34
unsigned int size() const
number of shape functions
Definition: brezzidouglasmarini1simplex2dlocalbasis.hh:52
BDM1Simplex2DLocalBasis(std::bitset< 3 > s)
Make set number s, where 0 <= s < 8.
Definition: brezzidouglasmarini1simplex2dlocalbasis.hh:45
First order Brezzi-Douglas-Marini shape functions on the reference triangle.
Definition: brezzidouglasmarini1simplex2dlocalbasis.hh:26
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: brezzidouglasmarini1simplex2dlocalbasis.hh:63
Type traits for LocalBasisVirtualInterface.
Definition: localbasis.hh:37
unsigned int order() const
Polynomial order of the shape functions.
Definition: brezzidouglasmarini1simplex2dlocalbasis.hh:125
D DomainType
domain type
Definition: localbasis.hh:49