5#ifndef DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS_RAVIARTTHOMAS03D_RAVIARTTHOMAS03DLOCALBASIS_HH
6#define DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS_RAVIARTTHOMAS03D_RAVIARTTHOMAS03DLOCALBASIS_HH
10#include <dune/common/fmatrix.hh>
25 template<
class D,
class R>
35 for (
int i=0; i<4; i++)
36 sign_[i] = s[i] ? -1.0 : 1.0;
47 std::vector<typename Traits::RangeType>& out)
const
50 auto c = std::sqrt(2.0);
51 out[0] = {sign_[0]*c* in[0], sign_[0]*c* in[1], sign_[0]*c*(in[2]-D(1))};
52 out[1] = {sign_[1]*c* in[0], sign_[1]*c*(in[1]-D(1)), sign_[1]*c* in[2] };
53 out[2] = {sign_[2]*c*(in[0]-D(1)), sign_[2]*c* in[1], sign_[2]*c* in[2] };
54 out[3] = {sign_[3]*c* in[0], sign_[3]*c* in[1], sign_[3]*c* in[2] };
60 std::vector<typename Traits::JacobianType>& out)
const
63 for (
int i=0; i<4; i++)
65 auto c = std::sqrt(2.0);
66 out[i][0] = {c*sign_[i], 0, 0};
67 out[i][1] = { 0,c*sign_[i], 0};
68 out[i][2] = { 0, 0,c*sign_[i]};
75 std::vector<typename Traits::RangeType>& out)
const
77 auto totalOrder = std::accumulate(
order.begin(),
order.end(), 0);
78 if (totalOrder == 0) {
80 }
else if (totalOrder == 1) {
81 auto const direction = std::distance(
order.begin(), std::find(
order.begin(),
order.end(), 1));
84 for (
int i=0; i<
size(); i++)
86 out[i][direction] = sign_[i]* std::sqrt(2.0) ;
87 out[i][(direction+1)%3] = 0;
88 out[i][(direction+2)%3] = 0;
92 for (std::size_t i = 0; i <
size(); ++i)
93 for (std::size_t j = 0; j < 3; ++j)
108 std::array<R,4> sign_;
Definition bdfmcube.hh:18
Type traits for LocalBasisVirtualInterface.
Definition common/localbasis.hh:35
D DomainType
domain type
Definition common/localbasis.hh:43
Definition raviartthomas03dlocalbasis.hh:27
RT03DLocalBasis(std::bitset< 4 > s=0)
Make set number s, where 0 <= s < 16.
Definition raviartthomas03dlocalbasis.hh:33
void partial(const std::array< unsigned int, 3 > &order, const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate partial derivatives of all shape functions.
Definition raviartthomas03dlocalbasis.hh:73
unsigned int order() const
Polynomial order of the shape functions.
Definition raviartthomas03dlocalbasis.hh:100
unsigned int size() const
number of shape functions
Definition raviartthomas03dlocalbasis.hh:40
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition raviartthomas03dlocalbasis.hh:59
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition raviartthomas03dlocalbasis.hh:46
LocalBasisTraits< D, 3, Dune::FieldVector< D, 3 >, R, 3, Dune::FieldVector< R, 3 >, Dune::FieldMatrix< R, 3, 3 > > Traits
Definition raviartthomas03dlocalbasis.hh:30