3 #ifndef DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI1_SIMPLEX2D_LOCALBASIS_HH 4 #define DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI1_SIMPLEX2D_LOCALBASIS_HH 11 #include <dune/common/fmatrix.hh> 13 #include "../../common/localbasis.hh" 26 template<
class D,
class R>
37 for (
size_t i=0; i<3; i++)
48 for (
size_t i=0; i<3; i++)
49 sign_[i] = s[i] ? -1.0 : 1.0;
65 std::vector<typename Traits::RangeType>& out)
const 69 out[0][0] = sign_[0]*in[0];
70 out[0][1] = sign_[0]*(in[1] - 1.0);
71 out[1][0] = sign_[1]*(in[0] - 1.0);
72 out[1][1] = sign_[1]*in[1];
73 out[2][0] = sign_[2]*in[0];
74 out[2][1] = sign_[2]*in[1];
75 out[3][0] = 3.0*in[0];
76 out[3][1] = 3.0 - 6.0*in[0] - 3.0*in[1];
77 out[4][0] = -3.0 + 3.0*in[0] + 6.0*in[1];
78 out[4][1] = -3.0*in[1];
79 out[5][0] = -3.0*in[0];
80 out[5][1] = 3.0*in[1];
90 std::vector<typename Traits::JacobianType>& out)
const 94 out[0][0][0] = sign_[0];
97 out[0][1][1] = sign_[0];
99 out[1][0][0] = sign_[1];
102 out[1][1][1] = sign_[1];
104 out[2][0][0] = sign_[2];
107 out[2][1][1] = sign_[2];
128 std::vector<typename Traits::RangeType>& out)
const 130 auto totalOrder = std::accumulate(
order.begin(),
order.end(), 0);
131 if (totalOrder == 0) {
133 }
else if (totalOrder == 1) {
135 auto const direction = std::distance(
order.begin(), std::find(
order.begin(),
order.end(), 1));
139 out[0][0] = sign_[0];
142 out[1][0] = sign_[1];
145 out[2][0] = sign_[2];
159 out[0][1] = sign_[0];
162 out[1][1] = sign_[1];
165 out[2][1] = sign_[2];
177 DUNE_THROW(RangeError,
"Component out of range.");
180 DUNE_THROW(NotImplemented,
"Desired derivative order is not implemented");
191 std::array<R,3> sign_;
194 #endif // DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI1_SIMPLEX2D_LOCALBASIS_HH BDM1Simplex2DLocalBasis(std::bitset< 3 > s)
Make set number s, where 0 <= s < 8.
Definition: brezzidouglasmarini1simplex2dlocalbasis.hh:46
LocalBasisTraits< D, 2, Dune::FieldVector< D, 2 >, R, 2, Dune::FieldVector< R, 2 >, Dune::FieldMatrix< R, 2, 2 > > Traits
Definition: brezzidouglasmarini1simplex2dlocalbasis.hh:32
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: brezzidouglasmarini1simplex2dlocalbasis.hh:64
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:15
BDM1Simplex2DLocalBasis()
Standard constructor.
Definition: brezzidouglasmarini1simplex2dlocalbasis.hh:35
D DomainType
domain type
Definition: localbasis.hh:49
unsigned int size() const
number of shape functions
Definition: brezzidouglasmarini1simplex2dlocalbasis.hh:53
First order Brezzi-Douglas-Marini shape functions on the reference triangle.
Definition: brezzidouglasmarini1simplex2dlocalbasis.hh:27
unsigned int order() const
Polynomial order of the shape functions.
Definition: brezzidouglasmarini1simplex2dlocalbasis.hh:185
Type traits for LocalBasisVirtualInterface.
Definition: localbasis.hh:37
void partial(const std::array< unsigned int, 2 > &order, const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate partial derivatives of all shape functions.
Definition: brezzidouglasmarini1simplex2dlocalbasis.hh:126
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: brezzidouglasmarini1simplex2dlocalbasis.hh:89