3 #ifndef DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS1_CUBE2D_LOCALBASIS_HH
4 #define DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS1_CUBE2D_LOCALBASIS_HH
8 #include <dune/common/fmatrix.hh>
10 #include "../../common/localbasis.hh"
23 template<
class D,
class R>
34 sign0 = sign1 = sign2 = sign3 = 1.0;
44 sign0 = sign1 = sign2 = sign3 = 1.0;
76 std::vector<typename Traits::RangeType>& out)
const
80 out[0][0] = sign0*(-1.0 + 4.0*in[0]-3*in[0]*in[0]);
82 out[1][0] = 3.0 - 12.0*in[0] - 6.0*in[1] + 24.0*in[0]*in[1]+9*in[0]*in[0] - 18.0*in[0]*in[0]*in[1];
84 out[2][0] = sign1*(-2.0*in[0] + 3.0*in[0]*in[0]);
86 out[3][0] = -6.0*in[0] + 12.0*in[0]*in[1] + 9.0*in[0]*in[0] - 18.0*in[0]*in[0]*in[1];
89 out[4][1] = sign2*(-1.0 + 4.0*in[1] - 3.0*in[1]*in[1]);
91 out[5][1] = -3.0 + 6.0*in[0] + 12.0*in[1] - 24.0*in[0]*in[1] - 9.0*in[1]*in[1] + 18.0*in[0]*in[1]*in[1];
93 out[6][1] = sign3*(-2.0*in[1] + 3.0*in[1]*in[1]);
95 out[7][1] = 6.0*in[1] - 12.0*in[0]*in[1] - 9.0*in[1]*in[1] + 18.0*in[0]*in[1]*in[1];
96 out[8][0] = 24.0*in[0] - 36.0*in[0]*in[1] - 24.0*in[0]*in[0] + 36.0*in[0]*in[0]*in[1];
99 out[9][1] = 24.0*in[1] - 36.0*in[0]*in[1] - 24.0*in[1]*in[1] + 36.0*in[0]*in[1]*in[1];
100 out[10][0] = -36.0*in[0] + 72.0*in[0]*in[1] + 36.0*in[0]*in[0] - 72.0*in[0]*in[0]*in[1];
103 out[11][1] = -36.0*in[1] + 72.0*in[0]*in[1] + 36*in[1]*in[1] - 72.0*in[0]*in[1]*in[1];
113 std::vector<typename Traits::JacobianType>& out)
const
117 out[0][0][0] = sign0*(4.0 - 6.0*in[0]);
122 out[1][0][0] = -12.0 + 24.0*in[1] + 18.0*in[0] - 36.0*in[0]*in[1];
123 out[1][0][1] = -6 + 24.0*in[0] - 18.0*in[0]*in[0];
127 out[2][0][0] = sign1*(-2.0 + 6.0*in[0]);
132 out[3][0][0] = -6.0 + 12.0*in[1] + 18.0*in[0] - 36.0*in[0]*in[1];
133 out[3][0][1] = 12.0*in[0] - 18.0*in[0]*in[0];
140 out[4][1][1] = sign2*(4.0 - 6.0*in[1]);
144 out[5][1][0] = 6.0 - 24.0*in[1] + 18.0*in[1]*in[1];
145 out[5][1][1] = 12.0 - 24.0*in[0] - 18.0*in[1] + 36.0*in[0]*in[1];
150 out[6][1][1] = sign3*(-2.0 + 6.0*in[1]);
154 out[7][1][0] = -12.0*in[1] + 18.0*in[1]*in[1];
155 out[7][1][1] = 6.0 - 12.0*in[0] - 18.0*in[1] + 36.0*in[1]*in[0];
157 out[8][0][0] = 24.0 - 36.0*in[1] - 48.0*in[0] + 72.0*in[0]*in[1];
158 out[8][0][1] = -36.0*in[0] + 36.0*in[0]*in[0];
164 out[9][1][0] = -36.0*in[1] + 36.0*in[1]*in[1];
165 out[9][1][1] = 24.0 - 36.0*in[0] - 48.0*in[1] + 72.0*in[0]*in[1];
167 out[10][0][0] = -36.0 + 72.0*in[1] + 72.0*in[0] - 144.0*in[0]*in[1];
168 out[10][0][1] = 72.0*in[0] - 72.0*in[0]*in[0];
174 out[11][1][0] = 72.0*in[1] - 72.0*in[1]*in[1];
175 out[11][1][1] = -36.0 + 72.0*in[0] + 72.0*in[1] - 144.0*in[0]*in[1];
185 R sign0, sign1, sign2, sign3;
188 #endif // DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS1_CUBE2D_LOCALBASIS_HH
unsigned int size() const
number of shape functions
Definition: raviartthomas1cube2dlocalbasis.hh:64
First order Raviart-Thomas shape functions on the reference quadrilateral.
Definition: raviartthomas1cube2dlocalbasis.hh:24
Type traits for LocalBasisVirtualInterface.
Definition: localbasis.hh:37
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:14
unsigned int order() const
Polynomial order of the shape functions.
Definition: raviartthomas1cube2dlocalbasis.hh:179
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: raviartthomas1cube2dlocalbasis.hh:75
D DomainType
domain type
Definition: localbasis.hh:49
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: raviartthomas1cube2dlocalbasis.hh:112
RT1Cube2DLocalBasis()
Standard constructor.
Definition: raviartthomas1cube2dlocalbasis.hh:32
RT1Cube2DLocalBasis(unsigned int s)
Make set number s, where 0 <= s < 16.
Definition: raviartthomas1cube2dlocalbasis.hh:42
LocalBasisTraits< D, 2, Dune::FieldVector< D, 2 >, R, 2, Dune::FieldVector< R, 2 >, Dune::FieldMatrix< R, 2, 2 > > Traits
Definition: raviartthomas1cube2dlocalbasis.hh:29