Dune Core Modules (2.7.1)

brezzidouglasmarini1simplex2dlocalinterpolation.hh
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 
9 #include <dune/localfunctions/common/localinterpolation.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  sign0 = sign1 = sign2 = 1.0;
40  if (s & 1)
41  {
42  sign0 = -1.0;
43  }
44  if (s & 2)
45  {
46  sign1 = -1.0;
47  }
48  if (s & 4)
49  {
50  sign2 = -1.0;
51  }
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 
78  auto&& f = Impl::makeFunctionWithCallOperator<typename LB::Traits::DomainType>(ff);
79 
80  out.resize(6);
81  fill(out.begin(), out.end(), 0.0);
82 
83  const int qOrder = 4;
85 
86  for (typename Dune::QuadratureRule<Scalar,1>::const_iterator it=rule.begin(); it!=rule.end(); ++it)
87  {
88  Scalar qPos = it->position();
89  typename LB::Traits::DomainType localPos;
90 
91  localPos[0] = qPos;
92  localPos[1] = 0.0;
93  auto y = f(localPos);
94  out[0] += (y[0]*n0[0] + y[1]*n0[1])*it->weight()*sign0/c0;
95  out[3] += (y[0]*n0[0] + y[1]*n0[1])*(2.0*qPos - 1.0)*it->weight()/c0;
96 
97  localPos[0] = 0.0;
98  localPos[1] = qPos;
99  y = f(localPos);
100  out[1] += (y[0]*n1[0] + y[1]*n1[1])*it->weight()*sign1/c1;
101  out[4] += (y[0]*n1[0] + y[1]*n1[1])*(1.0 - 2.0*qPos)*it->weight()/c1;
102 
103  localPos[0] = 1.0 - qPos;
104  localPos[1] = qPos;
105  y = f(localPos);
106  out[2] += (y[0]*n2[0] + y[1]*n2[1])*it->weight()*sign2/c2;
107  out[5] += (y[0]*n2[0] + y[1]*n2[1])*(2.0*qPos - 1.0)*it->weight()/c2;
108  }
109  }
110 
111  private:
112  typename LB::Traits::RangeFieldType sign0,sign1,sign2;
113  typename LB::Traits::DomainType n0,n1,n2;
114  typename LB::Traits::RangeFieldType c0,c1,c2;
115  };
116 }
117 
118 #endif // DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI1_SIMPLEX2D_LOCALINTERPOLATION_HH
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:73
BDM1Simplex2DLocalInterpolation(unsigned int s)
Make set number s, where 0 <= s < 8.
Definition: brezzidouglasmarini1simplex2dlocalinterpolation.hh:37
BDM1Simplex2DLocalInterpolation()
Standard constructor.
Definition: brezzidouglasmarini1simplex2dlocalinterpolation.hh:27
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
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (May 16, 22:29, 2024)