DUNE-FEM (unstable)

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