Dune Core Modules (2.9.0)

brezzidouglasmarini2simplex2dlocalinterpolation.hh
1 // -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2 // vi: set et ts=4 sw=2 sts=2:
3 // SPDX-FileCopyrightInfo: Copyright (C) DUNE Project contributors, see file LICENSE.md in module root
4 // SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
5 #ifndef DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI2_SIMPLEX2D_LOCALINTERPOLATION_HH
6 #define DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI2_SIMPLEX2D_LOCALINTERPOLATION_HH
7 
8 #include <vector>
9 
11 #include <dune/localfunctions/common/localinterpolation.hh>
12 
13 namespace Dune
14 {
15 
24  template<class LB>
26  {
27 
28  public:
31  {
32  sign0 = sign1 = sign2 = 1.0;
33  }
34 
41  {
42  sign0 = sign1 = sign2 = 1.0;
43  if (s & 1)
44  {
45  sign0 = -1.0;
46  }
47  if (s & 2)
48  {
49  sign1 = -1.0;
50  }
51  if (s & 4)
52  {
53  sign2 = -1.0;
54  }
55 
56  m0[0] = 0.5;
57  m0[1] = 0.0;
58  m1[0] = 0.0;
59  m1[1] = 0.5;
60  m2[0] = 0.5;
61  m2[1] = 0.5;
62  n0[0] = 0.0;
63  n0[1] = -1.0;
64  n1[0] = -1.0;
65  n1[1] = 0.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];
71  }
72 
81  template<typename F, typename C>
82  void interpolate(const F& ff, std::vector<C>& out) const
83  {
84  // f gives v*outer normal at a point on the edge!
85  typedef typename LB::Traits::RangeFieldType Scalar;
86  typedef typename LB::Traits::DomainFieldType Vector;
87 
88  auto&& f = Impl::makeFunctionWithCallOperator<typename LB::Traits::DomainType>(ff);
89 
90  out.resize(12);
91  fill(out.begin(), out.end(), 0.0);
92 
93  const int qOrder = 4;
95 
96  for (typename Dune::QuadratureRule<Scalar,1>::const_iterator it=rule.begin(); it!=rule.end(); ++it)
97  {
98  Scalar qPos = it->position();
99 
100  typename LB::Traits::DomainType localPos;
101 
102  localPos[0] = qPos;
103  localPos[1] = 0.0;
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;
108 
109  localPos[0] = 0.0;
110  localPos[1] = qPos;
111  y = f(localPos);
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;
115 
116  localPos[0] = 1.0 - qPos;
117  localPos[1] = qPos;
118  y = f(localPos);
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;
122  }
123 
124  // a volume part is needed here for dofs: 9 10 11
126 
127  for (typename QuadratureRule<Vector,2>::const_iterator it=rule2.begin(); it!=rule2.end(); ++it)
128  {
129  typename LB::Traits::DomainType localPos = it->position();
130  auto y = f(localPos);
131 
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();
136  }
137  }
138 
139  private:
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;
144  };
145 } // end namespace Dune
146 #endif // DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI2_SIMPLEX2D_LOCALINTERPOLATION_HH
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
Abstract base class for quadrature rules.
Definition: quadraturerules.hh:154
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:266
constexpr GeometryType simplex(unsigned int dim)
Returns a GeometryType representing a simplex of dimension dim.
Definition: type.hh:463
typename Overloads::ScalarType< std::decay_t< V > >::type Scalar
Element type of some SIMD type.
Definition: interface.hh:235
Dune namespace.
Definition: alignedallocator.hh:13
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (May 2, 22:35, 2024)