Dune Core Modules (2.9.0)

raviartthomas0pyramidlocalinterpolation.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_RAVIARTTHOMAS0_PYRAMID_LOCALINTERPOLATION_HH
6 #define DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS0_PYRAMID_LOCALINTERPOLATION_HH
7 
8 #include <vector>
9 
10 #include <dune/localfunctions/common/localinterpolation.hh>
11 
12 namespace Dune
13 {
22  template<class LB>
24  {
25 
26  public:
27 
33  RT0PyramidLocalInterpolation (std::bitset<5> s = 0)
34  {
35  typedef typename LB::Traits::RangeFieldType Scalar;
36 
37  for (size_t i=0; i<5; i++)
38  sign[i] = (s[i]) ? -1.0 : 1.0;
39 
40  Scalar r = 1/std::sqrt(2);
41 
42  n[0] = { 0.0, 0.0, -1.0};
43  n[1] = {-1.0, 0.0, 0.0};
44  n[2] = { r, 0.0, r};
45  n[3] = { 0.0, -1.0, 0.0};
46  n[4] = { 0.0, r, r};
47 
48  c[0] = 1.0;
49  c[1] = 1/2.0;
50  c[2] = 1/2.0 * std::sqrt(2);
51  c[3] = 1/2.0;
52  c[4] = 1/2.0 * std::sqrt(2);
53 
54  m[0] = { 0.5, 0.5, 0.0};
55  m[1] = { 0.0, 1/3.0, 1/3.0};
56  m[2] = { 2/3.0, 1/3.0, 1/3.0};
57  m[3] = { 1/3.0, 0.0, 1/3.0};
58  m[4] = { 1/3.0, 2/3.0, 1/3.0};
59  }
60 
69  template<class F, class C>
70  void interpolate (const F& ff, std::vector<C>& out) const
71  {
72  auto&& f = Impl::makeFunctionWithCallOperator<typename LB::Traits::DomainType>(ff);
73 
74  out.resize(5);
75  for(int i=0; i<5; i++)
76  out[i] = f(m[i]).dot(n[i]) * c[i] * sign[i];
77  }
78 
79  private:
80  // Facet orientations
81  std::array<typename LB::Traits::RangeFieldType, 5> sign;
82  // Facet area
83  std::array<typename LB::Traits::RangeFieldType, 5> c;
84 
85  // Facet normals
86  std::array<typename LB::Traits::DomainType, 5> n;
87  // Facet midpoints
88  std::array<typename LB::Traits::DomainType, 5> m;
89  };
90 }
91 #endif // DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS0_PYRAMID_LOCALINTERPOLATION_HH
First order Raviart-Thomas shape functions on the reference hexahedron.
Definition: raviartthomas0pyramidlocalinterpolation.hh:24
void interpolate(const F &ff, std::vector< C > &out) const
Interpolate a given function with shape functions.
Definition: raviartthomas0pyramidlocalinterpolation.hh:70
RT0PyramidLocalInterpolation(std::bitset< 5 > s=0)
Make set number s, where 0 <= s < 32.
Definition: raviartthomas0pyramidlocalinterpolation.hh:33
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 (Apr 26, 22:29, 2024)