3 #ifndef DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS12DLOCALBASIS_HH
4 #define DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS12DLOCALBASIS_HH
8 #include <dune/common/fmatrix.hh>
10 #include "../../common/localbasis.hh"
23 template<
class D,
class R>
33 sign0 = sign1 = sign2 = 1.0;
43 sign0 = sign1 = sign2 = 1.0;
71 std::vector<typename Traits::RangeType>& out)
const
74 out[0][0] = sign0*(in[0] - 4.0*in[0]*in[1]);
75 out[0][1] = sign0*(-1.0 + 5.0*in[1] - 4.0*in[1]*in[1]);
76 out[1][0] = sign1*(-1.0 + 5.0*in[0] - 4.0*in[0]*in[0]);
77 out[1][1] = sign1*(in[1] - 4.0*in[0]*in[1]);
78 out[2][0] = sign2*(-3.0*in[0] + 4.0*in[0]*in[0] + 4.0*in[1]*in[0]);
79 out[2][1] = sign2*(-3.0*in[1] + 4.0*in[0]*in[1] + 4.0*in[1]*in[1]);
80 out[3][0] = -5.0*in[0] + 8.0*in[0]*in[0] + 4.0*in[1]*in[0];
81 out[3][1] = 3.0 - 6.0*in[0] - 7.0*in[1] + 8.0*in[0]*in[1] + 4.0*in[1]*in[1];
82 out[4][0] = -3.0 + 7.0*in[0] + 6.0*in[1] - 4.0*in[0]*in[0] - 8.0*in[1]*in[0];
83 out[4][1] = 5.0*in[1] - 4.0*in[0]*in[1] - 8.0*in[1]*in[1];
84 out[5][0] = in[0] - 4.0*in[0]*in[0] + 4.0*in[1]*in[0];
85 out[5][1] = -1.0*in[1] - 4.0*in[0]*in[1] + 4.0*in[1]*in[1];
86 out[6][0] = 16.0*in[0] - 16.0*in[0]*in[0] - 8.0*in[1]*in[0];
87 out[6][1] = 8.0*in[1] - 16.0*in[0]*in[1] - 8.0*in[1]*in[1];
88 out[7][0] = 8.0*in[0] - 8.0*in[0]*in[0] - 16.0*in[1]*in[0];
89 out[7][1] = 16.0*in[1] - 8.0*in[0]*in[1] - 16.0*in[1]*in[1];
99 std::vector<typename Traits::JacobianType>& out)
const
103 out[0][0][0] = sign0*(1.0 - 4.0*in[1]);
104 out[0][0][1] = sign0*(-4.0*in[0]);
106 out[0][1][1] = sign0*(5.0 - 8.0*in[1]);
108 out[1][0][0] = sign1*(5.0 - 8.0*in[0]);
110 out[1][1][0] = sign1*(-4.0*in[1]);
111 out[1][1][1] = sign1*(1.0 - 4.0*in[0]);
113 out[2][0][0] = sign2*(-3.0 + 8.0*in[0] + 4.0*in[1]);
114 out[2][0][1] = sign2*(4.0*in[0]);
115 out[2][1][0] = sign2*(4.0*in[1]);
116 out[2][1][1] = sign2*(-3.0 + 4.0*in[0] + 8.0*in[1]);
118 out[3][0][0] = -5.0 + 16.0*in[0] + 4.0*in[1];
119 out[3][0][1] = 4.0*in[0];
120 out[3][1][0] = -6.0 + 8.0*in[1];
121 out[3][1][1] = -7.0 + 8.0*in[0] + 8.0*in[1];
123 out[4][0][0] = 7.0 - 8.0*in[0] - 8.0*in[1];
124 out[4][0][1] = 6.0 - 8.0*in[0];
125 out[4][1][0] = -4.0*in[1];
126 out[4][1][1] = 5.0 - 4.0*in[0] - 16.0*in[1];
128 out[5][0][0] = 1.0 - 8.0*in[0] + 4*in[1];
129 out[5][0][1] = 4.0*in[0];
130 out[5][1][0] = -4.0*in[1];
131 out[5][1][1] = -1.0 - 4.0*in[0] + 8.0*in[1];
133 out[6][0][0] = 16.0 - 32.0*in[0] - 8.0*in[1];
134 out[6][0][1] = -8.0*in[0];
135 out[6][1][0] = -16.0*in[1];
136 out[6][1][1] = 8.0 - 16.0*in[0] - 16.0*in[1];
138 out[7][0][0] = 8.0 - 16.0*in[0] - 16.0*in[1];
139 out[7][0][1] = -16.0*in[0];
140 out[7][1][0] = -8.0*in[1];
141 out[7][1][1] = 16.0 - 8.0*in[0] - 32.0*in[1];
151 R sign0, sign1, sign2;
154 #endif // DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS12DLOCALBASIS_HH
LocalBasisTraits< D, 2, Dune::FieldVector< D, 2 >, R, 2, Dune::FieldVector< R, 2 >, Dune::FieldMatrix< R, 2, 2 > > Traits
Definition: raviartthomas12dlocalbasis.hh:29
First order Raviart-Thomas shape functions on the reference triangle.
Definition: raviartthomas12dlocalbasis.hh:24
unsigned int order() const
Polynomial order of the shape functions.
Definition: raviartthomas12dlocalbasis.hh:145
RT12DLocalBasis(unsigned int s)
Make set number s, where 0 <= s < 8.
Definition: raviartthomas12dlocalbasis.hh:41
Type traits for LocalBasisVirtualInterface.
Definition: localbasis.hh:37
RT12DLocalBasis()
Standard constructor.
Definition: raviartthomas12dlocalbasis.hh:31
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: raviartthomas12dlocalbasis.hh:98
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: raviartthomas12dlocalbasis.hh:70
D DomainType
domain type
Definition: localbasis.hh:49
unsigned int size() const
number of shape functions
Definition: raviartthomas12dlocalbasis.hh:59