Dune Core Modules (2.7.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#ifndef DUNE_REFINED_SIMPLEX_LOCALBASIS_HH
4#define DUNE_REFINED_SIMPLEX_LOCALBASIS_HH
5
12#include <dune/localfunctions/common/localbasis.hh>
13
14namespace Dune
15{
16 template<class D, int dim>
17 class RefinedSimplexLocalBasis
18 {
19 protected:
20 RefinedSimplexLocalBasis()
21 {
22 DUNE_THROW(Dune::NotImplemented,"RefinedSimplexLocalBasis not implemented for dim > 3.");
23 }
24 };
25
33 template<class D>
34 class RefinedSimplexLocalBasis<D,1>
35 {
36 protected:
37
40
51 static int getSubElement(const FieldVector<D,1>& global)
52 {
53 if (global[0] <= 0.5)
54 return 0;
55 else if (global[0] <= 1.0)
56 return 1;
57
58 DUNE_THROW(InvalidStateException, "no subelement defined");
59 }
60
67 static void getSubElement(const FieldVector<D,1>& global,
68 int& subElement,
69 FieldVector<D,1>& local)
70 {
71 if (global[0] <= 0.5) {
72 subElement = 0;
73 local[0] = 2.0 * global[0];
74 return;
75 }
76
77 subElement = 1;
78 local[0] = 2.0 * global[0] - 1.0;
79 }
80
81 };
82
83
94 template<class D>
95 class RefinedSimplexLocalBasis<D,2>
96 {
97 protected:
98
101
117 static int getSubElement(const FieldVector<D,2>& global)
118 {
119 if (global[0] + global[1] <= 0.5)
120 return 0;
121 else if (global[0] >= 0.5)
122 return 1;
123 else if (global[1] >= 0.5)
124 return 2;
125
126 return 3;
127 }
128
135 static void getSubElement(const FieldVector<D,2>& global,
136 int& subElement,
137 FieldVector<D,2>& local)
138 {
139 if (global[0] + global[1] <= 0.5) {
140 subElement = 0;
141 local[0] = 2*global[0];
142 local[1] = 2*global[1];
143 return;
144 } else if (global[0] >= 0.5) {
145 subElement = 1;
146 local[0] = 2*global[0]-1;
147 local[1] = 2*global[1];
148 return;
149 } else if (global[1] >= 0.5) {
150 subElement = 2;
151 local[0] = 2*global[0];
152 local[1] = 2*global[1]-1;
153 return;
154 }
155
156 subElement = 3;
157 local[0] = -2 * global[0] + 1;
158 local[1] = -2 * global[1] + 1;
159
160 }
161
162
163 };
164
175 template<class D>
176 class RefinedSimplexLocalBasis<D,3>
177 {
178 protected:
179
182
213 static int getSubElement(const FieldVector<D,3>& global)
214 {
215 if (global[0] + global[1] + global[2] <= 0.5)
216 return 0;
217 else if (global[0] >= 0.5)
218 return 1;
219 else if (global[1] >= 0.5)
220 return 2;
221 else if (global[2] >= 0.5)
222 return 3;
223 else if ((global[0] + global[1] <= 0.5)and (global[1] + global[2] <= 0.5))
224 return 4;
225 else if ((global[0] + global[1] >= 0.5)and (global[1] + global[2] <= 0.5))
226 return 5;
227 else if ((global[0] + global[1] <= 0.5)and (global[1] + global[2] >= 0.5))
228 return 6;
229 else if ((global[0] + global[1] >= 0.5)and (global[1] + global[2] >= 0.5))
230 return 7;
231
232 DUNE_THROW(InvalidStateException, "no subelement defined");
233
234 }
241 static void getSubElement(const FieldVector<D,3>& global,
242 int& subElement,
243 FieldVector<D,3>& local)
244 {
245 if (global[0] + global[1] + global[2] <= 0.5) {
246 subElement = 0;
247 local = global;
248 local *= 2.0;
249 return;
250 } else if (global[0] >= 0.5) {
251 subElement = 1;
252 local = global;
253 local[0] -= 0.5;
254 local *= 2.0;
255 return;
256 } else if (global[1] >= 0.5) {
257 subElement = 2;
258 local = global;
259 local[1] -= 0.5;
260 local *= 2.0;
261 return;
262 } else if (global[2] >= 0.5) {
263 subElement = 3;
264 local = global;
265 local[2] -= 0.5;
266 local *= 2.0;
267 return;
268 } else if ((global[0] + global[1] <= 0.5)and (global[1] + global[2] <= 0.5)) {
269 subElement = 4;
270 local[0] = 2.0 * global[1];
271 local[1] = 2.0 * (0.5 - global[0] - global[1]);
272 local[2] = 2.0 * (-0.5 + global[0] + global[1] + global[2]);
273 // Dune::FieldMatrix<double,3,3> A(0.0);
274 // A[0][1] = 2.0;
275 // A[1][0] = -2.0;
276 // A[1][1] = -2.0;
277 // A[2][0] = 2.0;
278 // A[2][1] = 2.0;
279 // A[2][2] = 2.0;
280 // A.mv(global,local);
281 // local[1] += 1.0;
282 // local[2] -= 1.0;
283 return;
284 } else if ((global[0] + global[1] >= 0.5)and (global[1] + global[2] <= 0.5)) {
285 subElement = 5;
286 local[0] = 2.0 * (0.5 - global[0]);
287 local[1] = 2.0 * (0.5 - global[1] - global[2]);
288 local[2] = 2.0 * global[2];
289 // Dune::FieldMatrix<double,3,3> A(0.0);
290 // A[0][0] = -2.0;
291 // A[1][1] = -2.0;
292 // A[1][2] = -2.0;
293 // A[2][2] = 2.0;
294 // A.mv(global,local);
295 // local[0] += 1.0;
296 // local[1] += 1.0;
297 return;
298 } else if ((global[0] + global[1] <= 0.5)and (global[1] + global[2] >= 0.5)) {
299 subElement = 6;
300 local[0] = 2.0 * (0.5 - global[0] - global[1]);
301 local[1] = 2.0 * global[0];
302 local[2] = 2.0 * (-0.5 + global[1] + global[2]);
303 // Dune::FieldMatrix<double,3,3> A(0.0);
304 // A[0][0] = -2.0;
305 // A[0][1] = -2.0;
306 // A[1][0] = 2.0;
307 // A[2][1] = 2.0;
308 // A[2][2] = 2.0;
309 // A.mv(global,local);
310 // local[0] += 1.0;
311 // local[2] -= 1.0;
312 return;
313 } else if ((global[0] + global[1] >= 0.5)and (global[1] + global[2] >= 0.5)) {
314 subElement = 7;
315 local[0] = 2.0 * (-0.5 + global[1] + global[2]);
316 local[1] = 2.0 * (0.5 - global[1]);
317 local[2] = 2.0 * (-0.5 + global[0] + global[1]);
318 // Dune::FieldMatrix<double,3,3> A(0.0);
319 // A[0][1] = 2.0;
320 // A[0][2] = 2.0;
321 // A[1][1] = -2.0;
322 // A[2][0] = 2.0;
323 // A[2][1] = 2.0;
324 // A.mv(global,local);
325 // local[0] -= 1.0;
326 // local[1] += 1.0;
327 // local[2] -= 1.0;
328 return;
329 }
330
331 DUNE_THROW(InvalidStateException, "no subelement defined");
332
333 }
334
335 };
336
337
338}
339
340#endif
vector space out of a tensor product of fields.
Definition: fvector.hh:96
Default exception if a function was called while the object is not in a valid state for that function...
Definition: exceptions.hh:279
Default exception for dummy implementations.
Definition: exceptions.hh:261
RefinedSimplexLocalBasis()
Protected default constructor so this class can only be instantiated as a base class.
Definition: refinedsimplexlocalbasis.hh:39
static void getSubElement(const FieldVector< D, 1 > &global, int &subElement, FieldVector< D, 1 > &local)
Get local coordinates in the subelement.
Definition: refinedsimplexlocalbasis.hh:67
static int getSubElement(const FieldVector< D, 1 > &global)
Get the number of the subelement containing a given point.
Definition: refinedsimplexlocalbasis.hh:51
RefinedSimplexLocalBasis()
Protected default constructor so this class can only be instantiated as a base class.
Definition: refinedsimplexlocalbasis.hh:100
static int getSubElement(const FieldVector< D, 2 > &global)
Get the number of the subtriangle containing a given point.
Definition: refinedsimplexlocalbasis.hh:117
static void getSubElement(const FieldVector< D, 2 > &global, int &subElement, FieldVector< D, 2 > &local)
Get local coordinates in the subtriangle.
Definition: refinedsimplexlocalbasis.hh:135
static void getSubElement(const FieldVector< D, 3 > &global, int &subElement, FieldVector< D, 3 > &local)
Get local coordinates in the subsimplex.
Definition: refinedsimplexlocalbasis.hh:241
RefinedSimplexLocalBasis()
Protected default constructor so this class can only be instantiated as a base class.
Definition: refinedsimplexlocalbasis.hh:181
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:213
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:216
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)