3 #ifndef DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS1_CUBE2D_LOCALBASIS_HH 4 #define DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS1_CUBE2D_LOCALBASIS_HH 9 #include <dune/common/fmatrix.hh> 11 #include "../../common/localbasis.hh" 24 template<
class D,
class R>
35 sign0 = sign1 = sign2 = sign3 = 1.0;
45 sign0 = sign1 = sign2 = sign3 = 1.0;
77 std::vector<typename Traits::RangeType>& out)
const 81 out[0][0] = sign0*(-1.0 + 4.0*in[0]-3*in[0]*in[0]);
83 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];
85 out[2][0] = sign1*(-2.0*in[0] + 3.0*in[0]*in[0]);
87 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];
90 out[4][1] = sign2*(-1.0 + 4.0*in[1] - 3.0*in[1]*in[1]);
92 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];
94 out[6][1] = sign3*(-2.0*in[1] + 3.0*in[1]*in[1]);
96 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];
97 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];
100 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];
101 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];
104 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];
114 std::vector<typename Traits::JacobianType>& out)
const 118 out[0][0][0] = sign0*(4.0 - 6.0*in[0]);
123 out[1][0][0] = -12.0 + 24.0*in[1] + 18.0*in[0] - 36.0*in[0]*in[1];
124 out[1][0][1] = -6 + 24.0*in[0] - 18.0*in[0]*in[0];
128 out[2][0][0] = sign1*(-2.0 + 6.0*in[0]);
133 out[3][0][0] = -6.0 + 12.0*in[1] + 18.0*in[0] - 36.0*in[0]*in[1];
134 out[3][0][1] = 12.0*in[0] - 18.0*in[0]*in[0];
141 out[4][1][1] = sign2*(4.0 - 6.0*in[1]);
145 out[5][1][0] = 6.0 - 24.0*in[1] + 18.0*in[1]*in[1];
146 out[5][1][1] = 12.0 - 24.0*in[0] - 18.0*in[1] + 36.0*in[0]*in[1];
151 out[6][1][1] = sign3*(-2.0 + 6.0*in[1]);
155 out[7][1][0] = -12.0*in[1] + 18.0*in[1]*in[1];
156 out[7][1][1] = 6.0 - 12.0*in[0] - 18.0*in[1] + 36.0*in[1]*in[0];
158 out[8][0][0] = 24.0 - 36.0*in[1] - 48.0*in[0] + 72.0*in[0]*in[1];
159 out[8][0][1] = -36.0*in[0] + 36.0*in[0]*in[0];
165 out[9][1][0] = -36.0*in[1] + 36.0*in[1]*in[1];
166 out[9][1][1] = 24.0 - 36.0*in[0] - 48.0*in[1] + 72.0*in[0]*in[1];
168 out[10][0][0] = -36.0 + 72.0*in[1] + 72.0*in[0] - 144.0*in[0]*in[1];
169 out[10][0][1] = 72.0*in[0] - 72.0*in[0]*in[0];
175 out[11][1][0] = 72.0*in[1] - 72.0*in[1]*in[1];
176 out[11][1][1] = -36.0 + 72.0*in[0] + 72.0*in[1] - 144.0*in[0]*in[1];
182 std::vector<typename Traits::RangeType>& out)
const 184 auto totalOrder = std::accumulate(
order.begin(),
order.end(), 0);
185 if (totalOrder == 0) {
188 DUNE_THROW(NotImplemented,
"Desired derivative order is not implemented");
199 R sign0, sign1, sign2, sign3;
202 #endif // DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS1_CUBE2D_LOCALBASIS_HH LocalBasisTraits< D, 2, Dune::FieldVector< D, 2 >, R, 2, Dune::FieldVector< R, 2 >, Dune::FieldMatrix< R, 2, 2 > > Traits
Definition: raviartthomas1cube2dlocalbasis.hh:30
unsigned int order() const
Polynomial order of the shape functions.
Definition: raviartthomas1cube2dlocalbasis.hh:193
RT1Cube2DLocalBasis()
Standard constructor.
Definition: raviartthomas1cube2dlocalbasis.hh:33
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: raviartthomas1cube2dlocalbasis.hh:180
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:15
RT1Cube2DLocalBasis(unsigned int s)
Make set number s, where 0 <= s < 16.
Definition: raviartthomas1cube2dlocalbasis.hh:43
First order Raviart-Thomas shape functions on the reference quadrilateral.
Definition: raviartthomas1cube2dlocalbasis.hh:25
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: raviartthomas1cube2dlocalbasis.hh:113
D DomainType
domain type
Definition: localbasis.hh:49
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: raviartthomas1cube2dlocalbasis.hh:76
unsigned int size() const
number of shape functions
Definition: raviartthomas1cube2dlocalbasis.hh:65
Type traits for LocalBasisVirtualInterface.
Definition: localbasis.hh:37