dune-localfunctions  2.7.0
raviartthomas12dlocalinterpolation.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_RAVIARTTHOMAS12DLOCALINTERPOLATION_HH
4 #define DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS12DLOCALINTERPOLATION_HH
5 
6 #include <vector>
7 
8 #include <dune/geometry/quadraturerules.hh>
10 
11 namespace Dune
12 {
13 
22  template<class LB>
24  {
25 
26  public:
29  {
30  sign0 = sign1 = sign2 = 1.0;
31  }
32 
38  RT12DLocalInterpolation (unsigned int s)
39  {
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  n0[0] = 0.0;
54  n0[1] = -1.0;
55  n1[0] = -1.0;
56  n1[1] = 0.0;
57  n2[0] = 1.0/sqrt(2.0);
58  n2[1] = 1.0/sqrt(2.0);
59  c0 = 0.5*n0[0] - 1.0*n0[1];
60  c1 = -1.0*n1[0] + 0.5*n1[1];
61  c2 = 0.5*n2[0] + 0.5*n2[1];
62  }
63 
72  template<typename F, typename C>
73  void interpolate (const F& ff, std::vector<C>& out) const
74  {
75  // f gives v*outer normal at a point on the edge!
76  typedef typename LB::Traits::RangeFieldType Scalar;
77  typedef typename LB::Traits::DomainFieldType Vector;
78 
79  auto&& f = Impl::makeFunctionWithCallOperator<typename LB::Traits::DomainType>(ff);
80 
81  out.resize(8);
82  fill(out.begin(), out.end(), 0.0);
83 
84  const int qOrder1 = 4;
85  const Dune::QuadratureRule<Scalar,1>& rule1 = Dune::QuadratureRules<Scalar,1>::rule(Dune::GeometryTypes::simplex(1), qOrder1);
86 
87  for (typename Dune::QuadratureRule<Scalar,1>::const_iterator it = rule1.begin();
88  it != rule1.end(); ++it)
89  {
90  Scalar qPos = it->position();
91  typename LB::Traits::DomainType localPos;
92 
93  localPos[0] = qPos;
94  localPos[1] = 0.0;
95  auto y = f(localPos);
96  out[0] += (y[0]*n0[0] + y[1]*n0[1])*it->weight()*sign0/c0;
97  out[3] += (y[0]*n0[0] + y[1]*n0[1])*(2.0*qPos - 1.0)*it->weight()/c0;
98 
99  localPos[0] = 0.0;
100  localPos[1] = qPos;
101  y = f(localPos);
102  out[1] += (y[0]*n1[0] + y[1]*n1[1])*it->weight()*sign1/c1;
103  out[4] += (y[0]*n1[0] + y[1]*n1[1])*(1.0 - 2.0*qPos)*it->weight()/c1;
104 
105  localPos[0] = 1.0 - qPos;
106  localPos[1] = qPos;
107  y = f(localPos);
108  out[2] += (y[0]*n2[0] + y[1]*n2[1])*it->weight()*sign2/c2;
109  out[5] += (y[0]*n2[0] + y[1]*n2[1])*(2.0*qPos - 1.0)*it->weight()/c2;
110  }
111 
112  const int qOrder2 = 8;
113  const Dune::QuadratureRule<Vector,2>& rule2 = Dune::QuadratureRules<Vector,2>::rule(Dune::GeometryTypes::simplex(2), qOrder2);
114 
115  for (typename Dune::QuadratureRule<Vector,2>::const_iterator it = rule2.begin();
116  it != rule2.end(); ++it)
117  {
118  Dune::FieldVector<double,2> qPos = it->position();
119 
120  auto y = f(qPos);
121  out[6] += y[0]*it->weight();
122  out[7] += y[1]*it->weight();
123  }
124  }
125 
126  private:
127  typename LB::Traits::RangeFieldType sign0,sign1,sign2;
128  typename LB::Traits::DomainType n0,n1,n2;
129  typename LB::Traits::RangeFieldType c0,c1,c2;
130  };
131 }
132 #endif // DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS12DLOCALINTERPOLATION_HH
localinterpolation.hh
Dune
Definition: bdfmcube.hh:15
Dune::RT12DLocalInterpolation::interpolate
void interpolate(const F &ff, std::vector< C > &out) const
Interpolate a given function with shape functions.
Definition: raviartthomas12dlocalinterpolation.hh:73
Dune::RT12DLocalInterpolation::RT12DLocalInterpolation
RT12DLocalInterpolation(unsigned int s)
Make set number s, where 0 <= s < 8.
Definition: raviartthomas12dlocalinterpolation.hh:38
Dune::RT12DLocalInterpolation
First order Raviart-Thomas shape functions on the reference quadrilateral.
Definition: raviartthomas12dlocalinterpolation.hh:23
Dune::RT12DLocalInterpolation::RT12DLocalInterpolation
RT12DLocalInterpolation()
Standard constructor.
Definition: raviartthomas12dlocalinterpolation.hh:28