3#ifndef DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI2_SIMPLEX2D_LOCALINTERPOLATION_HH
4#define DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI2_SIMPLEX2D_LOCALINTERPOLATION_HH
9#include <dune/localfunctions/common/localinterpolation.hh>
30 sign0 = sign1 = sign2 = 1.0;
40 sign0 = sign1 = sign2 = 1.0;
64 n2[0] = 1.0/sqrt(2.0);
65 n2[1] = 1.0/sqrt(2.0);
66 c0 = 0.5*n0[0] - 1.0*n0[1];
67 c1 = -1.0*n1[0] + 0.5*n1[1];
68 c2 = 0.5*n2[0] + 0.5*n2[1];
79 template<
typename F,
typename C>
83 typedef typename LB::Traits::RangeFieldType
Scalar;
84 typedef typename LB::Traits::DomainFieldType Vector;
86 auto&& f = Impl::makeFunctionWithCallOperator<typename LB::Traits::DomainType>(ff);
89 fill(out.begin(), out.end(), 0.0);
96 Scalar qPos = it->position();
98 typename LB::Traits::DomainType localPos;
102 auto y = f(localPos);
103 out[0] += (y[0]*n0[0] + y[1]*n0[1])*it->weight()*sign0/c0;
104 out[1] += (y[0]*n0[0] + y[1]*n0[1])*(1.0 - 2.0*qPos)*it->weight()/c0;
105 out[2] += (y[0]*n0[0] + y[1]*n0[1])*(6.0*qPos*qPos - 6.0*qPos + 1.0)*it->weight()*sign0/c0;
110 out[3] += (y[0]*n1[0]+y[1]*n1[1])*it->weight()*sign1/c1;
111 out[4] += (y[0]*n1[0]+y[1]*n1[1])*(2.0*qPos-1.0)*it->weight()/c1;
112 out[5] += (y[0]*n1[0]+y[1]*n1[1])*(6.0*qPos*qPos - 6.0*qPos + 1.0)*it->weight()*sign1/c1;
114 localPos[0] = 1.0 - qPos;
117 out[6] += (y[0]*n2[0] + y[1]*n2[1])*it->weight()*sign2/c2;
118 out[7] += (y[0]*n2[0] + y[1]*n2[1])*(1.0 - 2.0*qPos)*it->weight()/c2;
119 out[8] += (y[0]*n2[0] + y[1]*n2[1])*(6.0*qPos*qPos - 6.0*qPos + 1.0)*it->weight()*sign2/c2;
127 typename LB::Traits::DomainType localPos = it->position();
128 auto y = f(localPos);
130 out[9] += y[0]*it->weight();
131 out[10] += y[1]*it->weight();
132 out[11] += (y[0]*(localPos[0]-2.0*localPos[0]*localPos[1]-localPos[0]*localPos[0])
133 +y[1]*(-localPos[1]+2.0*localPos[0]*localPos[1]+localPos[1]*localPos[1]))*it->weight();
138 typename LB::Traits::RangeFieldType sign0, sign1, sign2;
139 typename LB::Traits::DomainType m0, m1, m2;
140 typename LB::Traits::DomainType n0, n1, n2;
141 typename LB::Traits::RangeFieldType c0, c1, c2;
First order Brezzi-Douglas-Marini shape functions on triangles.
Definition: brezzidouglasmarini2simplex2dlocalinterpolation.hh:24
BDM2Simplex2DLocalInterpolation()
Standard constructor.
Definition: brezzidouglasmarini2simplex2dlocalinterpolation.hh:28
void interpolate(const F &ff, std::vector< C > &out) const
Interpolate a given function with shape functions.
Definition: brezzidouglasmarini2simplex2dlocalinterpolation.hh:80
BDM2Simplex2DLocalInterpolation(unsigned int s)
Make set number s, where 0 <= s < 8.
Definition: brezzidouglasmarini2simplex2dlocalinterpolation.hh:38
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
constexpr GeometryType simplex(unsigned int dim)
Returns a GeometryType representing a simplex of dimension dim.
Definition: type.hh:766
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