Dune Core Modules (2.10.0)

refinedsimplexlocalbasis.hh
Go to the documentation of this file.
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_REFINED_SIMPLEX_LOCALBASIS_HH
6#define DUNE_REFINED_SIMPLEX_LOCALBASIS_HH
7
14#include <dune/localfunctions/common/localbasis.hh>
15
16namespace Dune
17{
18 template<class D, int dim>
19 class RefinedSimplexLocalBasis
20 {
21 protected:
22 RefinedSimplexLocalBasis()
23 {
24 DUNE_THROW(Dune::NotImplemented,"RefinedSimplexLocalBasis not implemented for dim > 3.");
25 }
26 };
27
36 template<class D>
37 class RefinedSimplexLocalBasis<D,1>
38 {
39 protected:
40
43
54 static int getSubElement(const FieldVector<D,1>& global)
55 {
56 if (global[0] <= 0.5)
57 return 0;
58 else if (global[0] <= 1.0)
59 return 1;
60
61 DUNE_THROW(InvalidStateException, "no subelement defined");
62 }
63
70 static void getSubElement(const FieldVector<D,1>& global,
71 int& subElement,
72 FieldVector<D,1>& local)
73 {
74 if (global[0] <= 0.5) {
75 subElement = 0;
76 local[0] = 2.0 * global[0];
77 return;
78 }
79
80 subElement = 1;
81 local[0] = 2.0 * global[0] - 1.0;
82 }
83
84 };
85
86
98 template<class D>
99 class RefinedSimplexLocalBasis<D,2>
100 {
101 protected:
102
105
121 static int getSubElement(const FieldVector<D,2>& global)
122 {
123 if (global[0] + global[1] <= 0.5)
124 return 0;
125 else if (global[0] >= 0.5)
126 return 1;
127 else if (global[1] >= 0.5)
128 return 2;
129
130 return 3;
131 }
132
139 static void getSubElement(const FieldVector<D,2>& global,
140 int& subElement,
141 FieldVector<D,2>& local)
142 {
143 if (global[0] + global[1] <= 0.5) {
144 subElement = 0;
145 local[0] = 2*global[0];
146 local[1] = 2*global[1];
147 return;
148 } else if (global[0] >= 0.5) {
149 subElement = 1;
150 local[0] = 2*global[0]-1;
151 local[1] = 2*global[1];
152 return;
153 } else if (global[1] >= 0.5) {
154 subElement = 2;
155 local[0] = 2*global[0];
156 local[1] = 2*global[1]-1;
157 return;
158 }
159
160 subElement = 3;
161 local[0] = -2 * global[0] + 1;
162 local[1] = -2 * global[1] + 1;
163
164 }
165
166
167 };
168
180 template<class D>
181 class RefinedSimplexLocalBasis<D,3>
182 {
183 protected:
184
187
218 static int getSubElement(const FieldVector<D,3>& global)
219 {
220 if (global[0] + global[1] + global[2] <= 0.5)
221 return 0;
222 else if (global[0] >= 0.5)
223 return 1;
224 else if (global[1] >= 0.5)
225 return 2;
226 else if (global[2] >= 0.5)
227 return 3;
228 else if ((global[0] + global[1] <= 0.5)and (global[1] + global[2] <= 0.5))
229 return 4;
230 else if ((global[0] + global[1] >= 0.5)and (global[1] + global[2] <= 0.5))
231 return 5;
232 else if ((global[0] + global[1] <= 0.5)and (global[1] + global[2] >= 0.5))
233 return 6;
234 else if ((global[0] + global[1] >= 0.5)and (global[1] + global[2] >= 0.5))
235 return 7;
236
237 DUNE_THROW(InvalidStateException, "no subelement defined");
238
239 }
246 static void getSubElement(const FieldVector<D,3>& global,
247 int& subElement,
248 FieldVector<D,3>& local)
249 {
250 if (global[0] + global[1] + global[2] <= 0.5) {
251 subElement = 0;
252 local = global;
253 local *= 2.0;
254 return;
255 } else if (global[0] >= 0.5) {
256 subElement = 1;
257 local = global;
258 local[0] -= 0.5;
259 local *= 2.0;
260 return;
261 } else if (global[1] >= 0.5) {
262 subElement = 2;
263 local = global;
264 local[1] -= 0.5;
265 local *= 2.0;
266 return;
267 } else if (global[2] >= 0.5) {
268 subElement = 3;
269 local = global;
270 local[2] -= 0.5;
271 local *= 2.0;
272 return;
273 } else if ((global[0] + global[1] <= 0.5)and (global[1] + global[2] <= 0.5)) {
274 subElement = 4;
275 local[0] = 2.0 * global[1];
276 local[1] = 2.0 * (0.5 - global[0] - global[1]);
277 local[2] = 2.0 * (-0.5 + global[0] + global[1] + global[2]);
278 // Dune::FieldMatrix<double,3,3> A(0.0);
279 // A[0][1] = 2.0;
280 // A[1][0] = -2.0;
281 // A[1][1] = -2.0;
282 // A[2][0] = 2.0;
283 // A[2][1] = 2.0;
284 // A[2][2] = 2.0;
285 // A.mv(global,local);
286 // local[1] += 1.0;
287 // local[2] -= 1.0;
288 return;
289 } else if ((global[0] + global[1] >= 0.5)and (global[1] + global[2] <= 0.5)) {
290 subElement = 5;
291 local[0] = 2.0 * (0.5 - global[0]);
292 local[1] = 2.0 * (0.5 - global[1] - global[2]);
293 local[2] = 2.0 * global[2];
294 // Dune::FieldMatrix<double,3,3> A(0.0);
295 // A[0][0] = -2.0;
296 // A[1][1] = -2.0;
297 // A[1][2] = -2.0;
298 // A[2][2] = 2.0;
299 // A.mv(global,local);
300 // local[0] += 1.0;
301 // local[1] += 1.0;
302 return;
303 } else if ((global[0] + global[1] <= 0.5)and (global[1] + global[2] >= 0.5)) {
304 subElement = 6;
305 local[0] = 2.0 * (0.5 - global[0] - global[1]);
306 local[1] = 2.0 * global[0];
307 local[2] = 2.0 * (-0.5 + global[1] + global[2]);
308 // Dune::FieldMatrix<double,3,3> A(0.0);
309 // A[0][0] = -2.0;
310 // A[0][1] = -2.0;
311 // A[1][0] = 2.0;
312 // A[2][1] = 2.0;
313 // A[2][2] = 2.0;
314 // A.mv(global,local);
315 // local[0] += 1.0;
316 // local[2] -= 1.0;
317 return;
318 } else if ((global[0] + global[1] >= 0.5)and (global[1] + global[2] >= 0.5)) {
319 subElement = 7;
320 local[0] = 2.0 * (-0.5 + global[1] + global[2]);
321 local[1] = 2.0 * (0.5 - global[1]);
322 local[2] = 2.0 * (-0.5 + global[0] + global[1]);
323 // Dune::FieldMatrix<double,3,3> A(0.0);
324 // A[0][1] = 2.0;
325 // A[0][2] = 2.0;
326 // A[1][1] = -2.0;
327 // A[2][0] = 2.0;
328 // A[2][1] = 2.0;
329 // A.mv(global,local);
330 // local[0] -= 1.0;
331 // local[1] += 1.0;
332 // local[2] -= 1.0;
333 return;
334 }
335
336 DUNE_THROW(InvalidStateException, "no subelement defined");
337
338 }
339
340 };
341
342
343}
344
345#endif
vector space out of a tensor product of fields.
Definition: fvector.hh:91
Default exception if a function was called while the object is not in a valid state for that function...
Definition: exceptions.hh:281
Default exception for dummy implementations.
Definition: exceptions.hh:263
RefinedSimplexLocalBasis()
Protected default constructor so this class can only be instantiated as a base class.
Definition: refinedsimplexlocalbasis.hh:42
static void getSubElement(const FieldVector< D, 1 > &global, int &subElement, FieldVector< D, 1 > &local)
Get local coordinates in the subelement.
Definition: refinedsimplexlocalbasis.hh:70
static int getSubElement(const FieldVector< D, 1 > &global)
Get the number of the subelement containing a given point.
Definition: refinedsimplexlocalbasis.hh:54
RefinedSimplexLocalBasis()
Protected default constructor so this class can only be instantiated as a base class.
Definition: refinedsimplexlocalbasis.hh:104
static int getSubElement(const FieldVector< D, 2 > &global)
Get the number of the subtriangle containing a given point.
Definition: refinedsimplexlocalbasis.hh:121
static void getSubElement(const FieldVector< D, 2 > &global, int &subElement, FieldVector< D, 2 > &local)
Get local coordinates in the subtriangle.
Definition: refinedsimplexlocalbasis.hh:139
static void getSubElement(const FieldVector< D, 3 > &global, int &subElement, FieldVector< D, 3 > &local)
Get local coordinates in the subsimplex.
Definition: refinedsimplexlocalbasis.hh:246
RefinedSimplexLocalBasis()
Protected default constructor so this class can only be instantiated as a base class.
Definition: refinedsimplexlocalbasis.hh:186
static int getSubElement(const FieldVector< D, 3 > &global)
Get the number of the subsimplex containing a given point in the reference element.
Definition: refinedsimplexlocalbasis.hh:218
A few common exception classes.
Implements a vector constructed from a given type representing a field and a compile-time given size.
#define DUNE_THROW(E, m)
Definition: exceptions.hh:218
Dune namespace.
Definition: alignedallocator.hh:13
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Dec 27, 23:30, 2024)