Dune Core Modules (2.9.0)

raviartthomas3cube2dlocalinterpolation.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_RAVIARTTHOMAS3_CUBE2D_LOCALINTERPOLATION_HH
6 #define DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS3_CUBE2D_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:
29 
35  RT3Cube2DLocalInterpolation (std::bitset<4> s = 0)
36  {
37  for (size_t i=0; i<4; i++)
38  sign_[i] = (s[i]) ? -1.0 : 1.0;
39 
40  n_[0] = {-1.0, 0.0};
41  n_[1] = { 1.0, 0.0};
42  n_[2] = { 0.0, -1.0};
43  n_[3] = { 0.0, 1.0};
44  }
45 
54  template<typename F, typename C>
55  void interpolate (const F& ff, std::vector<C>& out) const
56  {
57  // f gives v*outer normal at a point on the edge!
58  typedef typename LB::Traits::RangeFieldType Scalar;
59  typedef typename LB::Traits::DomainFieldType Vector;
60 
61  auto&& f = Impl::makeFunctionWithCallOperator<typename LB::Traits::DomainType>(ff);
62 
63  out.resize(40);
64  fill(out.begin(), out.end(), 0.0);
65 
66  const int qOrder = 9;
67  const auto& rule1 = QuadratureRules<Scalar,1>::rule(GeometryTypes::cube(1), qOrder);
68 
69  for (auto&& qp : rule1)
70  {
71  Scalar qPos = qp.position();
72  typename LB::Traits::DomainType localPos;
73 
74  localPos = {0.0, qPos};
75  auto y = f(localPos);
76  out[0] += (y[0]*n_[0][0] + y[1]*n_[0][1])*qp.weight()*sign_[0];
77  out[1] += (y[0]*n_[0][0] + y[1]*n_[0][1])*(2.0*qPos - 1.0)*qp.weight();
78  out[2] += (y[0]*n_[0][0] + y[1]*n_[0][1])*(6.0*qPos*qPos - 6.0*qPos + 1.0)*qp.weight()*sign_[0];
79  out[3] += (y[0]*n_[0][0] + y[1]*n_[0][1])*(20.0*qPos*qPos*qPos - 30.0*qPos*qPos + 12.0*qPos - 1.0)*qp.weight();
80 
81  localPos = {1.0, qPos};
82  y = f(localPos);
83  out[4] += (y[0]*n_[1][0] + y[1]*n_[1][1])*qp.weight()*sign_[1];
84  out[5] += (y[0]*n_[1][0] + y[1]*n_[1][1])*(1.0 - 2.0*qPos)*qp.weight();
85  out[6] += (y[0]*n_[1][0] + y[1]*n_[1][1])*(6.0*qPos*qPos - 6.0*qPos + 1.0)*qp.weight()*sign_[1];
86  out[7] += (y[0]*n_[1][0] + y[1]*n_[1][1])*(-20.0*qPos*qPos*qPos + 30.0*qPos*qPos - 12.0*qPos + 1.0)*qp.weight();
87 
88  localPos = {qPos, 0.0};
89  y = f(localPos);
90  out[8] += (y[0]*n_[2][0] + y[1]*n_[2][1])*qp.weight()*sign_[2];
91  out[9] += (y[0]*n_[2][0] + y[1]*n_[2][1])*(1.0 - 2.0*qPos)*qp.weight();
92  out[10] += (y[0]*n_[2][0] + y[1]*n_[2][1])*(6.0*qPos*qPos - 6.0*qPos + 1.0)*qp.weight()*sign_[2];
93  out[11] += (y[0]*n_[2][0] + y[1]*n_[2][1])*(-20.0*qPos*qPos*qPos + 30.0*qPos*qPos - 12.0*qPos + 1.0)*qp.weight();
94 
95  localPos = {qPos, 1.0};
96  y = f(localPos);
97  out[12] += (y[0]*n_[3][0] + y[1]*n_[3][1])*qp.weight()*sign_[3];
98  out[13] += (y[0]*n_[3][0] + y[1]*n_[3][1])*(2.0*qPos - 1.0)*qp.weight();
99  out[14] += (y[0]*n_[3][0] + y[1]*n_[3][1])*(6.0*qPos*qPos - 6.0*qPos + 1.0)*qp.weight()*sign_[3];
100  out[15] += (y[0]*n_[3][0] + y[1]*n_[3][1])*(20.0*qPos*qPos*qPos - 30.0*qPos*qPos + 12.0*qPos - 1.0)*qp.weight();
101  }
102 
103  const auto& rule2 = QuadratureRules<Vector,2>::rule(GeometryTypes::cube(2), qOrder);
104 
105  for (auto&& qp : rule2)
106  {
107  auto qPos = qp.position();
108 
109  auto y = f(qPos);
110  double l0_x=1.0;
111  double l1_x=2.0*qPos[0]-1.0;
112  double l2_x=6.0*qPos[0]*qPos[0]-6.0*qPos[0]+1.0;
113  double l3_x=20.0*qPos[0]*qPos[0]*qPos[0] - 30.0*qPos[0]*qPos[0] + 12.0*qPos[0] - 1.0;
114  double l0_y=1.0;
115  double l1_y=2.0*qPos[1]-1.0;
116  double l2_y=6.0*qPos[1]*qPos[1]-6.0*qPos[1]+1.0;
117  double l3_y=20.0*qPos[1]*qPos[1]*qPos[1] - 30.0*qPos[1]*qPos[1] + 12.0*qPos[1] - 1.0;
118 
119  out[16] += y[0]*l0_x*l0_y*qp.weight();
120  out[17] += y[0]*l0_x*l1_y*qp.weight();
121  out[18] += y[0]*l0_x*l2_y*qp.weight();
122  out[19] += y[0]*l0_x*l3_y*qp.weight();
123  out[20] += y[0]*l1_x*l0_y*qp.weight();
124  out[21] += y[0]*l1_x*l1_y*qp.weight();
125  out[22] += y[0]*l1_x*l2_y*qp.weight();
126  out[23] += y[0]*l1_x*l3_y*qp.weight();
127  out[24] += y[0]*l2_x*l0_y*qp.weight();
128  out[25] += y[0]*l2_x*l1_y*qp.weight();
129  out[26] += y[0]*l2_x*l2_y*qp.weight();
130  out[27] += y[0]*l2_x*l3_y*qp.weight();
131 
132  out[28] += y[1]*l0_x*l0_y*qp.weight();
133  out[29] += y[1]*l0_x*l1_y*qp.weight();
134  out[30] += y[1]*l0_x*l2_y*qp.weight();
135  out[31] += y[1]*l1_x*l0_y*qp.weight();
136  out[32] += y[1]*l1_x*l1_y*qp.weight();
137  out[33] += y[1]*l1_x*l2_y*qp.weight();
138  out[34] += y[1]*l2_x*l0_y*qp.weight();
139  out[35] += y[1]*l2_x*l1_y*qp.weight();
140  out[36] += y[1]*l2_x*l2_y*qp.weight();
141  out[37] += y[1]*l3_x*l0_y*qp.weight();
142  out[38] += y[1]*l3_x*l1_y*qp.weight();
143  out[39] += y[1]*l3_x*l2_y*qp.weight();
144  }
145  }
146 
147  private:
148  // Edge orientations
149  std::array<typename LB::Traits::RangeFieldType, 4> sign_;
150 
151  // Edge normals
152  std::array<typename LB::Traits::DomainType, 4> n_;
153  };
154 }
155 
156 #endif // DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS3_CUBE2D_LOCALINTERPOLATION_HH
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
Second order Raviart-Thomas shape functions on the reference quadrilateral.
Definition: raviartthomas3cube2dlocalinterpolation.hh:26
RT3Cube2DLocalInterpolation(std::bitset< 4 > s=0)
Make set number s, where 0 <= s < 16.
Definition: raviartthomas3cube2dlocalinterpolation.hh:35
void interpolate(const F &ff, std::vector< C > &out) const
Interpolate a given function with shape functions.
Definition: raviartthomas3cube2dlocalinterpolation.hh:55
constexpr GeometryType cube(unsigned int dim)
Returns a GeometryType representing a hypercube of dimension dim.
Definition: type.hh:472
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)