3#ifndef DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS1_CUBE3D_LOCALINTERPOLATION_HH
4#define DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS1_CUBE3D_LOCALINTERPOLATION_HH
9#include <dune/localfunctions/common/localinterpolation.hh>
34 sign0 = sign1 = sign2 = sign3 = sign4 = sign5 = 1.0;
88 template<
class F,
class C>
92 typedef typename LB::Traits::RangeFieldType
Scalar;
93 typedef typename LB::Traits::DomainFieldType Vector;
95 auto&& f = Impl::makeFunctionWithCallOperator<typename LB::Traits::DomainType>(ff);
98 fill(out.begin(), out.end(), 0.0);
100 const int qOrder = 3;
104 it != rule1.end(); ++it)
107 typename LB::Traits::DomainType localPos;
110 localPos[1] = qPos[0];
111 localPos[2] = qPos[1];
112 auto y = f(localPos);
113 out[0] += (y[0]*n0[0] + y[1]*n0[1] + y[2]*n0[2])*it->weight()*sign0;
114 out[6] += (y[0]*n0[0] + y[1]*n0[1] + y[2]*n0[2])*(2.0*qPos[0] - 1.0)*it->weight();
115 out[12] += (y[0]*n0[0] + y[1]*n0[1] + y[2]*n0[2])*(2.0*qPos[1] - 1.0)*it->weight();
116 out[18] += (y[0]*n0[0] + y[1]*n0[1] + y[2]*n0[2])*(2.0*qPos[0] - 1.0)*(2.0*qPos[1] - 1.0)*it->weight();
119 localPos[1] = qPos[0];
120 localPos[2] = qPos[1];
122 out[1] += (y[0]*n1[0] + y[1]*n1[1] + y[2]*n1[2])*it->weight()*sign1;
123 out[7] += (y[0]*n1[0] + y[1]*n1[1] + y[2]*n1[2])*(1.0 - 2.0*qPos[0])*it->weight();
124 out[13] += (y[0]*n1[0] + y[1]*n1[1] + y[2]*n1[2])*(1.0 - 2.0*qPos[1])*it->weight();
125 out[19] += (y[0]*n1[0] + y[1]*n1[1] + y[2]*n1[2])*(1.0 - 2.0*qPos[0])*(2.0*qPos[1] - 1.0)*it->weight();
127 localPos[0] = qPos[0];
129 localPos[2] = qPos[1];
131 out[2] += (y[0]*n2[0] + y[1]*n2[1] + y[2]*n2[2])*it->weight()*sign2;
132 out[8] += (y[0]*n2[0] + y[1]*n2[1] + y[2]*n2[2])*(1.0 - 2.0*qPos[0])*it->weight();
133 out[14] += (y[0]*n2[0] + y[1]*n2[1] + y[2]*n2[2])*(2.0*qPos[1] - 1.0)*it->weight();
134 out[20] += (y[0]*n2[0] + y[1]*n2[1] + y[2]*n2[2])*(1.0 - 2.0*qPos[0])*(2.0*qPos[1] - 1.0)*it->weight();
136 localPos[0] = qPos[0];
138 localPos[2] = qPos[1];
140 out[3] += (y[0]*n3[0] + y[1]*n3[1] + y[2]*n3[2])*it->weight()*sign3;
141 out[9] += (y[0]*n3[0] + y[1]*n3[1] + y[2]*n3[2])*(2.0*qPos[0] - 1.0)*it->weight();
142 out[15] += (y[0]*n3[0] + y[1]*n3[1] + y[2]*n3[2])*(1.0 - 2.0*qPos[1])*it->weight();
143 out[21] += (y[0]*n3[0] + y[1]*n3[1] + y[2]*n3[2])*(2.0*qPos[0] - 1.0)*(2.0*qPos[1] - 1.0)*it->weight();
145 localPos[0] = qPos[0];
146 localPos[1] = qPos[1];
149 out[4] += (y[0]*n4[0] + y[1]*n4[1] + y[2]*n4[2])*it->weight()*sign4;
150 out[10] += (y[0]*n4[0] + y[1]*n4[1] + y[2]*n4[2])*(1.0 - 2.0*qPos[0])*it->weight();
151 out[16] += (y[0]*n4[0] + y[1]*n4[1] + y[2]*n4[2])*(1.0 - 2.0*qPos[1])*it->weight();
152 out[22] += (y[0]*n4[0] + y[1]*n4[1] + y[2]*n4[2])*(1.0 - 2.0*qPos[0])*(2.0*qPos[1] - 1.0)*it->weight();
154 localPos[0] = qPos[0];
155 localPos[1] = qPos[1];
158 out[5] += (y[0]*n5[0] + y[1]*n5[1] + y[2]*n5[2])*it->weight()*sign5;
159 out[11] += (y[0]*n5[0] + y[1]*n5[1] + y[2]*n5[2])*(2.0*qPos[0] - 1.0)*it->weight();
160 out[17] += (y[0]*n5[0] + y[1]*n5[1] + y[2]*n5[2])*(2.0*qPos[1] - 1.0)*it->weight();
161 out[23] += (y[0]*n5[0] + y[1]*n5[1] + y[2]*n5[2])*(2.0*qPos[0] - 1.0)*(2.0*qPos[1] - 1.0)*it->weight();
166 it != rule2.end(); ++it)
171 out[24] += y[0]*it->weight();
172 out[25] += y[1]*it->weight();
173 out[26] += y[2]*it->weight();
174 out[27] += y[0]*qPos[1]*it->weight();
175 out[28] += y[0]*qPos[2]*it->weight();
176 out[29] += y[1]*qPos[0]*it->weight();
177 out[30] += y[1]*qPos[2]*it->weight();
178 out[31] += y[2]*qPos[0]*it->weight();
179 out[32] += y[2]*qPos[1]*it->weight();
180 out[33] += y[0]*qPos[1]*qPos[2]*it->weight();
181 out[34] += y[1]*qPos[0]*qPos[2]*it->weight();
182 out[35] += y[2]*qPos[0]*qPos[1]*it->weight();
187 typename LB::Traits::RangeFieldType sign0, sign1, sign2, sign3, sign4, sign5;
188 typename LB::Traits::DomainType n0, n1, n2, n3, n4, n5;
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
First order Raviart-Thomas shape functions on the reference hexahedron.
Definition: raviartthomas1cube3dlocalinterpolation.hh:23
RT1Cube3DLocalInterpolation(unsigned int s=0)
Make set number s, where 0 <= s < 64.
Definition: raviartthomas1cube3dlocalinterpolation.hh:32
void interpolate(const F &ff, std::vector< C > &out) const
Interpolate a given function with shape functions.
Definition: raviartthomas1cube3dlocalinterpolation.hh:89
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