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 out[0] = {sign_[0]*2.0* in[0], sign_[0]*2.0* in[1], sign_[0]*2.0*(in[2]-D(1))};
51 out[1] = {sign_[1]*2.0* in[0], sign_[1]*2.0*(in[1]-D(1)), sign_[1]*2.0* in[2] };
52 out[2] = {sign_[2]*2.0*(in[0]-D(1)), sign_[2]*2.0* in[1], sign_[2]*2.0* in[2] };
53 out[3] = {sign_[3]*2.0* in[0], sign_[3]*2.0* in[1], sign_[3]*2.0* in[2] };
59 std::vector<typename Traits::JacobianType>& out)
const
62 for (
int i=0; i<4; i++)
64 out[i][0] = {2.0*sign_[i], 0, 0};
65 out[i][1] = { 0,2.0*sign_[i], 0};
66 out[i][2] = { 0, 0,2.0*sign_[i]};
73 std::vector<typename Traits::RangeType>& out)
const
75 auto totalOrder = std::accumulate(
order.begin(),
order.end(), 0);
76 if (totalOrder == 0) {
78 }
else if (totalOrder == 1) {
79 auto const direction = std::distance(
order.begin(), std::find(
order.begin(),
order.end(), 1));
82 for (std::size_t i=0; i<
size(); i++)
84 out[i][direction] = sign_[i]*2.0;
85 out[i][(direction+1)%3] = 0;
86 out[i][(direction+2)%3] = 0;
90 for (std::size_t i = 0; i <
size(); ++i)
91 for (std::size_t j = 0; j < 3; ++j)
106 std::array<R,4> sign_;
Definition bdfmcube.hh:18
Type traits for LocalBasisVirtualInterface.
Definition common/localbasis.hh:35
Dune::LocalBasisTraits< D, 3, Dune::FieldVector< D, 3 >, R, 3, Dune::FieldVector< R, 3 >, Dune::FieldMatrix< R, 3, 3 > >::DomainType D DomainType
Definition common/localbasis.hh:43
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:71
unsigned int order() const
Polynomial order of the shape functions.
Definition raviartthomas03dlocalbasis.hh:98
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:58
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