5#ifndef DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI1_SIMPLEX2D_LOCALBASIS_HH
6#define DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI1_SIMPLEX2D_LOCALBASIS_HH
13#include <dune/common/fmatrix.hh>
15#include "../../common/localbasis.hh"
28 template<
class D,
class R>
39 for (
size_t i=0; i<3; i++)
50 for (
size_t i=0; i<3; i++)
51 sign_[i] = s[i] ? -1.0 : 1.0;
67 std::vector<typename Traits::RangeType>& out)
const
71 out[0][0] = sign_[0]*in[0];
72 out[0][1] = sign_[0]*(in[1] - 1.0);
73 out[1][0] = sign_[1]*(in[0] - 1.0);
74 out[1][1] = sign_[1]*in[1];
75 out[2][0] = sign_[2]*in[0];
76 out[2][1] = sign_[2]*in[1];
77 out[3][0] = 3.0*in[0];
78 out[3][1] = 3.0 - 6.0*in[0] - 3.0*in[1];
79 out[4][0] = -3.0 + 3.0*in[0] + 6.0*in[1];
80 out[4][1] = -3.0*in[1];
81 out[5][0] = -3.0*in[0];
82 out[5][1] = 3.0*in[1];
92 std::vector<typename Traits::JacobianType>& out)
const
96 out[0][0][0] = sign_[0];
99 out[0][1][1] = sign_[0];
101 out[1][0][0] = sign_[1];
104 out[1][1][1] = sign_[1];
106 out[2][0][0] = sign_[2];
109 out[2][1][1] = sign_[2];
130 std::vector<typename Traits::RangeType>& out)
const
132 auto totalOrder = std::accumulate(
order.begin(),
order.end(), 0);
133 if (totalOrder == 0) {
135 }
else if (totalOrder == 1) {
137 auto const direction = std::distance(
order.begin(), std::find(
order.begin(),
order.end(), 1));
141 out[0][0] = sign_[0];
144 out[1][0] = sign_[1];
147 out[2][0] = sign_[2];
161 out[0][1] = sign_[0];
164 out[1][1] = sign_[1];
167 out[2][1] = sign_[2];
179 DUNE_THROW(RangeError,
"Component out of range.");
182 DUNE_THROW(NotImplemented,
"Desired derivative order is not implemented");
193 std::array<R,3> sign_;
Definition bdfmcube.hh:18
First order Brezzi-Douglas-Marini shape functions on the reference triangle.
Definition brezzidouglasmarini1simplex2dlocalbasis.hh:30
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition brezzidouglasmarini1simplex2dlocalbasis.hh:91
unsigned int size() const
number of shape functions
Definition brezzidouglasmarini1simplex2dlocalbasis.hh:55
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition brezzidouglasmarini1simplex2dlocalbasis.hh:66
LocalBasisTraits< D, 2, Dune::FieldVector< D, 2 >, R, 2, Dune::FieldVector< R, 2 >, Dune::FieldMatrix< R, 2, 2 > > Traits
Definition brezzidouglasmarini1simplex2dlocalbasis.hh:34
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:128
BDM1Simplex2DLocalBasis()
Standard constructor.
Definition brezzidouglasmarini1simplex2dlocalbasis.hh:37
unsigned int order() const
Polynomial order of the shape functions.
Definition brezzidouglasmarini1simplex2dlocalbasis.hh:187
BDM1Simplex2DLocalBasis(std::bitset< 3 > s)
Make set number s, where 0 <= s < 8.
Definition brezzidouglasmarini1simplex2dlocalbasis.hh:48
Type traits for LocalBasisVirtualInterface.
Definition common/localbasis.hh:35
D DomainType
domain type
Definition common/localbasis.hh:43