dune-localfunctions  2.8.0
brezzidouglasmarini1simplex2dlocalinterpolation.hh
Go to the documentation of this file.
1 // -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2 // vi: set et ts=4 sw=2 sts=2:
3 #ifndef DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI1_SIMPLEX2D_LOCALINTERPOLATION_HH
4 #define DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI1_SIMPLEX2D_LOCALINTERPOLATION_HH
5 
6 #include <vector>
7 
8 #include <dune/geometry/quadraturerules.hh>
10 
11 namespace Dune
12 {
21  template<class LB>
23  {
24 
25  public:
28  {
29  sign0 = sign1 = sign2 = 1.0;
30  }
31 
38  {
39  using std::sqrt;
40  sign0 = sign1 = sign2 = 1.0;
41  if (s & 1)
42  {
43  sign0 = -1.0;
44  }
45  if (s & 2)
46  {
47  sign1 = -1.0;
48  }
49  if (s & 4)
50  {
51  sign2 = -1.0;
52  }
53 
54  n0[0] = 0.0;
55  n0[1] = -1.0;
56  n1[0] = -1.0;
57  n1[1] = 0.0;
58  n2[0] = 1.0/sqrt(2.0);
59  n2[1] = 1.0/sqrt(2.0);
60  c0 = 0.5*n0[0] - 1.0*n0[1];
61  c1 = -1.0*n1[0] + 0.5*n1[1];
62  c2 = 0.5*n2[0] + 0.5*n2[1];
63  }
64 
73  template<typename F, typename C>
74  void interpolate (const F& ff, std::vector<C>& out) const
75  {
76  // f gives v*outer normal at a point on the edge!
77  typedef typename LB::Traits::RangeFieldType Scalar;
78 
79  auto&& f = Impl::makeFunctionWithCallOperator<typename LB::Traits::DomainType>(ff);
80 
81  out.resize(6);
82  fill(out.begin(), out.end(), 0.0);
83 
84  const int qOrder = 4;
85  const Dune::QuadratureRule<Scalar,1>& rule = Dune::QuadratureRules<Scalar,1>::rule(Dune::GeometryTypes::simplex(1), qOrder);
86 
87  for (typename Dune::QuadratureRule<Scalar,1>::const_iterator it=rule.begin(); it!=rule.end(); ++it)
88  {
89  Scalar qPos = it->position();
90  typename LB::Traits::DomainType localPos;
91 
92  localPos[0] = qPos;
93  localPos[1] = 0.0;
94  auto y = f(localPos);
95  out[0] += (y[0]*n0[0] + y[1]*n0[1])*it->weight()*sign0/c0;
96  out[3] += (y[0]*n0[0] + y[1]*n0[1])*(2.0*qPos - 1.0)*it->weight()/c0;
97 
98  localPos[0] = 0.0;
99  localPos[1] = qPos;
100  y = f(localPos);
101  out[1] += (y[0]*n1[0] + y[1]*n1[1])*it->weight()*sign1/c1;
102  out[4] += (y[0]*n1[0] + y[1]*n1[1])*(1.0 - 2.0*qPos)*it->weight()/c1;
103 
104  localPos[0] = 1.0 - qPos;
105  localPos[1] = qPos;
106  y = f(localPos);
107  out[2] += (y[0]*n2[0] + y[1]*n2[1])*it->weight()*sign2/c2;
108  out[5] += (y[0]*n2[0] + y[1]*n2[1])*(2.0*qPos - 1.0)*it->weight()/c2;
109  }
110  }
111 
112  private:
113  typename LB::Traits::RangeFieldType sign0,sign1,sign2;
114  typename LB::Traits::DomainType n0,n1,n2;
115  typename LB::Traits::RangeFieldType c0,c1,c2;
116  };
117 }
118 
119 #endif // DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI1_SIMPLEX2D_LOCALINTERPOLATION_HH
Definition: bdfmcube.hh:16
First order Brezzi-Douglas-Marini shape functions on the reference triangle.
Definition: brezzidouglasmarini1simplex2dlocalinterpolation.hh:23
void interpolate(const F &ff, std::vector< C > &out) const
Interpolate a given function with shape functions.
Definition: brezzidouglasmarini1simplex2dlocalinterpolation.hh:74
BDM1Simplex2DLocalInterpolation(unsigned int s)
Make set number s, where 0 <= s < 8.
Definition: brezzidouglasmarini1simplex2dlocalinterpolation.hh:37
BDM1Simplex2DLocalInterpolation()
Standard constructor.
Definition: brezzidouglasmarini1simplex2dlocalinterpolation.hh:27