Dune Core Modules (2.9.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 (C) 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
35 template<class D>
36 class RefinedSimplexLocalBasis<D,1>
37 {
38 protected:
39
42
53 static int getSubElement(const FieldVector<D,1>& global)
54 {
55 if (global[0] <= 0.5)
56 return 0;
57 else if (global[0] <= 1.0)
58 return 1;
59
60 DUNE_THROW(InvalidStateException, "no subelement defined");
61 }
62
69 static void getSubElement(const FieldVector<D,1>& global,
70 int& subElement,
71 FieldVector<D,1>& local)
72 {
73 if (global[0] <= 0.5) {
74 subElement = 0;
75 local[0] = 2.0 * global[0];
76 return;
77 }
78
79 subElement = 1;
80 local[0] = 2.0 * global[0] - 1.0;
81 }
82
83 };
84
85
96 template<class D>
97 class RefinedSimplexLocalBasis<D,2>
98 {
99 protected:
100
103
119 static int getSubElement(const FieldVector<D,2>& global)
120 {
121 if (global[0] + global[1] <= 0.5)
122 return 0;
123 else if (global[0] >= 0.5)
124 return 1;
125 else if (global[1] >= 0.5)
126 return 2;
127
128 return 3;
129 }
130
137 static void getSubElement(const FieldVector<D,2>& global,
138 int& subElement,
139 FieldVector<D,2>& local)
140 {
141 if (global[0] + global[1] <= 0.5) {
142 subElement = 0;
143 local[0] = 2*global[0];
144 local[1] = 2*global[1];
145 return;
146 } else if (global[0] >= 0.5) {
147 subElement = 1;
148 local[0] = 2*global[0]-1;
149 local[1] = 2*global[1];
150 return;
151 } else if (global[1] >= 0.5) {
152 subElement = 2;
153 local[0] = 2*global[0];
154 local[1] = 2*global[1]-1;
155 return;
156 }
157
158 subElement = 3;
159 local[0] = -2 * global[0] + 1;
160 local[1] = -2 * global[1] + 1;
161
162 }
163
164
165 };
166
177 template<class D>
178 class RefinedSimplexLocalBasis<D,3>
179 {
180 protected:
181
184
215 static int getSubElement(const FieldVector<D,3>& global)
216 {
217 if (global[0] + global[1] + global[2] <= 0.5)
218 return 0;
219 else if (global[0] >= 0.5)
220 return 1;
221 else if (global[1] >= 0.5)
222 return 2;
223 else if (global[2] >= 0.5)
224 return 3;
225 else if ((global[0] + global[1] <= 0.5)and (global[1] + global[2] <= 0.5))
226 return 4;
227 else if ((global[0] + global[1] >= 0.5)and (global[1] + global[2] <= 0.5))
228 return 5;
229 else if ((global[0] + global[1] <= 0.5)and (global[1] + global[2] >= 0.5))
230 return 6;
231 else if ((global[0] + global[1] >= 0.5)and (global[1] + global[2] >= 0.5))
232 return 7;
233
234 DUNE_THROW(InvalidStateException, "no subelement defined");
235
236 }
243 static void getSubElement(const FieldVector<D,3>& global,
244 int& subElement,
245 FieldVector<D,3>& local)
246 {
247 if (global[0] + global[1] + global[2] <= 0.5) {
248 subElement = 0;
249 local = global;
250 local *= 2.0;
251 return;
252 } else if (global[0] >= 0.5) {
253 subElement = 1;
254 local = global;
255 local[0] -= 0.5;
256 local *= 2.0;
257 return;
258 } else if (global[1] >= 0.5) {
259 subElement = 2;
260 local = global;
261 local[1] -= 0.5;
262 local *= 2.0;
263 return;
264 } else if (global[2] >= 0.5) {
265 subElement = 3;
266 local = global;
267 local[2] -= 0.5;
268 local *= 2.0;
269 return;
270 } else if ((global[0] + global[1] <= 0.5)and (global[1] + global[2] <= 0.5)) {
271 subElement = 4;
272 local[0] = 2.0 * global[1];
273 local[1] = 2.0 * (0.5 - global[0] - global[1]);
274 local[2] = 2.0 * (-0.5 + global[0] + global[1] + global[2]);
275 // Dune::FieldMatrix<double,3,3> A(0.0);
276 // A[0][1] = 2.0;
277 // A[1][0] = -2.0;
278 // A[1][1] = -2.0;
279 // A[2][0] = 2.0;
280 // A[2][1] = 2.0;
281 // A[2][2] = 2.0;
282 // A.mv(global,local);
283 // local[1] += 1.0;
284 // local[2] -= 1.0;
285 return;
286 } else if ((global[0] + global[1] >= 0.5)and (global[1] + global[2] <= 0.5)) {
287 subElement = 5;
288 local[0] = 2.0 * (0.5 - global[0]);
289 local[1] = 2.0 * (0.5 - global[1] - global[2]);
290 local[2] = 2.0 * global[2];
291 // Dune::FieldMatrix<double,3,3> A(0.0);
292 // A[0][0] = -2.0;
293 // A[1][1] = -2.0;
294 // A[1][2] = -2.0;
295 // A[2][2] = 2.0;
296 // A.mv(global,local);
297 // local[0] += 1.0;
298 // local[1] += 1.0;
299 return;
300 } else if ((global[0] + global[1] <= 0.5)and (global[1] + global[2] >= 0.5)) {
301 subElement = 6;
302 local[0] = 2.0 * (0.5 - global[0] - global[1]);
303 local[1] = 2.0 * global[0];
304 local[2] = 2.0 * (-0.5 + global[1] + global[2]);
305 // Dune::FieldMatrix<double,3,3> A(0.0);
306 // A[0][0] = -2.0;
307 // A[0][1] = -2.0;
308 // A[1][0] = 2.0;
309 // A[2][1] = 2.0;
310 // A[2][2] = 2.0;
311 // A.mv(global,local);
312 // local[0] += 1.0;
313 // local[2] -= 1.0;
314 return;
315 } else if ((global[0] + global[1] >= 0.5)and (global[1] + global[2] >= 0.5)) {
316 subElement = 7;
317 local[0] = 2.0 * (-0.5 + global[1] + global[2]);
318 local[1] = 2.0 * (0.5 - global[1]);
319 local[2] = 2.0 * (-0.5 + global[0] + global[1]);
320 // Dune::FieldMatrix<double,3,3> A(0.0);
321 // A[0][1] = 2.0;
322 // A[0][2] = 2.0;
323 // A[1][1] = -2.0;
324 // A[2][0] = 2.0;
325 // A[2][1] = 2.0;
326 // A.mv(global,local);
327 // local[0] -= 1.0;
328 // local[1] += 1.0;
329 // local[2] -= 1.0;
330 return;
331 }
332
333 DUNE_THROW(InvalidStateException, "no subelement defined");
334
335 }
336
337 };
338
339
340}
341
342#endif
vector space out of a tensor product of fields.
Definition: fvector.hh:95
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:41
static void getSubElement(const FieldVector< D, 1 > &global, int &subElement, FieldVector< D, 1 > &local)
Get local coordinates in the subelement.
Definition: refinedsimplexlocalbasis.hh:69
static int getSubElement(const FieldVector< D, 1 > &global)
Get the number of the subelement containing a given point.
Definition: refinedsimplexlocalbasis.hh:53
RefinedSimplexLocalBasis()
Protected default constructor so this class can only be instantiated as a base class.
Definition: refinedsimplexlocalbasis.hh:102
static int getSubElement(const FieldVector< D, 2 > &global)
Get the number of the subtriangle containing a given point.
Definition: refinedsimplexlocalbasis.hh:119
static void getSubElement(const FieldVector< D, 2 > &global, int &subElement, FieldVector< D, 2 > &local)
Get local coordinates in the subtriangle.
Definition: refinedsimplexlocalbasis.hh:137
static void getSubElement(const FieldVector< D, 3 > &global, int &subElement, FieldVector< D, 3 > &local)
Get local coordinates in the subsimplex.
Definition: refinedsimplexlocalbasis.hh:243
RefinedSimplexLocalBasis()
Protected default constructor so this class can only be instantiated as a base class.
Definition: refinedsimplexlocalbasis.hh:183
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:215
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 21, 23:30, 2024)