Dune Core Modules (2.7.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#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:
26
32 RT1Cube3DLocalInterpolation (unsigned int s = 0)
33 {
34 sign0 = sign1 = sign2 = sign3 = sign4 = sign5 = 1.0;
35 if (s & 1)
36 {
37 sign0 = -1.0;
38 }
39 if (s & 2)
40 {
41 sign1 = -1.0;
42 }
43 if (s & 4)
44 {
45 sign2 = -1.0;
46 }
47 if (s & 8)
48 {
49 sign3 = -1.0;
50 }
51 if (s & 16)
52 {
53 sign4 = -1.0;
54 }
55 if (s & 32)
56 {
57 sign5 = -1.0;
58 }
59
60 n0[0] = -1.0;
61 n0[1] = 0.0;
62 n0[2] = 0.0;
63 n1[0] = 1.0;
64 n1[1] = 0.0;
65 n1[2] = 0.0;
66 n2[0] = 0.0;
67 n2[1] = -1.0;
68 n2[2] = 0.0;
69 n3[0] = 0.0;
70 n3[1] = 1.0;
71 n3[2] = 0.0;
72 n4[0] = 0.0;
73 n4[1] = 0.0;
74 n4[2] = -1.0;
75 n5[0] = 0.0;
76 n5[1] = 0.0;
77 n5[2] = 1.0;
78 }
79
88 template<class F, class C>
89 void interpolate (const F& ff, std::vector<C>& out) const
90 {
91 // f gives v*outer normal at a point on the edge!
92 typedef typename LB::Traits::RangeFieldType Scalar;
93 typedef typename LB::Traits::DomainFieldType Vector;
94
95 auto&& f = Impl::makeFunctionWithCallOperator<typename LB::Traits::DomainType>(ff);
96
97 out.resize(36);
98 fill(out.begin(), out.end(), 0.0);
99
100 const int qOrder = 3;
102
103 for (typename QuadratureRule<Scalar,2>::const_iterator it = rule1.begin();
104 it != rule1.end(); ++it)
105 {
106 Dune::FieldVector<Scalar,2> qPos = it->position();
107 typename LB::Traits::DomainType localPos;
108
109 localPos[0] = 0.0;
110 localPos[1] = qPos[0];
111 localPos[2] = qPos[1];
112 auto y = f(localPos);
113 out[0] += (y[0]*n0[0] + y[1]*n0[1] + y[2]*n0[2])*it->weight()*sign0;
114 out[6] += (y[0]*n0[0] + y[1]*n0[1] + y[2]*n0[2])*(2.0*qPos[0] - 1.0)*it->weight();
115 out[12] += (y[0]*n0[0] + y[1]*n0[1] + y[2]*n0[2])*(2.0*qPos[1] - 1.0)*it->weight();
116 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();
117
118 localPos[0] = 1.0;
119 localPos[1] = qPos[0];
120 localPos[2] = qPos[1];
121 y = f(localPos);
122 out[1] += (y[0]*n1[0] + y[1]*n1[1] + y[2]*n1[2])*it->weight()*sign1;
123 out[7] += (y[0]*n1[0] + y[1]*n1[1] + y[2]*n1[2])*(1.0 - 2.0*qPos[0])*it->weight();
124 out[13] += (y[0]*n1[0] + y[1]*n1[1] + y[2]*n1[2])*(1.0 - 2.0*qPos[1])*it->weight();
125 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();
126
127 localPos[0] = qPos[0];
128 localPos[1] = 0.0;
129 localPos[2] = qPos[1];
130 y = f(localPos);
131 out[2] += (y[0]*n2[0] + y[1]*n2[1] + y[2]*n2[2])*it->weight()*sign2;
132 out[8] += (y[0]*n2[0] + y[1]*n2[1] + y[2]*n2[2])*(1.0 - 2.0*qPos[0])*it->weight();
133 out[14] += (y[0]*n2[0] + y[1]*n2[1] + y[2]*n2[2])*(2.0*qPos[1] - 1.0)*it->weight();
134 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();
135
136 localPos[0] = qPos[0];
137 localPos[1] = 1.0;
138 localPos[2] = qPos[1];
139 y = f(localPos);
140 out[3] += (y[0]*n3[0] + y[1]*n3[1] + y[2]*n3[2])*it->weight()*sign3;
141 out[9] += (y[0]*n3[0] + y[1]*n3[1] + y[2]*n3[2])*(2.0*qPos[0] - 1.0)*it->weight();
142 out[15] += (y[0]*n3[0] + y[1]*n3[1] + y[2]*n3[2])*(1.0 - 2.0*qPos[1])*it->weight();
143 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();
144
145 localPos[0] = qPos[0];
146 localPos[1] = qPos[1];
147 localPos[2] = 0.0;
148 y = f(localPos);
149 out[4] += (y[0]*n4[0] + y[1]*n4[1] + y[2]*n4[2])*it->weight()*sign4;
150 out[10] += (y[0]*n4[0] + y[1]*n4[1] + y[2]*n4[2])*(1.0 - 2.0*qPos[0])*it->weight();
151 out[16] += (y[0]*n4[0] + y[1]*n4[1] + y[2]*n4[2])*(1.0 - 2.0*qPos[1])*it->weight();
152 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();
153
154 localPos[0] = qPos[0];
155 localPos[1] = qPos[1];
156 localPos[2] = 1.0;
157 y = f(localPos);
158 out[5] += (y[0]*n5[0] + y[1]*n5[1] + y[2]*n5[2])*it->weight()*sign5;
159 out[11] += (y[0]*n5[0] + y[1]*n5[1] + y[2]*n5[2])*(2.0*qPos[0] - 1.0)*it->weight();
160 out[17] += (y[0]*n5[0] + y[1]*n5[1] + y[2]*n5[2])*(2.0*qPos[1] - 1.0)*it->weight();
161 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();
162 }
163
165 for (typename QuadratureRule<Vector,3>::const_iterator it = rule2.begin();
166 it != rule2.end(); ++it)
167 {
168 FieldVector<double,3> qPos = it->position();
169
170 auto y = f(qPos);
171 out[24] += y[0]*it->weight();
172 out[25] += y[1]*it->weight();
173 out[26] += y[2]*it->weight();
174 out[27] += y[0]*qPos[1]*it->weight();
175 out[28] += y[0]*qPos[2]*it->weight();
176 out[29] += y[1]*qPos[0]*it->weight();
177 out[30] += y[1]*qPos[2]*it->weight();
178 out[31] += y[2]*qPos[0]*it->weight();
179 out[32] += y[2]*qPos[1]*it->weight();
180 out[33] += y[0]*qPos[1]*qPos[2]*it->weight();
181 out[34] += y[1]*qPos[0]*qPos[2]*it->weight();
182 out[35] += y[2]*qPos[0]*qPos[1]*it->weight();
183 }
184 }
185
186 private:
187 typename LB::Traits::RangeFieldType sign0, sign1, sign2, sign3, sign4, sign5;
188 typename LB::Traits::DomainType n0, n1, n2, n3, n4, n5;
189 };
190}
191#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=0)
Make set number s, where 0 <= s < 64.
Definition: raviartthomas1cube3dlocalinterpolation.hh:32
void interpolate(const F &ff, std::vector< C > &out) const
Interpolate a given function with shape functions.
Definition: raviartthomas1cube3dlocalinterpolation.hh:89
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 (Nov 13, 23:29, 2024)