3 #ifndef DUNE_RT0TRIANGLELOCALBASIS_HH 4 #define DUNE_RT0TRIANGLELOCALBASIS_HH 8 #include <dune/common/fmatrix.hh> 22 template<
class D,
class R>
32 std::fill(sign_.begin(), sign_.end(), 1.0);
38 for (
int i=0; i<3; i++)
39 sign_[i] = s[i] ? -1.0 : 1.0;
50 std::vector<typename Traits::RangeType>& out)
const 53 out[0] = {sign_[0]*in[0], sign_[0]*(in[1]-D(1))};
54 out[1] = {sign_[1]*(in[0]-D(1)), sign_[1]*in[1]};
55 out[2] = {sign_[2]*in[0], sign_[2]*in[1]};
61 std::vector<typename Traits::JacobianType>& out)
const 64 for (
int i=0; i<3; i++)
66 out[i][0] = {sign_[i], 0};
67 out[i][1] = { 0, sign_[i]};
74 std::vector<typename Traits::RangeType>& out)
const 76 auto totalOrder = std::accumulate(
order.begin(),
order.end(), 0);
77 if (totalOrder == 0) {
79 }
else if (totalOrder == 1) {
80 auto const direction = std::distance(
order.begin(), std::find(
order.begin(),
order.end(), 1));
83 for (
int i=0; i<3; i++)
85 out[i][direction] = sign_[i];
86 out[i][1-direction] = 0;
90 for (std::size_t i = 0; i <
size(); ++i)
91 for (std::size_t j = 0; j < 2; ++j)
106 std::array<R,3> sign_;
RT02DLocalBasis(std::bitset< 3 > s)
Make set numer s, where 0<=s<4.
Definition: raviartthomas02dlocalbasis.hh:36
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: raviartthomas02dlocalbasis.hh:49
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: raviartthomas02dlocalbasis.hh:72
Lowest order Raviart-Thomas shape functions on the reference triangle.
Definition: raviartthomas02dlocalbasis.hh:23
Type traits for LocalBasisVirtualInterface.
Definition: localbasis.hh:31
D DomainType
domain type
Definition: localbasis.hh:43
unsigned int size() const
number of shape functions
Definition: raviartthomas02dlocalbasis.hh:43
RT02DLocalBasis()
Standard constructor.
Definition: raviartthomas02dlocalbasis.hh:30
LocalBasisTraits< D, 2, Dune::FieldVector< D, 2 >, R, 2, Dune::FieldVector< R, 2 >, Dune::FieldMatrix< R, 2, 2 > > Traits
Definition: raviartthomas02dlocalbasis.hh:27
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:15
unsigned int order() const
Polynomial order of the shape functions.
Definition: raviartthomas02dlocalbasis.hh:98
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: raviartthomas02dlocalbasis.hh:60