Dune Core Modules (2.9.1)

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