DUNE-FEM (unstable)

brezzidouglasmarini2simplex2dlocalbasis.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_BREZZIDOUGLASMARINI2_SIMPLEX2D_LOCALBASIS_HH
6#define DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI2_SIMPLEX2D_LOCALBASIS_HH
7
8#include <array>
9#include <bitset>
10#include <numeric>
11#include <vector>
12
14
15#include "../../common/localbasis.hh"
16
17namespace Dune
18{
28 template<class D, class R>
30 {
31
32 public:
36
39 {
40 for (size_t i=0; i<3; i++)
41 sign_[i] = 1.0;
42 }
43
49 BDM2Simplex2DLocalBasis(std::bitset<3> s)
50 {
51 for (size_t i=0; i<3; i++)
52 sign_[i] = s[i] ? -1.0 : 1.0;
53 }
54
56 unsigned int size() const
57 {
58 return 12;
59 }
60
67 inline void evaluateFunction(const typename Traits::DomainType& in,
68 std::vector<typename Traits::RangeType>& out) const
69 {
70 out.resize(size());
71
72 out[0][0] = sign_[0]*(-2*in[0]*in[1] + in[0]*in[0]);
73 out[0][1] = sign_[0]*(-1 + 6*in[1] -2*in[0]*in[1] - 5*in[1]*in[1]);
74
75 out[1][0] = 1.5*in[0] + 3*in[0]*in[1] - 4.5*in[0]*in[0];
76 out[1][1] = -3 + 6*in[0] + 10.5*in[1] - 15*in[0]*in[1] - 7.5*in[1]*in[1];
77
78 out[2][0] = sign_[0]*(-7.5*in[0] + 5*in[0]*in[1] + 12.5*in[0]*in[0]);
79 out[2][1] = sign_[0]*(-5 + 30*in[0] + 7.5*in[1] - 25*in[0]*in[1] - 30*in[0]*in[0] - 2.5*in[1]*in[1]);
80
81
82
83 out[3][0] = sign_[1]*(-1 + 6*in[0] - 2*in[0]*in[1] - 5*in[0]*in[0]);
84 out[3][1] = sign_[1]*(-2*in[0]*in[1] + in[1]*in[1]);
85
86 out[4][0] = 3 - 10.5*in[0] - 6*in[1] + 15*in[0]*in[1] + 7.5*in[0]*in[0];
87 out[4][1] = -1.5*in[1] - 3*in[0]*in[1] + 4.5*in[1]*in[1];
88
89 out[5][0] = sign_[1]*(-5 + 7.5*in[0] + 30*in[1] - 25*in[0]*in[1] - 2.5*in[0]*in[0] - 30*in[1]*in[1]);
90 out[5][1] = sign_[1]*(-7.5*in[1] + 5*in[0]*in[1] + 12.5*in[1]*in[1]);
91
92
93
94 out[6][0] = sign_[2]*(-3*in[0] + 4*in[0]*in[1] + 4*in[0]*in[0]);
95 out[6][1] = sign_[2]*(-3*in[1] + 4*in[0]*in[1] + 4*in[1]*in[1]);
96
97 out[7][0] = -3*in[0] + 6*in[0]*in[0];
98 out[7][1] = 3*in[1] - 6*in[1]*in[1];
99
100 out[8][0] = sign_[2]*(-10*in[0]*in[1] + 5*in[0]*in[0]);
101 out[8][1] = sign_[2]*(-10*in[0]*in[1] + 5*in[1]*in[1]);
102
103
104
105 out[9][0] = 18*in[0] - 12*in[0]*in[1] - 18*in[0]*in[0];
106 out[9][1] = 6*in[1] - 12*in[0]*in[1] - 6*in[1]*in[1];
107
108 out[10][0] = 6*in[0] - 12*in[0]*in[1] - 6*in[0]*in[0];
109 out[10][1] = 18*in[1] - 12*in[0]*in[1] - 18*in[1]*in[1];
110
111 out[11][0] = 90*in[0] - 180*in[0]*in[1] - 90*in[0]*in[0];
112 out[11][1] = -90*in[1] + 180*in[0]*in[1] + 90*in[1]*in[1];
113 }
114
121 inline void evaluateJacobian(const typename Traits::DomainType& in,
122 std::vector<typename Traits::JacobianType>& out) const
123 {
124 out.resize(size());
125
126 out[0][0][0] = sign_[0]*(-2*in[1] + 2*in[0]);
127 out[0][0][1] = sign_[0]*(-2*in[0]);
128
129 out[0][1][0] = sign_[0]*(-2*in[1]);
130 out[0][1][1] = sign_[0]*(6 -2*in[0] - 10*in[1]);
131
132
133 out[1][0][0] = 1.5 + 3*in[1] - 9*in[0];
134 out[1][0][1] = 3*in[0];
135
136 out[1][1][0] = 6 - 15*in[1];
137 out[1][1][1] = 10.5 - 15*in[0] - 15*in[1];
138
139
140 out[2][0][0] = sign_[0]*(-7.5 + 5*in[1] + 25*in[0]);
141 out[2][0][1] = sign_[0]*(5*in[0]);
142
143 out[2][1][0] = sign_[0]*(30 - 25*in[1] - 60*in[0]);
144 out[2][1][1] = sign_[0]*(7.5 - 25*in[0] - 5*in[1]);
145
146
147
148 out[3][0][0] = sign_[1]*(6 - 2*in[1] - 10*in[0]);
149 out[3][0][1] = sign_[1]*(-2*in[0]);
150
151 out[3][1][0] = sign_[1]*(-2*in[1]);
152 out[3][1][1] = sign_[1]*(-2*in[0] + 2*in[1]);
153
154
155 out[4][0][0] = -10.5 + 15*in[1] + 15*in[0];
156 out[4][0][1] = -6 + 15*in[0];
157
158 out[4][1][0] = -3*in[1];
159 out[4][1][1] = -1.5 - 3*in[0] + 9*in[1];
160
161
162 out[5][0][0] = sign_[1]*(7.5 - 25*in[1] - 5*in[0]);
163 out[5][0][1] = sign_[1]*(30 - 25*in[0] - 60*in[1]);
164
165 out[5][1][0] = sign_[1]*(5*in[1]);
166 out[5][1][1] = sign_[1]*(-7.5 + 5*in[0] + 25*in[1]);
167
168
169
170 out[6][0][0] = sign_[2]*(-3 + 4*in[1] + 8*in[0]);
171 out[6][0][1] = sign_[2]*(4*in[0]);
172
173 out[6][1][0] = sign_[2]*(4*in[1]);
174 out[6][1][1] = sign_[2]*(-3 + 4*in[0] + 8*in[1]);
175
176
177 out[7][0][0] = -3 + 12*in[0];
178 out[7][0][1] = 0;
179
180 out[7][1][0] = 0;
181 out[7][1][1] = 3 - 12*in[1];
182
183
184 out[8][0][0] = sign_[2]*(-10*in[1] + 10*in[0]);
185 out[8][0][1] = sign_[2]*(-10*in[0]);
186
187 out[8][1][0] = sign_[2]*(-10*in[1]);
188 out[8][1][1] = sign_[2]*(-10*in[0] + 10*in[1]);
189
190
191 out[9][0][0] = 18 - 12*in[1] - 36*in[0];
192 out[9][0][1] = -12*in[0];
193
194 out[9][1][0] = -12*in[1];
195 out[9][1][1] = 6 - 12*in[0] - 12*in[1];
196
197 out[10][0][0] = 6 - 12*in[1] - 12*in[0];
198 out[10][0][1] = -12*in[0];
199
200 out[10][1][0] = -12*in[1];
201 out[10][1][1] = 18 - 12*in[0] - 36*in[1];
202
203 out[11][0][0] = 90 - 180*in[1] - 180*in[0];
204 out[11][0][1] = -180*in[0];
205
206 out[11][1][0] = 180*in[1];
207 out[11][1][1] = -90 + 180*in[0] + 180*in[1];
208 }
209
211 void partial (const std::array<unsigned int, 2>& order,
212 const typename Traits::DomainType& in, // position
213 std::vector<typename Traits::RangeType>& out) const // return value
214 {
215 auto totalOrder = std::accumulate(order.begin(), order.end(), 0);
216 if (totalOrder == 0) {
217 evaluateFunction(in, out);
218 } else if (totalOrder == 1) {
219 out.resize(size());
220 auto const direction = std::distance(order.begin(), std::find(order.begin(), order.end(), 1));
221
222 switch (direction) {
223 case 0:
224 out[0][0] = sign_[0]*(-2*in[1] + 2*in[0]);
225 out[0][1] = sign_[0]*(-2*in[1]);
226
227 out[1][0] = 1.5 + 3*in[1] - 9*in[0];
228 out[1][1] = 6 - 15*in[1];
229
230 out[2][0] = sign_[0]*(-7.5 + 5*in[1] + 25*in[0]);
231 out[2][1] = sign_[0]*(30 - 25*in[1] - 60*in[0]);
232
233 out[3][0] = sign_[1]*(6 - 2*in[1] - 10*in[0]);
234 out[3][1] = sign_[1]*(-2*in[1]);
235
236 out[4][0] = -10.5 + 15*in[1] + 15*in[0];
237 out[4][1] = -3*in[1];
238
239 out[5][0] = sign_[1]*(7.5 - 25*in[1] - 5*in[0]);
240 out[5][1] = sign_[1]*(5*in[1]);
241
242 out[6][0] = sign_[2]*(-3 + 4*in[1] + 8*in[0]);
243 out[6][1] = sign_[2]*(4*in[1]);
244
245 out[7][0] = -3 + 12*in[0];
246 out[7][1] = 0;
247
248 out[8][0] = sign_[2]*(-10*in[1] + 10*in[0]);
249 out[8][1] = sign_[2]*(-10*in[1]);
250
251 out[9][0] = 18 - 12*in[1] - 36*in[0];
252 out[9][1] = -12*in[1];
253
254 out[10][0] = 6 - 12*in[1] - 12*in[0];
255 out[10][1] = -12*in[1];
256
257 out[11][0] = 90 - 180*in[1] - 180*in[0];
258 out[11][1] = 180*in[1];
259 break;
260 case 1:
261 out[0][0] = sign_[0]*(-2*in[0]);
262 out[0][1] = sign_[0]*(6 -2*in[0] - 10*in[1]);
263
264 out[1][0] = 3*in[0];
265 out[1][1] = 10.5 - 15*in[0] - 15*in[1];
266
267 out[2][0] = sign_[0]*(5*in[0]);
268 out[2][1] = sign_[0]*(7.5 - 25*in[0] - 5*in[1]);
269
270 out[3][0] = sign_[1]*(-2*in[0]);
271 out[3][1] = sign_[1]*(-2*in[0] + 2*in[1]);
272
273 out[4][0] = -6 + 15*in[0];
274 out[4][1] = -1.5 - 3*in[0] + 9*in[1];
275
276 out[5][0] = sign_[1]*(30 - 25*in[0] - 60*in[1]);
277 out[5][1] = sign_[1]*(-7.5 + 5*in[0] + 25*in[1]);
278
279 out[6][0] = sign_[2]*(4*in[0]);
280 out[6][1] = sign_[2]*(-3 + 4*in[0] + 8*in[1]);
281
282 out[7][0] = 0;
283 out[7][1] = 3 - 12*in[1];
284
285 out[8][0] = sign_[2]*(-10*in[0]);
286 out[8][1] = sign_[2]*(-10*in[0] + 10*in[1]);
287
288 out[9][0] = -12*in[0];
289 out[9][1] = 6 - 12*in[0] - 12*in[1];
290
291 out[10][0] = -12*in[0];
292 out[10][1] = 18 - 12*in[0] - 36*in[1];
293
294 out[11][0] = -180*in[0];
295 out[11][1] = -90 + 180*in[0] + 180*in[1];
296 break;
297 default:
298 DUNE_THROW(RangeError, "Component out of range.");
299 }
300 } else {
301 DUNE_THROW(NotImplemented, "Desired derivative order is not implemented");
302 }
303 }
304
306 unsigned int order() const
307 {
308 return 2; // TODO: check whether this is not order 3
309 }
310
311 private:
312 std::array<R,3> sign_;
313 };
314} // end namespace Dune
315#endif // DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI2_SIMPLEX2D_LOCALBASIS_HH
First order Brezzi-Douglas-Marini shape functions on quadrilaterals.
Definition: brezzidouglasmarini2simplex2dlocalbasis.hh:30
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: brezzidouglasmarini2simplex2dlocalbasis.hh:121
BDM2Simplex2DLocalBasis(std::bitset< 3 > s)
Make set number s, where 0 <= s < 8.
Definition: brezzidouglasmarini2simplex2dlocalbasis.hh:49
void partial(const std::array< unsigned int, 2 > &order, const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate partial derivatives of all shape functions.
Definition: brezzidouglasmarini2simplex2dlocalbasis.hh:211
unsigned int size() const
number of shape functions
Definition: brezzidouglasmarini2simplex2dlocalbasis.hh:56
BDM2Simplex2DLocalBasis()
Standard constructor.
Definition: brezzidouglasmarini2simplex2dlocalbasis.hh:38
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: brezzidouglasmarini2simplex2dlocalbasis.hh:67
unsigned int order() const
Polynomial order of the shape functions.
Definition: brezzidouglasmarini2simplex2dlocalbasis.hh:306
A dense n x m matrix.
Definition: fmatrix.hh:117
vector space out of a tensor product of fields.
Definition: fvector.hh:91
Default exception for dummy implementations.
Definition: exceptions.hh:263
Default exception class for range errors.
Definition: exceptions.hh:254
Implements a matrix constructed from a given type representing a field and compile-time given number ...
#define DUNE_THROW(E, m)
Definition: exceptions.hh:218
constexpr T accumulate(Range &&range, T value, F &&f)
Accumulate values.
Definition: hybridutilities.hh:279
Dune namespace.
Definition: alignedallocator.hh:13
Type traits for LocalBasisVirtualInterface.
Definition: localbasis.hh:35
D DomainType
domain type
Definition: localbasis.hh:43
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 12, 23:30, 2024)