Dune Core Modules (2.7.0)

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#ifndef DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS1_CUBE3D_LOCALINTERPOLATION_HH
4#define DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS1_CUBE3D_LOCALINTERPOLATION_HH
5
6#include <vector>
7
9#include <dune/localfunctions/common/localinterpolation.hh>
10
11namespace Dune
12{
21 template<class LB>
23 {
24
25 public:
28 {
29 sign0 = sign1 = sign2 = sign3 = sign4 = sign5 = 1.0;
30 }
31
38 {
39 sign0 = sign1 = sign2 = sign3 = sign4 = sign5 = 1.0;
40 if (s & 1)
41 {
42 sign0 = -1.0;
43 }
44 if (s & 2)
45 {
46 sign1 = -1.0;
47 }
48 if (s & 4)
49 {
50 sign2 = -1.0;
51 }
52 if (s & 8)
53 {
54 sign3 = -1.0;
55 }
56 if (s & 16)
57 {
58 sign4 = -1.0;
59 }
60 if (s & 32)
61 {
62 sign5 = -1.0;
63 }
64
65 n0[0] = -1.0;
66 n0[1] = 0.0;
67 n0[2] = 0.0;
68 n1[0] = 1.0;
69 n1[1] = 0.0;
70 n1[2] = 0.0;
71 n2[0] = 0.0;
72 n2[1] = -1.0;
73 n2[2] = 0.0;
74 n3[0] = 0.0;
75 n3[1] = 1.0;
76 n3[2] = 0.0;
77 n4[0] = 0.0;
78 n4[1] = 0.0;
79 n4[2] = -1.0;
80 n5[0] = 0.0;
81 n5[1] = 0.0;
82 n5[2] = 1.0;
83 }
84
93 template<class F, class C>
94 void interpolate (const F& ff, std::vector<C>& out) const
95 {
96 // f gives v*outer normal at a point on the edge!
97 typedef typename LB::Traits::RangeFieldType Scalar;
98 typedef typename LB::Traits::DomainFieldType Vector;
99
100 auto&& f = Impl::makeFunctionWithCallOperator<typename LB::Traits::DomainType>(ff);
101
102 out.resize(36);
103 fill(out.begin(), out.end(), 0.0);
104
105 const int qOrder = 3;
107
108 for (typename QuadratureRule<Scalar,2>::const_iterator it = rule1.begin();
109 it != rule1.end(); ++it)
110 {
111 Dune::FieldVector<Scalar,2> qPos = it->position();
112 typename LB::Traits::DomainType localPos;
113
114 localPos[0] = 0.0;
115 localPos[1] = qPos[0];
116 localPos[2] = qPos[1];
117 auto y = f(localPos);
118 out[0] += (y[0]*n0[0] + y[1]*n0[1] + y[2]*n0[2])*it->weight()*sign0;
119 out[6] += (y[0]*n0[0] + y[1]*n0[1] + y[2]*n0[2])*(2.0*qPos[0] - 1.0)*it->weight();
120 out[12] += (y[0]*n0[0] + y[1]*n0[1] + y[2]*n0[2])*(2.0*qPos[1] - 1.0)*it->weight();
121 out[18] += (y[0]*n0[0] + y[1]*n0[1] + y[2]*n0[2])*(2.0*qPos[0] - 1.0)*(2.0*qPos[1] - 1.0)*it->weight();
122
123 localPos[0] = 1.0;
124 localPos[1] = qPos[0];
125 localPos[2] = qPos[1];
126 y = f(localPos);
127 out[1] += (y[0]*n1[0] + y[1]*n1[1] + y[2]*n1[2])*it->weight()*sign1;
128 out[7] += (y[0]*n1[0] + y[1]*n1[1] + y[2]*n1[2])*(1.0 - 2.0*qPos[0])*it->weight();
129 out[13] += (y[0]*n1[0] + y[1]*n1[1] + y[2]*n1[2])*(1.0 - 2.0*qPos[1])*it->weight();
130 out[19] += (y[0]*n1[0] + y[1]*n1[1] + y[2]*n1[2])*(1.0 - 2.0*qPos[0])*(2.0*qPos[1] - 1.0)*it->weight();
131
132 localPos[0] = qPos[0];
133 localPos[1] = 0.0;
134 localPos[2] = qPos[1];
135 y = f(localPos);
136 out[2] += (y[0]*n2[0] + y[1]*n2[1] + y[2]*n2[2])*it->weight()*sign2;
137 out[8] += (y[0]*n2[0] + y[1]*n2[1] + y[2]*n2[2])*(1.0 - 2.0*qPos[0])*it->weight();
138 out[14] += (y[0]*n2[0] + y[1]*n2[1] + y[2]*n2[2])*(2.0*qPos[1] - 1.0)*it->weight();
139 out[20] += (y[0]*n2[0] + y[1]*n2[1] + y[2]*n2[2])*(1.0 - 2.0*qPos[0])*(2.0*qPos[1] - 1.0)*it->weight();
140
141 localPos[0] = qPos[0];
142 localPos[1] = 1.0;
143 localPos[2] = qPos[1];
144 y = f(localPos);
145 out[3] += (y[0]*n3[0] + y[1]*n3[1] + y[2]*n3[2])*it->weight()*sign3;
146 out[9] += (y[0]*n3[0] + y[1]*n3[1] + y[2]*n3[2])*(2.0*qPos[0] - 1.0)*it->weight();
147 out[15] += (y[0]*n3[0] + y[1]*n3[1] + y[2]*n3[2])*(1.0 - 2.0*qPos[1])*it->weight();
148 out[21] += (y[0]*n3[0] + y[1]*n3[1] + y[2]*n3[2])*(2.0*qPos[0] - 1.0)*(2.0*qPos[1] - 1.0)*it->weight();
149
150 localPos[0] = qPos[0];
151 localPos[1] = qPos[1];
152 localPos[2] = 0.0;
153 y = f(localPos);
154 out[4] += (y[0]*n4[0] + y[1]*n4[1] + y[2]*n4[2])*it->weight()*sign4;
155 out[10] += (y[0]*n4[0] + y[1]*n4[1] + y[2]*n4[2])*(1.0 - 2.0*qPos[0])*it->weight();
156 out[16] += (y[0]*n4[0] + y[1]*n4[1] + y[2]*n4[2])*(1.0 - 2.0*qPos[1])*it->weight();
157 out[22] += (y[0]*n4[0] + y[1]*n4[1] + y[2]*n4[2])*(1.0 - 2.0*qPos[0])*(2.0*qPos[1] - 1.0)*it->weight();
158
159 localPos[0] = qPos[0];
160 localPos[1] = qPos[1];
161 localPos[2] = 1.0;
162 y = f(localPos);
163 out[5] += (y[0]*n5[0] + y[1]*n5[1] + y[2]*n5[2])*it->weight()*sign5;
164 out[11] += (y[0]*n5[0] + y[1]*n5[1] + y[2]*n5[2])*(2.0*qPos[0] - 1.0)*it->weight();
165 out[17] += (y[0]*n5[0] + y[1]*n5[1] + y[2]*n5[2])*(2.0*qPos[1] - 1.0)*it->weight();
166 out[23] += (y[0]*n5[0] + y[1]*n5[1] + y[2]*n5[2])*(2.0*qPos[0] - 1.0)*(2.0*qPos[1] - 1.0)*it->weight();
167 }
168
170 for (typename QuadratureRule<Vector,3>::const_iterator it = rule2.begin();
171 it != rule2.end(); ++it)
172 {
173 FieldVector<double,3> qPos = it->position();
174
175 auto y = f(qPos);
176 out[24] += y[0]*it->weight();
177 out[25] += y[1]*it->weight();
178 out[26] += y[2]*it->weight();
179 out[27] += y[0]*qPos[1]*it->weight();
180 out[28] += y[0]*qPos[2]*it->weight();
181 out[29] += y[1]*qPos[0]*it->weight();
182 out[30] += y[1]*qPos[2]*it->weight();
183 out[31] += y[2]*qPos[0]*it->weight();
184 out[32] += y[2]*qPos[1]*it->weight();
185 out[33] += y[0]*qPos[1]*qPos[2]*it->weight();
186 out[34] += y[1]*qPos[0]*qPos[2]*it->weight();
187 out[35] += y[2]*qPos[0]*qPos[1]*it->weight();
188 }
189 }
190
191 private:
192 typename LB::Traits::RangeFieldType sign0, sign1, sign2, sign3, sign4, sign5;
193 typename LB::Traits::DomainType n0, n1, n2, n3, n4, n5;
194 };
195}
196#endif // DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS1_CUBE3D_LOCALINTERPOLATION_HH
vector space out of a tensor product of fields.
Definition: fvector.hh:96
Abstract base class for quadrature rules.
Definition: quadraturerules.hh:126
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:254
First order Raviart-Thomas shape functions on the reference hexahedron.
Definition: raviartthomas1cube3dlocalinterpolation.hh:23
RT1Cube3DLocalInterpolation(unsigned int s)
Make set number s, where 0 <= s < 64.
Definition: raviartthomas1cube3dlocalinterpolation.hh:37
RT1Cube3DLocalInterpolation()
Standard constructor.
Definition: raviartthomas1cube3dlocalinterpolation.hh:27
void interpolate(const F &ff, std::vector< C > &out) const
Interpolate a given function with shape functions.
Definition: raviartthomas1cube3dlocalinterpolation.hh:94
constexpr GeometryType cube(unsigned int dim)
Returns a GeometryType representing a hypercube of dimension dim.
Definition: type.hh:775
typename Overloads::ScalarType< std::decay_t< V > >::type Scalar
Element type of some SIMD type.
Definition: interface.hh:233
Dune namespace.
Definition: alignedallocator.hh:14
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 15, 22:36, 2024)