Loading [MathJax]/extensions/tex2jax.js
ug4
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
mini.h
Go to the documentation of this file.
1/*
2 * Copyright (c) 2012-2016: G-CSC, Goethe University Frankfurt
3 * Author: Arne Nägel
4 *
5 * This file is part of UG4.
6 *
7 * UG4 is free software: you can redistribute it and/or modify it under the
8 * terms of the GNU Lesser General Public License version 3 (as published by the
9 * Free Software Foundation) with the following additional attribution
10 * requirements (according to LGPL/GPL v3 §7):
11 *
12 * (1) The following notice must be displayed in the Appropriate Legal Notices
13 * of covered and combined works: "Based on UG4 (www.ug4.org/license)".
14 *
15 * (2) The following notice must be displayed at a prominent place in the
16 * terminal output of covered works: "Based on UG4 (www.ug4.org/license)".
17 *
18 * (3) The following bibliography is recommended for citation and must be
19 * preserved in all covered files:
20 * "Reiter, S., Vogel, A., Heppner, I., Rupp, M., and Wittum, G. A massively
21 * parallel geometric multigrid solver on hierarchically distributed grids.
22 * Computing and visualization in science 16, 4 (2013), 151-164"
23 * "Vogel, A., Reiter, S., Rupp, M., Nägel, A., and Wittum, G. UG4 -- a novel
24 * flexible software system for simulating pde based models on high performance
25 * computers. Computing and visualization in science 16, 4 (2013), 165-179"
26 *
27 * This program is distributed in the hope that it will be useful,
28 * but WITHOUT ANY WARRANTY; without even the implied warranty of
29 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
30 * GNU Lesser General Public License for more details.
31 */
32
33#ifndef __H__UG__LIB_DISC__LOCAL_SHAPE_FUNCTION_SET__MINI__MINI_BUBBLE__
34#define __H__UG__LIB_DISC__LOCAL_SHAPE_FUNCTION_SET__MINI__MINI_BUBBLE__
35
36#include "../common/lagrange1d.h"
37#include "../local_finite_element_provider.h"
38#include "../local_dof_set.h"
46
47namespace ug{
48
50template <typename TRefElem>
52
53
55template <typename TRefElem>
57{
58 protected:
60 static const int refDim = TRefElem::dim;
61
62 public:
65 {
66 // get _the_ reference element
67 const TRefElem& rRefElem = Provider<TRefElem>::get();
68
69 if(refDim >= 2)
70 {
71 // face (or volume???)
72 // set local DoFs (located at vertices+bubble)
73 nsh = rRefElem.num(0)+1;
74
75 m_vLocalDoF.resize(nsh);
76 for(size_t i = 0; i< nsh-1; ++i)
77 m_vLocalDoF[i] = LocalDoF(0, i, 0);
78
79 m_vLocalDoF[nsh-1] = LocalDoF(refDim, nsh-1, 0); // bubble located at element
80 }
81 else
82 {
83 // edge or vertex
84 nsh = refDim+1;
85 m_vLocalDoF.resize(nsh);
86 for(size_t i = 0; i< nsh-1; ++i)
87 m_vLocalDoF[i] = LocalDoF(0, i, 0);
88
89 }
90 }
91
93 ReferenceObjectID roid() const {return TRefElem::REFERENCE_OBJECT_ID;}
94
96 size_t num_dof() const {return nsh;};
97
99 size_t num_dof(ReferenceObjectID type) const
100 {
101 const int d = ReferenceElementDimension(type);
102 if (d==0) return 1; // vertices
103 if (d == refDim) return 1; // element
104 return 0;
105 }
106
108 const LocalDoF& local_dof(size_t dof) const {return m_vLocalDoF[dof];}
109
111 bool exact_position_available() const {return true;};
112
113 protected:
115 size_t nsh;
116
119
121 std::vector<LocalDoF> m_vLocalDoF;
122};
123
125// ReferenceEdge
126//
127// function space span {1,x}
128//
130
131template <>
133: public MiniBubbleLDS<ReferenceEdge>,
134 public BaseLSFS<MiniBubbleLSFS<ReferenceEdge>, 1>
135{
136 public:
138 static const size_t order = 1;
139
141 static const int dim = 1;
142
144 static const size_t nsh = 3;
145
146 public:
149
152
155
156 public:
159
161 inline bool continuous() const {return true;}
162
164 inline size_t num_sh() const {return nsh;}
165
167 inline bool position(size_t i, MathVector<dim>& pos) const
168 {
169 switch(i)
170 {
171 case 0: pos[0] = 0;return true;
172 case 1: pos[0] = 1;return true;
173 case 2: pos[0] = 0.5;return true; // bubble
174 default: UG_THROW("MiniBubbleLSFS: shape function "<<i<<
175 " not found. Only "<<nsh<<" shapes present.");
176 }
177 }
178
180 inline number shape(const size_t i, const MathVector<dim>& x) const
181 {
182 // ReferenceEdge::check_position(x);
183
184 switch(i)
185 {
186 case 0: return 1-x[0];
187 case 1: return x[0];
188 case 2: return 2.0*x[0]*(1-x[0]);
189 default: UG_THROW("MiniBubbleLSFS: shape function "<<i<<
190 " not found. Only "<<nsh<<" shapes present.");
191 }
192 }
193
195 inline void grad(MathVector<dim>& g, const size_t i, const MathVector<dim>& x) const
196 {
197 // ReferenceEdge::check_position(x);
198
199 switch(i)
200 {
201 case 0: g[0] = -1.0;
202 case 1: g[0] = 1.0;
203 case 2: g[0] = 2.0-4.0*x[0];
204 default: UG_THROW("MiniBubbleLSFS: shape function "<<i<<
205 " not found. Only "<<nsh<<" shapes present.");
206 }
207 }
208};
209
211// ReferenceTriangle
212//
213// function space span {1,x,y}+bubble
214//
216
217template <>
219: public MiniBubbleLDS<ReferenceTriangle>,
220 public BaseLSFS<MiniBubbleLSFS<ReferenceTriangle>, 2>
221{
222 public:
224 static const size_t order = 2;
225
227 static const int dim = 2;
228
230 static const size_t nsh = 4;
231 public:
234
237
240
241 public:
244
246 inline bool continuous() const {return true;}
247
249 inline size_t num_sh() const {return nsh;}
250
252 inline bool position(size_t i, MathVector<dim>& pos) const
253 {
254 switch(i)
255 {
256 case 0: pos[0] = 0.0;
257 pos[1] = 0.0; return true;
258 case 1: pos[0] = 0.1;
259 pos[1] = 0.0; return true;
260 case 2: pos[0] = 0.0;
261 pos[1] = 0.1; return true;
262 case 3: pos[0] = 1.0/3.0;
263 pos[1] = 1.0/3.0; return true;
264 default: UG_THROW("MiniLSFS: shape function "<<i<<
265 " not found. Only "<<nsh<<" shapes present.");
266 }
267 }
268
270 inline number shape(const size_t i, const MathVector<dim>& x) const
271 {
272 // ReferenceTriangle::check_position(x);
273
274 switch(i)
275 {
276 case 0: return (1.0-x[0]-x[1]);
277 case 1: return x[0];
278 case 2: return x[1];
279 case 3: return 27.0*x[0]*x[1]*(1.0-x[0]-x[1]); // bubble
280 default: UG_THROW("MiniLSFS: shape function "<<i<<
281 " not found. Only "<<nsh<<" shapes present.");
282 }
283 }
284
286 inline void grad(MathVector<dim>& g, const size_t i, const MathVector<dim>& x) const
287 {
288 // ReferenceTriangle::check_position(x);
289
290 switch(i)
291 {
292 case 0: g[0] = -1.0;
293 g[1] = -1.0; return;
294 case 1: g[0] = 1.0;
295 g[1] = 0.0; return;
296 case 2: g[0] = 0.0;
297 g[1] = 1.0; return;
298 case 3: g[0] = 27*x[1]*(1.0-x[1]-2.0*x[0]);
299 g[1] = 27*x[0]*(1.0-x[0]-2.0*x[1]);
300 std::cout << "MiniGrad: " << g[0] << ","<< g[1] << "at"<< x[0]<< ","<< x[1]<< std::endl;
301 return;
302 default: UG_THROW("MiniLSFS: shape function "<<i<<
303 " not found. Only "<<nsh<<" shapes present.");
304 }
305 }
306};
307
308
310// ReferenceQuadrilateral
311//
312// function space: span {1,x,y,x^2-y^2} + 2 bubbles
313// ref: Wen Bai, CMAME 143 (1997), 41-47
315
316template <>
318: public MiniBubbleLDS<ReferenceQuadrilateral>,
319 public BaseLSFS<MiniBubbleLSFS<ReferenceQuadrilateral>, 2>
320{
321 protected:
322 static const double SQRT_FIVE;
323 static const double SQRT_FIVTH;
324
325 public:
327 static const size_t order = 2;
328
330 static const int dim = 2;
331
333 static const size_t nsh = 5;
334
335 public:
338
341
344
345 public:
348
350 inline bool continuous() const {return true;}
351
353 inline size_t num_sh() const {return nsh;}
354
356 inline bool position(size_t i, MathVector<dim>& pos) const
357 {
358 switch(i)
359 {
360 case 0: pos[0] = 0.0;
361 pos[1] = 0.0; return true;
362 case 1: pos[0] = 1.0;
363 pos[1] = 0.0; return true;
364 case 2: pos[0] = 1.0;
365 pos[1] = 1.0; return true;
366 case 3: pos[0] = 0.0;
367 pos[1] = 1.0; return true;
368 case 4: pos[0] = 0.5;
369 pos[1] = 0.5; return true;
370 case 5: pos[0] = SQRT_FIVTH;
371 pos[1] = SQRT_FIVTH; return true;
372
373 default: UG_THROW("MiniLSFS: shape function "<<i<<
374 " not found. Only "<<nsh<<" shapes present.");
375 }
376 }
377
379 inline number shape(const size_t i, const MathVector<dim>& x) const
380 {
381 // ReferenceQuadrilateral::check_position(x);
382
383 switch(i)
384 {
385 case 0: return (1.0-x[0])*(1.0-x[1]);
386 case 1: return x[0]*(1.0-x[1]);
387 case 2: return x[0]*x[1];
388 case 3: return x[1]*(1.0-x[0]);
389
390 // two bubbles (condensable)
391 case 4: return x[0]*(1.0-x[0])*x[1]*(1.0-x[1]);
392 case 5: return (1.0-x[0]*x[0])*(1.0-x[1]*x[1])*(x[0]+x[1])/(0.8*0.8*2.0)*SQRT_FIVTH; // max: sqrt(1/5) = sqrt(5)/5
393
394 default: UG_THROW("MiniLSFS: shape function "<<i<<
395 " not found. Only "<<nsh<<" shapes present.");
396 }
397 }
398
400 inline void grad(MathVector<dim>& g, const size_t i, const MathVector<dim>& x) const
401 {
402 // ReferenceQuadrilateral::check_position(x);
403
404 switch(i)
405 {
406 case 0: g[0] = -(1.0-x[1]);
407 g[1] = -(1.0-x[0]); return;
408 case 1: g[0] = (1.0-x[1]);
409 g[1] = -(0.0+x[0]); return;
410 case 2: g[0] = (0.0+x[1]);
411 g[1] = (0.0+x[0]); return;
412 case 3: g[0] = -(0.0+x[1]);
413 g[1] = (1.0-x[0]); return;
414 // bubble 1
415 case 4: g[0] = (1.0-2.0*x[0])*x[1]*(1.0-x[1]);
416 g[1] = x[0]*(1.0-x[0])*(1.0-2.0*x[1]); return;
417 // bubble 2, grad = 3*x^2*(y^2-1) + 2xy*(y^2-1) -(y^2-1) = (y^2-1)*(3*x^2+2xy-1)
418 case 5: g[0] =(x[1]*x[1]-1.0)*(3.0*x[0]*x[0]+2.0*x[0]*x[1]-1.0)/(0.8*0.8*2.0)*SQRT_FIVTH;
419 g[1] =(x[0]*x[0]-1.0)*(3.0*x[1]*x[1]+2.0*x[1]*x[0]-1.0)/(0.8*0.8*2.0)*SQRT_FIVTH; return;
420
421 default: UG_THROW("MiniLSFS: shape function "<<i<<
422 " not found. Only "<<nsh<<" shapes present.");
423 }
424 }
425};
426
428// ReferenceTetrahedron
429//
430// function space span {1,x,y,z} + bubble
431//
433
434template <>
436: public MiniBubbleLDS<ReferenceTetrahedron>,
437 public BaseLSFS<MiniBubbleLSFS<ReferenceTetrahedron>, 3>
438{
439 public:
441 static const size_t order = 1;
442
444 static const int dim = 3;
445
447 static const size_t nsh = 4+1;
448
449 public:
452
455
458
459 public:
462
464 inline bool continuous() const {return true;}
465
467 inline size_t num_sh() const {return nsh;}
468
470 inline bool position(size_t i, MathVector<dim>& pos) const
471 {
472 switch(i)
473 {
474 case 0: pos[0] = 0.0; pos[1] = 0.0; pos[2] = 0.0; return true;
475 case 1: pos[0] = 1.0; pos[1] = 0.0; pos[2] = 0.0; return true;
476 case 2: pos[0] = 0.0; pos[1] = 1.0; pos[2] = 0.0; return true;
477 case 3: pos[0] = 0.0; pos[1] = 0.0; pos[2] = 1.0; return true;
478 case 4: pos[0] = 0.25; pos[1] = 0.25; pos[2] = 0.25; return true;
479 default: UG_THROW("MiniLSFS: shape function "<<i<<
480 " not found. Only "<<nsh<<" shapes present.");
481 }
482 }
483
485 inline number shape(const size_t i, const MathVector<dim>& x) const
486 {
487 // ReferenceTetrahedron::check_position(x);
488 //number prod = x[0]*x[1]*x[2];
489 //number lambda4 = (1.0-x[0]-x[1]-x[2]);
490
491 switch(i)
492 {
493 case 0: return (1.0-x[0]-x[1]-x[2]);
494 case 1: return x[0];
495 case 2: return x[1];
496 case 3: return x[2];
497 case 4: return x[0]*x[1]*x[2]*(1.0-x[0]-x[1]-x[2]); // bubble
498 default: UG_THROW("MiniLSFS: shape function "<<i<<
499 " not found. Only "<<nsh<<" shapes present.");
500 }
501 }
502
504 inline void grad(MathVector<dim>& g, const size_t i, const MathVector<dim>& x) const
505 {
506 // ReferenceTetrahedron::check_position(x);
507 number prod = x[0]*x[1]*x[2];
508 number lambda4 = (1.0-x[0]-x[1]-x[2]);
509
510 switch(i)
511 {
512 case 0: g[0] = -1.0;
513 g[1] = -1.0;
514 g[1] = -1.0; return;
515 case 1: g[0] = 1.0;
516 g[1] = 0.0;
517 g[2] = 0.0; return;
518 case 2: g[0] = 0.0;
519 g[1] = 1.0;
520 g[2] = 0.0; return;
521 case 3: g[0] = 0.0;
522 g[1] = 0.0;
523 g[2] = 1.0; return;
524 case 4: g[0] = -prod + lambda4*x[1]*x[2];
525 g[1] = -prod + lambda4*x[0]*x[2];
526 g[2] = -prod + lambda4*x[0]*x[2]; return;
527
528 default: UG_THROW("MiniLSFS: shape function "<<i<<
529 " not found. Only "<<nsh<<" shapes present.");
530 }
531 }
532};
533
535// ReferenceHexahedron
536//
537// function space span {1,x} \times {1,y} \times {1,z} + bubble
538//
540template <>
542: public MiniBubbleLDS<ReferenceHexahedron>,
543 public BaseLSFS<MiniBubbleLSFS<ReferenceHexahedron>, 3>
544{
545 public:
547 static const size_t order = 6;
548
550 static const int dim = 3;
551
553 static const size_t nsh = 8+1;
554
555 public:
558
561
564
565 public:
568
570 inline bool continuous() const {return true;}
571
573 inline size_t num_sh() const {return nsh;}
574
575
577 inline bool position(size_t i, MathVector<dim>& pos) const
578 {
579 static const DimReferenceElement<3>& refElem
580 = ReferenceElementProvider::get<3>(ROID_HEXAHEDRON);
581
582 if (i<=7)
583 { // corner
584 pos = refElem.corner(i); return true;
585 } else if (i==8)
586 { // bubble
587 pos[0] = 0.5; pos[1] = 0.5; pos[2] = 0.5; return true;
588 } else
589 {
590 UG_THROW("MiniLSFS: shape function "<<i<< " not found. Only "<<nsh<<" shapes present.");
591 }
592 return false;
593 }
594
596 inline number shape(const size_t i, const MathVector<dim>& x) const
597 {
598
599 switch(i)
600 {
601 // corners 0..7
602 case 0: return((1.0-x[0])*(1.0-x[1])*(1.0-x[2]));
603 case 1: return((x[0])*(1.0-x[1])*(1.0-x[2]));
604 case 2: return((x[0])*(x[1])*(1.0-x[2]));
605 case 3: return((1.0-x[0])*(x[1])*(1.0-x[2]));
606 case 4: return((1.0-x[0])*(1.0-x[1])*(x[2]));
607 case 5: return((x[0])*(1.0-x[1])*(x[2]));
608 case 6: return((x[0])*(x[1])*(x[2]));
609 case 7: return((1.0-x[0])*(x[1])*(x[2]));
610 // bubble
611 case 8: return 64.0*x[0]*(1-x[0])*x[1]*(1-x[1])*x[2]*(1-x[2]);
612 default: UG_THROW("MiniLSFS: shape function "<<i<<
613 " not found. Only "<<nsh<<" shapes present.");
614 }
615
616 }
617
619 inline void grad(MathVector<dim>& value, const size_t i, const MathVector<dim>& x) const
620 {
621
622 switch(i)
623
624 {
625 case 0:
626 value[0] = -(1.0-x[1])*(1.0-x[2]);
627 value[1] = -(1.0-x[0])*(1.0-x[2]);
628 value[2] = -(1.0-x[0])*(1.0-x[1]);
629 return;
630 case 1:
631 value[0] = (1.0-x[1])*(1.0-x[2]);
632 value[1] = -(x[0])*(1.0-x[2]);
633 value[2] = -(x[0])*(1.0-x[1]);
634 return;
635 case 2:
636 value[0] = (x[1])*(1.0-x[2]);
637 value[1] = (x[0])*(1.0-x[2]);
638 value[2] = -x[0]*x[1];
639 return;
640 case 3:
641 value[0] = -(x[1])*(1.0-x[2]);
642 value[1] = (1.0-x[0])*(1.0-x[2]);
643 value[2] = -(1.0-x[0])*(x[1]);
644 return;
645 case 4:
646 value[0] = -(1.0-x[1])*(x[2]);
647 value[1] = -(1.0-x[0])*(x[2]);
648 value[2] = (1.0-x[0])*(1.0-x[1]);
649 return;
650 case 5:
651 value[0] = (1.0-x[1])*x[2];
652 value[1] = -(x[0])*x[2];
653 value[2] = (x[0])*(1.0-x[1]);
654 return;
655 case 6:
656 value[0] = (x[1])*x[2];
657 value[1] = (x[0])*x[2];
658 value[2] = x[0]*x[1];
659 return;
660 case 7:
661 value[0] = -(x[1])*x[2];
662 value[1] = (1.0-x[0])*x[2];
663 value[2] = (1.0-x[0])*x[1];
664 return;
665 case 8: // bubble
666 value[0] = 64.0*(1.0-2.0*x[0])*x[1]*(1-x[1])*x[2]*(1-x[2]);
667 value[1] = 64.0*(1.0-2.0*x[1])*x[0]*(1-x[0])*x[2]*(1-x[2]);
668 value[2] = 64.0*(1.0-2.0*x[2])*x[0]*(1-x[0])*x[1]*(1-x[1]);
669 return;
670
671 default: UG_THROW("MiniLSFS: Invalid shape fct index: "<<i);
672 }
673
674
675 }
676};
677
678} //namespace ug
679
680#endif /* __H__UG__LIB_DISC__LOCAL_SHAPE_FUNCTION_SET__MINI__MINI_BUBBLE__ */
681
static interface for trial spaces
Definition local_shape_function_set.h:176
dimension dependent base class for reference elements
Definition reference_element.h:183
const MathVector< dim > & corner(size_t i) const
coordinates of reference corner (i = 0 ... num(0))
Definition reference_element.h:192
Definition local_dof_set.h:58
Definition local_dof_set.h:120
a mathematical Vector with N entries.
Definition math_vector.h:97
MiniBubble Set (2D only!)
Definition mini.h:57
size_t nsh
number of shapes
Definition mini.h:115
MiniBubbleLDS()
constructor
Definition mini.h:64
const LocalDoF & local_dof(size_t dof) const
returns the dof storage
Definition mini.h:108
size_t num_dof(ReferenceObjectID type) const
returns the number of DoFs on a sub-geometric object type
Definition mini.h:99
size_t num_dof() const
returns the total number of DoFs on the finite element
Definition mini.h:96
bool exact_position_available() const
returns if the local dof position are exact
Definition mini.h:111
std::vector< LocalDoF > m_vLocalDoF
association to elements
Definition mini.h:121
ReferenceObjectID roid() const
returns the type of reference element
Definition mini.h:93
static const int refDim
dimension of reference element
Definition mini.h:60
bool position(size_t i, MathVector< dim > &pos) const
local position of DoF i
Definition mini.h:167
void grad(MathVector< dim > &g, const size_t i, const MathVector< dim > &x) const
evaluates the gradient of the shape function
Definition mini.h:195
MathVector< dim > grad_type
Gradient type.
Definition mini.h:151
number shape_type
Shape type.
Definition mini.h:148
MiniBubbleLSFS()
Constructor.
Definition mini.h:158
size_t num_sh() const
Definition mini.h:164
number shape(const size_t i, const MathVector< dim > &x) const
evaluates the shape function
Definition mini.h:180
bool continuous() const
returns if space constructs continuous functions
Definition mini.h:161
ReferenceEdge reference_element_type
Reference Element type.
Definition mini.h:154
MathVector< dim > grad_type
Gradient type.
Definition mini.h:560
void grad(MathVector< dim > &value, const size_t i, const MathVector< dim > &x) const
evaluates the gradient of the shape function
Definition mini.h:619
MiniBubbleLSFS()
Constructor.
Definition mini.h:567
bool position(size_t i, MathVector< dim > &pos) const
local position of DoF i
Definition mini.h:577
size_t num_sh() const
Definition mini.h:573
number shape(const size_t i, const MathVector< dim > &x) const
evaluates the shape function
Definition mini.h:596
number shape_type
Shape type.
Definition mini.h:557
ReferenceHexahedron reference_element_type
Reference Element type.
Definition mini.h:563
bool continuous() const
returns if space constructs continuous functions
Definition mini.h:570
static const double SQRT_FIVTH
Definition mini.h:323
MiniBubbleLSFS()
Constructor.
Definition mini.h:347
bool continuous() const
returns if space constructs continuous functions
Definition mini.h:350
ReferenceQuadrilateral reference_element_type
Reference Element type.
Definition mini.h:343
static const double SQRT_FIVE
Definition mini.h:322
number shape(const size_t i, const MathVector< dim > &x) const
evaluates the shape function
Definition mini.h:379
number shape_type
Shape type.
Definition mini.h:337
size_t num_sh() const
Definition mini.h:353
bool position(size_t i, MathVector< dim > &pos) const
local position of DoF i
Definition mini.h:356
void grad(MathVector< dim > &g, const size_t i, const MathVector< dim > &x) const
evaluates the gradient of the shape function
Definition mini.h:400
MathVector< dim > grad_type
Gradient type.
Definition mini.h:340
number shape(const size_t i, const MathVector< dim > &x) const
evaluates the shape function
Definition mini.h:485
bool position(size_t i, MathVector< dim > &pos) const
local position of DoF i
Definition mini.h:470
ReferenceTetrahedron reference_element_type
Reference Element type.
Definition mini.h:457
number shape_type
Shape type.
Definition mini.h:451
MiniBubbleLSFS()
Constructor.
Definition mini.h:461
bool continuous() const
returns if space constructs continuous functions
Definition mini.h:464
MathVector< dim > grad_type
Gradient type.
Definition mini.h:454
size_t num_sh() const
Definition mini.h:467
void grad(MathVector< dim > &g, const size_t i, const MathVector< dim > &x) const
evaluates the gradient of the shape function
Definition mini.h:504
MathVector< dim > grad_type
Gradient type.
Definition mini.h:236
number shape_type
Shape type.
Definition mini.h:233
bool continuous() const
returns if space constructs continuous functions
Definition mini.h:246
bool position(size_t i, MathVector< dim > &pos) const
local position of DoF i
Definition mini.h:252
size_t num_sh() const
Definition mini.h:249
void grad(MathVector< dim > &g, const size_t i, const MathVector< dim > &x) const
evaluates the gradient of the shape function
Definition mini.h:286
ReferenceTriangle reference_element_type
Reference Element type.
Definition mini.h:239
MiniBubbleLSFS()
Constructor.
Definition mini.h:243
number shape(const size_t i, const MathVector< dim > &x) const
evaluates the shape function
Definition mini.h:270
Lagrange Shape Function Set without virtual functions and fixed order.
Definition mini.h:51
static TClass & get()
returns a singleton based on the identifier
Definition provider.h:54
Definition reference_element.h:345
reference element for a hexahedron
Definition reference_element.h:648
Definition reference_element.h:445
Definition reference_element.h:494
Definition reference_element.h:394
#define UG_THROW(msg)
Definition error.h:57
double number
Definition types.h:124
the ug namespace
ReferenceObjectID
these ids are used to identify the shape of a geometric object.
Definition grid_base_objects.h:74
@ ROID_HEXAHEDRON
Definition grid_base_objects.h:81
int ReferenceElementDimension(ReferenceObjectID roid)
returns the reference element dimension at run-time
Definition reference_element_util.h:40