5 #ifndef DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI2_SIMPLEX2D_LOCALINTERPOLATION_HH
6 #define DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI2_SIMPLEX2D_LOCALINTERPOLATION_HH
10 #include <dune/geometry/quadraturerules.hh>
32 sign0 = sign1 = sign2 = 1.0;
42 sign0 = sign1 = sign2 = 1.0;
66 n2[0] = 1.0/sqrt(2.0);
67 n2[1] = 1.0/sqrt(2.0);
68 c0 = 0.5*n0[0] - 1.0*n0[1];
69 c1 = -1.0*n1[0] + 0.5*n1[1];
70 c2 = 0.5*n2[0] + 0.5*n2[1];
81 template<
typename F,
typename C>
85 typedef typename LB::Traits::RangeFieldType Scalar;
86 typedef typename LB::Traits::DomainFieldType Vector;
88 auto&& f = Impl::makeFunctionWithCallOperator<typename LB::Traits::DomainType>(ff);
91 fill(out.begin(), out.end(), 0.0);
94 const Dune::QuadratureRule<Scalar,1>& rule = Dune::QuadratureRules<Scalar,1>::rule(Dune::GeometryTypes::simplex(1), qOrder);
96 for (
typename Dune::QuadratureRule<Scalar,1>::const_iterator it=rule.begin(); it!=rule.end(); ++it)
98 Scalar qPos = it->position();
100 typename LB::Traits::DomainType localPos;
104 auto y = f(localPos);
105 out[0] += (y[0]*n0[0] + y[1]*n0[1])*it->weight()*sign0/c0;
106 out[1] += (y[0]*n0[0] + y[1]*n0[1])*(1.0 - 2.0*qPos)*it->weight()/c0;
107 out[2] += (y[0]*n0[0] + y[1]*n0[1])*(6.0*qPos*qPos - 6.0*qPos + 1.0)*it->weight()*sign0/c0;
112 out[3] += (y[0]*n1[0]+y[1]*n1[1])*it->weight()*sign1/c1;
113 out[4] += (y[0]*n1[0]+y[1]*n1[1])*(2.0*qPos-1.0)*it->weight()/c1;
114 out[5] += (y[0]*n1[0]+y[1]*n1[1])*(6.0*qPos*qPos - 6.0*qPos + 1.0)*it->weight()*sign1/c1;
116 localPos[0] = 1.0 - qPos;
119 out[6] += (y[0]*n2[0] + y[1]*n2[1])*it->weight()*sign2/c2;
120 out[7] += (y[0]*n2[0] + y[1]*n2[1])*(1.0 - 2.0*qPos)*it->weight()/c2;
121 out[8] += (y[0]*n2[0] + y[1]*n2[1])*(6.0*qPos*qPos - 6.0*qPos + 1.0)*it->weight()*sign2/c2;
125 const QuadratureRule<Vector,2>& rule2 = QuadratureRules<Vector,2>::rule(GeometryTypes::simplex(2), qOrder);
127 for (
typename QuadratureRule<Vector,2>::const_iterator it=rule2.begin(); it!=rule2.end(); ++it)
129 typename LB::Traits::DomainType localPos = it->position();
130 auto y = f(localPos);
132 out[9] += y[0]*it->weight();
133 out[10] += y[1]*it->weight();
134 out[11] += (y[0]*(localPos[0]-2.0*localPos[0]*localPos[1]-localPos[0]*localPos[0])
135 +y[1]*(-localPos[1]+2.0*localPos[0]*localPos[1]+localPos[1]*localPos[1]))*it->weight();
140 typename LB::Traits::RangeFieldType sign0, sign1, sign2;
141 typename LB::Traits::DomainType m0, m1, m2;
142 typename LB::Traits::DomainType n0, n1, n2;
143 typename LB::Traits::RangeFieldType c0, c1, c2;
Definition: bdfmcube.hh:18
First order Brezzi-Douglas-Marini shape functions on triangles.
Definition: brezzidouglasmarini2simplex2dlocalinterpolation.hh:26
BDM2Simplex2DLocalInterpolation()
Standard constructor.
Definition: brezzidouglasmarini2simplex2dlocalinterpolation.hh:30
void interpolate(const F &ff, std::vector< C > &out) const
Interpolate a given function with shape functions.
Definition: brezzidouglasmarini2simplex2dlocalinterpolation.hh:82
BDM2Simplex2DLocalInterpolation(unsigned int s)
Make set number s, where 0 <= s < 8.
Definition: brezzidouglasmarini2simplex2dlocalinterpolation.hh:40