Dune Core Modules (2.9.0)

raviartthomas12dlocalinterpolation.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_RAVIARTTHOMAS12DLOCALINTERPOLATION_HH
6#define DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS12DLOCALINTERPOLATION_HH
7
8#include <vector>
9
11#include <dune/localfunctions/common/localinterpolation.hh>
12
13namespace Dune
14{
15
24 template<class LB>
26 {
27
28 public:
29
35 RT12DLocalInterpolation (std::bitset<3> s = 0)
36 {
37 using std::sqrt;
38 for (size_t i=0; i<3; i++)
39 sign_[i] = (s[i]) ? -1.0 : 1.0;
40
41 n_[0] = { 0.0, -1.0};
42 n_[1] = {-1.0, 0.0};
43 n_[2] = { 1.0/sqrt(2.0), 1.0/sqrt(2.0)};
44
45 c_ = { 0.5*n_[0][0] - 1.0*n_[0][1],
46 -1.0*n_[1][0] + 0.5*n_[1][1],
47 0.5*n_[2][0] + 0.5*n_[2][1]};
48 }
49
58 template<typename F, typename C>
59 void interpolate (const F& ff, std::vector<C>& out) const
60 {
61 // f gives v*outer normal at a point on the edge!
62 typedef typename LB::Traits::RangeFieldType Scalar;
63 typedef typename LB::Traits::DomainFieldType Vector;
64
65 auto&& f = Impl::makeFunctionWithCallOperator<typename LB::Traits::DomainType>(ff);
66
67 out.resize(8);
68 fill(out.begin(), out.end(), 0.0);
69
70 const int qOrder1 = 4;
72
73 for (auto&& qp : rule1)
74 {
75 Scalar qPos = qp.position();
76 typename LB::Traits::DomainType localPos;
77
78 localPos = {qPos, 0.0};
79 auto y = f(localPos);
80 out[0] += (y[0]*n_[0][0] + y[1]*n_[0][1])*qp.weight()*sign_[0]/c_[0];
81 out[3] += (y[0]*n_[0][0] + y[1]*n_[0][1])*(2.0*qPos - 1.0)*qp.weight()/c_[0];
82
83 localPos = {0.0, qPos};
84 y = f(localPos);
85 out[1] += (y[0]*n_[1][0] + y[1]*n_[1][1])*qp.weight()*sign_[1]/c_[1];
86 out[4] += (y[0]*n_[1][0] + y[1]*n_[1][1])*(1.0 - 2.0*qPos)*qp.weight()/c_[1];
87
88 localPos = {1.0 - qPos, qPos};
89 y = f(localPos);
90 out[2] += (y[0]*n_[2][0] + y[1]*n_[2][1])*qp.weight()*sign_[2]/c_[2];
91 out[5] += (y[0]*n_[2][0] + y[1]*n_[2][1])*(2.0*qPos - 1.0)*qp.weight()/c_[2];
92 }
93
94 const int qOrder2 = 8;
96
97 for (auto&& qp : rule2)
98 {
99 auto qPos = qp.position();
100
101 auto y = f(qPos);
102 out[6] += y[0]*qp.weight();
103 out[7] += y[1]*qp.weight();
104 }
105 }
106
107 private:
108 // Edge orientations
109 std::array<typename LB::Traits::RangeFieldType, 3> sign_;
110
111 // Edge normals
112 std::array<typename LB::Traits::DomainType, 3> n_;
113
114 std::array<typename LB::Traits::RangeFieldType, 3> c_;
115 };
116}
117#endif // DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS12DLOCALINTERPOLATION_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
First order Raviart-Thomas shape functions on the reference quadrilateral.
Definition: raviartthomas12dlocalinterpolation.hh:26
RT12DLocalInterpolation(std::bitset< 3 > s=0)
Make set number s, where 0 <= s < 8.
Definition: raviartthomas12dlocalinterpolation.hh:35
void interpolate(const F &ff, std::vector< C > &out) const
Interpolate a given function with shape functions.
Definition: raviartthomas12dlocalinterpolation.hh:59
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.111.3 (Dec 21, 23:30, 2024)