3#ifndef DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS3_CUBE2D_LOCALINTERPOLATION_HH
4#define DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS3_CUBE2D_LOCALINTERPOLATION_HH
9#include <dune/localfunctions/common/localinterpolation.hh>
35 sign0 = sign1 = sign2 = sign3 = 1.0;
71 template<
typename F,
typename C>
75 typedef typename LB::Traits::RangeFieldType
Scalar;
76 typedef typename LB::Traits::DomainFieldType Vector;
78 auto&& f = Impl::makeFunctionWithCallOperator<typename LB::Traits::DomainType>(ff);
81 fill(out.begin(), out.end(), 0.0);
88 Scalar qPos = it->position();
89 typename LB::Traits::DomainType localPos;
94 out[0] += (y[0]*n0[0] + y[1]*n0[1])*it->weight()*sign0;
95 out[1] += (y[0]*n0[0] + y[1]*n0[1])*(2.0*qPos - 1.0)*it->weight();
96 out[2] += (y[0]*n0[0] + y[1]*n0[1])*(6.0*qPos*qPos - 6.0*qPos + 1.0)*it->weight()*sign0;
97 out[3] += (y[0]*n0[0] + y[1]*n0[1])*(20.0*qPos*qPos*qPos - 30.0*qPos*qPos + 12.0*qPos - 1.0)*it->weight();
102 out[4] += (y[0]*n1[0] + y[1]*n1[1])*it->weight()*sign1;
103 out[5] += (y[0]*n1[0] + y[1]*n1[1])*(1.0 - 2.0*qPos)*it->weight();
104 out[6] += (y[0]*n1[0] + y[1]*n1[1])*(6.0*qPos*qPos - 6.0*qPos + 1.0)*it->weight()*sign1;
105 out[7] += (y[0]*n1[0] + y[1]*n1[1])*(-20.0*qPos*qPos*qPos + 30.0*qPos*qPos - 12.0*qPos + 1.0)*it->weight();
110 out[8] += (y[0]*n2[0] + y[1]*n2[1])*it->weight()*sign2;
111 out[9] += (y[0]*n2[0] + y[1]*n2[1])*(1.0 - 2.0*qPos)*it->weight();
112 out[10] += (y[0]*n2[0] + y[1]*n2[1])*(6.0*qPos*qPos - 6.0*qPos + 1.0)*it->weight()*sign2;
113 out[11] += (y[0]*n2[0] + y[1]*n2[1])*(-20.0*qPos*qPos*qPos + 30.0*qPos*qPos - 12.0*qPos + 1.0)*it->weight();
118 out[12] += (y[0]*n3[0] + y[1]*n3[1])*it->weight()*sign3;
119 out[13] += (y[0]*n3[0] + y[1]*n3[1])*(2.0*qPos - 1.0)*it->weight();
120 out[14] += (y[0]*n3[0] + y[1]*n3[1])*(6.0*qPos*qPos - 6.0*qPos + 1.0)*it->weight()*sign3;
121 out[15] += (y[0]*n3[0] + y[1]*n3[1])*(20.0*qPos*qPos*qPos - 30.0*qPos*qPos + 12.0*qPos - 1.0)*it->weight();
127 it != rule2.end(); ++it)
133 double l1_x=2.0*qPos[0]-1.0;
134 double l2_x=6.0*qPos[0]*qPos[0]-6.0*qPos[0]+1.0;
135 double l3_x=20.0*qPos[0]*qPos[0]*qPos[0] - 30.0*qPos[0]*qPos[0] + 12.0*qPos[0] - 1.0;
137 double l1_y=2.0*qPos[1]-1.0;
138 double l2_y=6.0*qPos[1]*qPos[1]-6.0*qPos[1]+1.0;
139 double l3_y=20.0*qPos[1]*qPos[1]*qPos[1] - 30.0*qPos[1]*qPos[1] + 12.0*qPos[1] - 1.0;
141 out[16] += y[0]*l0_x*l0_y*it->weight();
142 out[17] += y[0]*l0_x*l1_y*it->weight();
143 out[18] += y[0]*l0_x*l2_y*it->weight();
144 out[19] += y[0]*l0_x*l3_y*it->weight();
145 out[20] += y[0]*l1_x*l0_y*it->weight();
146 out[21] += y[0]*l1_x*l1_y*it->weight();
147 out[22] += y[0]*l1_x*l2_y*it->weight();
148 out[23] += y[0]*l1_x*l3_y*it->weight();
149 out[24] += y[0]*l2_x*l0_y*it->weight();
150 out[25] += y[0]*l2_x*l1_y*it->weight();
151 out[26] += y[0]*l2_x*l2_y*it->weight();
152 out[27] += y[0]*l2_x*l3_y*it->weight();
154 out[28] += y[1]*l0_x*l0_y*it->weight();
155 out[29] += y[1]*l0_x*l1_y*it->weight();
156 out[30] += y[1]*l0_x*l2_y*it->weight();
157 out[31] += y[1]*l1_x*l0_y*it->weight();
158 out[32] += y[1]*l1_x*l1_y*it->weight();
159 out[33] += y[1]*l1_x*l2_y*it->weight();
160 out[34] += y[1]*l2_x*l0_y*it->weight();
161 out[35] += y[1]*l2_x*l1_y*it->weight();
162 out[36] += y[1]*l2_x*l2_y*it->weight();
163 out[37] += y[1]*l3_x*l0_y*it->weight();
164 out[38] += y[1]*l3_x*l1_y*it->weight();
165 out[39] += y[1]*l3_x*l2_y*it->weight();
170 typename LB::Traits::RangeFieldType sign0, sign1, sign2, sign3;
171 typename LB::Traits::DomainType n0, n1, n2, n3;
vector space out of a tensor product of fields.
Definition: fvector.hh:96
Abstract base class for quadrature rules.
Definition: quadraturerules.hh:126
static const QuadratureRule & rule(const GeometryType &t, int p, QuadratureType::Enum qt=QuadratureType::GaussLegendre)
select the appropriate QuadratureRule for GeometryType t and order p
Definition: quadraturerules.hh:254
Second order Raviart-Thomas shape functions on the reference quadrilateral.
Definition: raviartthomas3cube2dlocalinterpolation.hh:24
RT3Cube2DLocalInterpolation(unsigned int s=0)
Make set number s, where 0 <= s < 16.
Definition: raviartthomas3cube2dlocalinterpolation.hh:33
void interpolate(const F &ff, std::vector< C > &out) const
Interpolate a given function with shape functions.
Definition: raviartthomas3cube2dlocalinterpolation.hh:72
constexpr GeometryType cube(unsigned int dim)
Returns a GeometryType representing a hypercube of dimension dim.
Definition: type.hh:775
typename Overloads::ScalarType< std::decay_t< V > >::type Scalar
Element type of some SIMD type.
Definition: interface.hh:233
Dune namespace.
Definition: alignedallocator.hh:14