ug4
crouzeix_raviart.h
Go to the documentation of this file.
1 /*
2  * Copyright (c) 2012-2015: G-CSC, Goethe University Frankfurt
3  * Author: Christian Wehner
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__CROUZEIX_RAVIART__
34 #define __H__UG__LIB_DISC__LOCAL_SHAPE_FUNCTION_SET__CROUZEIX_RAVIART__
35 
39 
40 namespace ug{
41 
42 
44 template <typename TRefElem>
46 {
47  public:
49  static const int dim = TRefElem::dim;
50 
52  static const size_t nsh = TRefElem::numSides;
53 
54  public:
57  {
58  for(size_t i = 0; i < nsh; ++i)
59  m_vLocalDoF[i] = LocalDoF(dim-1, i, 0);
60  }
61 
63  ReferenceObjectID roid() const {return TRefElem::REFERENCE_OBJECT_ID;}
64 
66  size_t num_dof() const {return nsh;};
67  size_t num_sh() const {return nsh;}
68 
70  size_t num_dof(ReferenceObjectID type) const
71  {
72  if(ReferenceElementDimension(type) == dim-1) return 1;
73  else return 0;
74  }
75 
77  const LocalDoF& local_dof(size_t dof) const {return m_vLocalDoF[dof];}
78 
80  bool exact_position_available() const {return true;};
81 
83  bool continuous() const {return false;}
84 
85  protected:
88 };
89 
91 template <typename TRefElem>
93 
95 // ReferenceTriangle
96 //
97 // function space span {1,x,y}
98 //
100 
101 template <>
103 : public CrouzeixRaviartBase<ReferenceTriangle>,
104  public BaseLSFS<CrouzeixRaviartLSFS<ReferenceTriangle>, 2>
105 {
106  public:
108  static const int dim = 2;
109 
110  public:
112  bool position(size_t i, MathVector<dim>& pos) const
113  {
114  switch(i)
115  {
116  case 0: pos[0] = 0.5;
117  pos[1] = 0.0; return true;
118  case 1: pos[0] = 0.5;
119  pos[1] = 0.5; return true;
120  case 2: pos[0] = 0.0;
121  pos[1] = 0.5; return true;
122  default: UG_THROW("CrouzeixRaviartLSFS: shape function "<<i<<
123  " not found. Only "<<nsh<<" shapes present.");
124  }
125  }
126 
128  number shape(const size_t i, const MathVector<dim>& x) const
129  {
130  switch(i)
131  {
132  case 0: return 1-2*x[1];
133  case 1: return -1+2*x[0]+2*x[1];
134  case 2: return 1-2*x[0];
135  default: UG_THROW("CrouzeixRaviartLSFS: shape function "<<i<<
136  " not found. Only "<<nsh<<" shapes present.");
137  }
138  }
139 
141  void grad(MathVector<dim>& g, const size_t i, const MathVector<dim>& x) const
142  {
143  switch(i)
144  {
145  case 0: g[0] = 0.0;
146  g[1] = -2.0; return;
147  case 1: g[0] = 2.0;
148  g[1] = 2.0; return;
149  case 2: g[0] = -2.0;
150  g[1] = 0.0; return;
151  default: UG_THROW("CrouzeixRaviartLSFS: shape function "<<i<<
152  " not found. Only "<<nsh<<" shapes present.");
153  }
154  }
155 };
156 
157 
159 // ReferenceQuadrilateral
160 //
161 // function space span {1,x,y,x^2-y^2}
162 //
164 
165 template <>
167 : public CrouzeixRaviartBase<ReferenceQuadrilateral>,
168  public BaseLSFS<CrouzeixRaviartLSFS<ReferenceQuadrilateral>, 2>
169 {
170  public:
172  static const int dim = 2;
173 
174  public:
176  bool position(size_t i, MathVector<dim>& pos) const
177  {
178  switch(i)
179  {
180  case 0: pos[0] = 0.5;
181  pos[1] = 0.0; return true;
182  case 1: pos[0] = 1.0;
183  pos[1] = 0.5; return true;
184  case 2: pos[0] = 0.5;
185  pos[1] = 1.0; return true;
186  case 3: pos[0] = 0.0;
187  pos[1] = 0.5; return true;
188  return true;
189  default: UG_THROW("CrouzeixRaviartLSFS: shape function "<<i<<
190  " not found. Only "<<nsh<<" shapes present.");
191  }
192  }
193 
195  number shape(const size_t i, const MathVector<dim>& x) const
196  {
197  switch(i)
198  {
199  case 0: return 0.75+x[0]-2*x[1]-x[0]*x[0]+x[1]*x[1];
200  case 1: return -0.25+x[1]+x[0]*x[0]-x[1]*x[1];
201  case 2: return -0.25+x[0]-x[0]*x[0]+x[1]*x[1];
202  case 3: return 0.75-2*x[0]+x[1]+x[0]*x[0]-x[1]*x[1];
203  default: UG_THROW("CrouzeixRaviartLSFS: shape function "<<i<<
204  " not found. Only "<<nsh<<" shapes present.");
205  }
206  }
207 
209  void grad(MathVector<dim>& g, const size_t i, const MathVector<dim>& x) const
210  {
211  switch(i)
212  {
213  case 0: g[0] = 1-2*x[0];
214  g[1] = -2.0+2*x[1]; return;
215  case 1: g[0] = 2.0*x[0];
216  g[1] = 1-2.0*x[1]; return;
217  case 2: g[0] = 1-2.0*x[0];
218  g[1] = 2*x[1]; return;
219  case 3: g[0] = -2+2*x[0];
220  g[1] = 1-2*x[1]; return;
221  default: UG_THROW("CrouzeixRaviartLSFS: shape function "<<i<<
222  " not found. Only "<<nsh<<" shapes present.");
223  }
224  }
225 };
226 
228 // ReferenceTetrahedron
229 //
230 // function space span {1,x,y,z}
231 //
233 
234 template <>
236 : public CrouzeixRaviartBase<ReferenceTetrahedron>,
237  public BaseLSFS<CrouzeixRaviartLSFS<ReferenceTetrahedron>, 3>
238 {
239  public:
241  static const int dim = ReferenceTetrahedron::dim;
242 
243  public:
245  bool position(size_t i, MathVector<dim>& pos) const
246  {
247  switch(i)
248  {
249  case 0: pos[0] = 1.0/3.0;
250  pos[1] = 1.0/3.0;
251  pos[2] = 0.0; return true;
252  case 1: pos[0] = 1.0/3.0;
253  pos[1] = 1.0/3.0;
254  pos[2] = 1.0/3.0; return true;
255  case 2: pos[0] = 0.0;
256  pos[1] = 1.0/3.0;
257  pos[2] = 1.0/3.0; return true;
258  case 3: pos[0] = 1.0/3.0;
259  pos[1] = 0.0;
260  pos[2] = 1.0/3.0; return true;
261  return true;
262  default: UG_THROW("CrouzeixRaviartLSFS: shape function "<<i<<
263  " not found. Only "<<nsh<<" shapes present.");
264  }
265  }
266 
268  number shape(const size_t i, const MathVector<dim>& x) const
269  {
270  switch(i)
271  {
272  case 0: return 1 - 3*x[2];;
273  case 1: return -2 + 3*x[0]+3*x[1]+3*x[2];
274  case 2: return 1 - 3*x[0];
275  case 3: return 1 - 3*x[1];
276  default: UG_THROW("CrouzeixRaviartLSFS: shape function "<<i<<
277  " not found. Only "<<nsh<<" shapes present.");
278  }
279  }
280 
282  void grad(MathVector<dim>& g, const size_t i, const MathVector<dim>& x) const
283  {
284  switch(i)
285  {
286  case 0: g[0] = 0.0;
287  g[1] = 0.0;
288  g[2] = -3.0;return;
289  case 1: g[0] = 3.0;
290  g[1] = 3.0;
291  g[2] = 3.0;return;
292  case 2: g[0] = -3.0;
293  g[1] = 0.0;
294  g[2] = 0.0;return;
295  case 3: g[0] = 0.0;
296  g[1] = -3;
297  g[2] = 0.0;return;
298  default: UG_THROW("CrouzeixRaviartLSFS: shape function "<<i<<
299  " not found. Only "<<nsh<<" shapes present.");
300  }
301  }
302 };
303 
304 
306 // ReferenceHexahedron
307 //
308 // function space span {1,x,y,z,x^2-y^2,y^2-z^2}
309 //
311 
312 template <>
314 : public CrouzeixRaviartBase<ReferenceHexahedron>,
315  public BaseLSFS<CrouzeixRaviartLSFS<ReferenceHexahedron>, 3>
316 {
317  public:
319  static const int dim = ReferenceHexahedron::dim;
320 
321  public:
323  bool position(size_t i, MathVector<dim>& pos) const
324  {
325  switch(i)
326  {
327  case 0: pos[0] = 0.5;
328  pos[1] = 0.5;
329  pos[2] = 0.0;return true;
330  case 1: pos[0] = 0.5;
331  pos[1] = 0.0;
332  pos[2] = 0.5; return true;
333  case 2: pos[0] = 1.0;
334  pos[1] = 0.5;
335  pos[2] = 0.5; return true;
336  case 3: pos[0] = 0.5;
337  pos[1] = 1.0;
338  pos[2] = 0.5; return true;
339  case 4: pos[0] = 0.0;
340  pos[1] = 0.5;
341  pos[2] = 0.5; return true;
342  case 5: pos[0] = 0.5;
343  pos[1] = 0.5;
344  pos[2] = 1.0; return true;
345  return true;
346  default: UG_THROW("CrouzeixRaviartLSFS: shape function "<<i<<
347  " not found. Only "<<nsh<<" shapes present.");
348  }
349  }
350 
352  number shape(const size_t i, const MathVector<dim>& x) const
353  {
354  switch(i)
355  {
356  case 0: return 1.0/3.0*(2 + 2*x[0] + 2*x[1] - 7*x[2]-2*x[0]*x[0]-2*x[1]*x[1]+4*x[2]*x[2]);
357  case 1: return 1.0/3.0*(2+2*x[0]-7*x[1]+2*x[2]-2*x[0]*x[0]+4*x[1]*x[1]-2*x[2]*x[2]);
358  case 2: return 1.0/3.0*(-1-x[0]+2*x[1]+2*x[2]+4*x[0]*x[0]-2*x[1]*x[1]-2*x[2]*x[2]);
359  case 3: return 1.0/3.0*(-1+2*x[0]-x[1]+2*x[2]-2*x[0]*x[0]+4*x[1]*x[1]-2*x[2]*x[2]);
360  case 4: return 1.0/3.0*(2-7*x[0]+2*x[1]+2*x[2]+4*x[0]*x[0]-2*x[1]*x[1]-2*x[2]*x[2]);
361  case 5: return 1.0/3.0*(-1+2*x[0]+2*x[1]-x[2]-2*x[0]*x[0]-2*x[1]*x[1]+4*x[2]*x[2]);
362  default: UG_THROW("CrouzeixRaviartLSFS: shape function "<<i<<
363  " not found. Only "<<nsh<<" shapes present.");
364  }
365  }
366 
368  void grad(MathVector<dim>& g, const size_t i, const MathVector<dim>& x) const
369  {
370  switch(i)
371  {
372  case 0: g[0] = 1.0/3.0*(2-4*x[0]);
373  g[1] = 1.0/3.0*(2-4*x[1]);
374  g[2] = 1.0/3.0*(-7+8*x[2]);return;
375  case 1: g[0] = 1.0/3.0*(2-4*x[0]);
376  g[1] = 1.0/3.0*(-7+8*x[1]);
377  g[2] = 1.0/3.0*(2-4*x[2]);return;
378  case 2: g[0] = 1.0/3.0*(-1+8*x[0]);
379  g[1] = 1.0/3.0*(2-4*x[1]);
380  g[2] = 1.0/3.0*(2-4*x[2]);return;
381  case 3: g[0] = 1.0/3.0*(2-4*x[0]);
382  g[1] = 1.0/3.0*(-1+8*x[1]);
383  g[2] = 1.0/3.0*(2-4*x[2]);return;
384  case 4: g[0] = 1.0/3.0*(-7+8*x[0]);
385  g[1] = 1.0/3.0*(2-4*x[1]);
386  g[2] = 1.0/3.0*(2-4*x[2]);return;
387  case 5: g[0] = 1.0/3.0*(2-4*x[0]);
388  g[1] = 1.0/3.0*(2-4*x[1]);
389  g[2] = 1.0/3.0*(-1+8*x[2]);return;
390  default: UG_THROW("CrouzeixRaviartLSFS: shape function "<<i<<
391  " not found. Only "<<nsh<<" shapes present.");
392  }
393  }
394 };
395 
397 // ReferencePrism
398 //
399 // function space span {1,x,y,z,x^2-y^2-z^2}
400 //
402 
403 template <>
405 : public CrouzeixRaviartBase<ReferencePrism>,
406  public BaseLSFS<CrouzeixRaviartLSFS<ReferencePrism>, 3>
407 {
408  public:
410  static const int dim = 3;
411 
412  public:
414  bool position(size_t i, MathVector<dim>& pos) const
415  {
416  switch(i)
417  {
418  case 0: pos[0] = 1.0/3.0;
419  pos[1] = 1.0/3.0;
420  pos[2] = 0.0;return true;
421  case 1: pos[0] = 0.5;
422  pos[1] = 0.0;
423  pos[2] = 0.5; return true;
424  case 2: pos[0] = 0.5;
425  pos[1] = 0.5;
426  pos[2] = 0.5; return true;
427  case 3: pos[0] = 0.0;
428  pos[1] = 0.5;
429  pos[2] = 0.5; return true;
430  case 4: pos[0] = 1.0/3.0;
431  pos[1] = 1.0/3.0;
432  pos[2] = 1.0; return true;
433  default: UG_THROW("CrouzeixRaviartLSFS: shape function "<<i<<
434  " not found. Only "<<nsh<<" shapes present.");
435  }
436  }
437 
439  number shape(const size_t i, const MathVector<dim>& x) const
440  {
441  switch(i)
442  {
443  case 0: return 1+x[0]-x[1]-3*x[2]-2*(x[0]*x[0]-x[1]*x[1]-x[2]*x[2]);
444  case 1: return 1.0/3.0*(2-2*x[0]-4*x[1]+4*x[2]+4*(x[0]*x[0]-x[1]*x[1]-x[2]*x[2]));
445  case 2: return 1.0/3.0*(-4+4*x[0]+8*x[1]+4*x[2]+4*(x[0]*x[0]-x[1]*x[1]-x[2]*x[2]));
446  case 3: return 1.0/3.0*(2-8*x[0]+2*x[1]+4*x[2]+4*(x[0]*x[0]-x[1]*x[1]-x[2]*x[2]));
447  case 4: return x[0]-x[1]-x[2]-2*(x[0]*x[0]-x[1]*x[1]-x[2]*x[2]);
448  default: UG_THROW("CrouzeixRaviartLSFS: shape function "<<i<<
449  " not found. Only "<<nsh<<" shapes present.");
450  }
451  }
452 
454  void grad(MathVector<dim>& g, const size_t i, const MathVector<dim>& x) const
455  {
456  switch(i)
457  {
458  case 0: g[0] = 1-4*x[0];
459  g[1] = -1+4*x[1];
460  g[2] = -3+4*x[2];return;
461  case 1: g[0] = 1.0/3.0*(-2+8*x[0]);
462  g[1] = 1.0/3.0*(-4-8*x[1]);
463  g[2] = 1.0/3.0*(4-8*x[2]);return;
464  case 2: g[0] = 1.0/3.0*(4+8*x[0]);
465  g[1] = 1.0/3.0*(8-8*x[1]);
466  g[2] = 1.0/3.0*(4-8*x[2]);return;
467  case 3: g[0] = 1.0/3.0*(-8+8*x[0]);
468  g[1] = 1.0/3.0*(2-8*x[1]);
469  g[2] = 1.0/3.0*(4-8*x[2]);return;
470  case 4: g[0] = 1-4*x[0];
471  g[1] = -1+4*x[1];
472  g[2] = -1+4*x[2];return;
473  default: UG_THROW("CrouzeixRaviartLSFS: shape function "<<i<<
474  " not found. Only "<<nsh<<" shapes present.");
475  }
476  }
477 };
478 
479 
481 // ReferencePyramid
482 //
483 // function space span {1,x,y,z,x^2-y^2-z^2}
484 //
486 
487 template <>
489 : public CrouzeixRaviartBase<ReferencePyramid>,
490  public BaseLSFS<CrouzeixRaviartLSFS<ReferencePyramid>, 3>
491 {
492  public:
494  static const int dim = 3;
495 
496  public:
498  bool position(size_t i, MathVector<dim>& pos) const
499  {
500  switch(i)
501  {
502  case 0: pos[0] = 0.5;
503  pos[1] = 0.5;
504  pos[2] = 0.0;return true;
505  case 1: pos[0] = 1.0/3.0;
506  pos[1] = 0.0;
507  pos[2] = 1.0/3.0; return true;
508  case 2: pos[0] = 2.0/3.0;
509  pos[1] = 1.0/3.0;
510  pos[2] = 1.0/3.0; return true;
511  case 3: pos[0] = 1.0/3.0;
512  pos[1] = 2.0/3.0;
513  pos[2] = 1.0/3.0; return true;
514  case 4: pos[0] = 0.0;
515  pos[1] = 1.0/3.0;
516  pos[2] = 1.0/3.0; return true;
517  default: UG_THROW("CrouzeixRaviartLSFS: shape function "<<i<<
518  " not found. Only "<<nsh<<" shapes present.");
519  }
520  }
521 
523  number shape(const size_t i, const MathVector<dim>& x) const
524  {
525  switch(i)
526  {
527  case 0: return 1-3*x[2];
528  case 1: return 0.75+1.5*x[0]-3*x[1]-0.75*x[2]-2.25*(x[0]*x[0]-x[1]*x[1]-x[2]*x[2]);
529  case 2: return -0.75+1.5*x[1]+2.25*x[2]+2.25*(x[0]*x[0]-x[1]*x[1]-x[2]*x[2]);
530  case 3: return -0.75+1.5*x[0]+0.75*x[2]-2.25*(x[0]*x[0]-x[1]*x[1]-x[2]*x[2]);
531  case 4: return 0.75-3*x[0]+1.5*x[1]+0.75*x[2]+2.25*(x[0]*x[0]-x[1]*x[1]-x[2]*x[2]);
532  default: UG_THROW("CrouzeixRaviartLSFS: shape function "<<i<<
533  " not found. Only "<<nsh<<" shapes present.");
534  }
535  }
536 
538  void grad(MathVector<dim>& g, const size_t i, const MathVector<dim>& x) const
539  {
540  switch(i)
541  {
542  case 0: g[0] = 0.0;
543  g[1] = 0.0;
544  g[2] = -3.0;return;
545  case 1: g[0] = 1.5-4.5*x[0];
546  g[1] = -3+4.5*x[1];
547  g[2] = -0.75+4.5*x[2];return;
548  case 2: g[0] = 4.5*x[0];
549  g[1] = 1.5-4.5*x[1];
550  g[2] = 2.25-4.5*x[2];return;
551  case 3: g[0] = 1.5-4.5*x[0];
552  g[1] = 4.5*x[1];
553  g[2] = 0.75+4.5*x[2];return;
554  case 4: g[0] = -3+4.5*x[0];
555  g[1] = 1.5-4.5*x[1];
556  g[2] = 0.75-4.5*x[2];return;
557  default: UG_THROW("CrouzeixRaviartLSFS: shape function "<<i<<
558  " not found. Only "<<nsh<<" shapes present.");
559  }
560  }
561 };
562 
563 
564 } //namespace ug
565 
566 #endif /* __H__UG__LIB_DISC__LOCAL_SHAPE_FUNCTION_SET__CROUZEIX_RAVIART__ */
567 
static interface for trial spaces
Definition: local_shape_function_set.h:176
Crouzeix - Raviart Set.
Definition: crouzeix_raviart.h:46
const LocalDoF & local_dof(size_t dof) const
returns the dof storage
Definition: crouzeix_raviart.h:77
static const size_t nsh
Number of shape functions.
Definition: crouzeix_raviart.h:52
bool exact_position_available() const
returns if the local dof position are exact
Definition: crouzeix_raviart.h:80
ReferenceObjectID roid() const
returns the type of reference element
Definition: crouzeix_raviart.h:63
LocalDoF m_vLocalDoF[nsh]
association to elements
Definition: crouzeix_raviart.h:87
CrouzeixRaviartBase()
constructor
Definition: crouzeix_raviart.h:56
static const int dim
dimension of reference element
Definition: crouzeix_raviart.h:49
bool continuous() const
returns if space constructs continuous functions
Definition: crouzeix_raviart.h:83
size_t num_dof() const
returns the total number of DoFs on the finite element
Definition: crouzeix_raviart.h:66
size_t num_dof(ReferenceObjectID type) const
returns the number of DoFs on a sub-geometric object type
Definition: crouzeix_raviart.h:70
size_t num_sh() const
Definition: crouzeix_raviart.h:67
number shape(const size_t i, const MathVector< dim > &x) const
evaluates the shape function
Definition: crouzeix_raviart.h:352
void grad(MathVector< dim > &g, const size_t i, const MathVector< dim > &x) const
evaluates the gradient of the shape function
Definition: crouzeix_raviart.h:368
bool position(size_t i, MathVector< dim > &pos) const
local position of DoF i
Definition: crouzeix_raviart.h:323
void grad(MathVector< dim > &g, const size_t i, const MathVector< dim > &x) const
evaluates the gradient of the shape function
Definition: crouzeix_raviart.h:454
number shape(const size_t i, const MathVector< dim > &x) const
evaluates the shape function
Definition: crouzeix_raviart.h:439
bool position(size_t i, MathVector< dim > &pos) const
local position of DoF i
Definition: crouzeix_raviart.h:414
bool position(size_t i, MathVector< dim > &pos) const
local position of DoF i
Definition: crouzeix_raviart.h:498
number shape(const size_t i, const MathVector< dim > &x) const
evaluates the shape function
Definition: crouzeix_raviart.h:523
void grad(MathVector< dim > &g, const size_t i, const MathVector< dim > &x) const
evaluates the gradient of the shape function
Definition: crouzeix_raviart.h:538
void grad(MathVector< dim > &g, const size_t i, const MathVector< dim > &x) const
evaluates the gradient of the shape function
Definition: crouzeix_raviart.h:209
bool position(size_t i, MathVector< dim > &pos) const
local position of DoF i
Definition: crouzeix_raviart.h:176
number shape(const size_t i, const MathVector< dim > &x) const
evaluates the shape function
Definition: crouzeix_raviart.h:195
void grad(MathVector< dim > &g, const size_t i, const MathVector< dim > &x) const
evaluates the gradient of the shape function
Definition: crouzeix_raviart.h:282
number shape(const size_t i, const MathVector< dim > &x) const
evaluates the shape function
Definition: crouzeix_raviart.h:268
bool position(size_t i, MathVector< dim > &pos) const
local position of DoF i
Definition: crouzeix_raviart.h:245
number shape(const size_t i, const MathVector< dim > &x) const
evaluates the shape function
Definition: crouzeix_raviart.h:128
void grad(MathVector< dim > &g, const size_t i, const MathVector< dim > &x) const
evaluates the gradient of the shape function
Definition: crouzeix_raviart.h:141
bool position(size_t i, MathVector< dim > &pos) const
local position of DoF i
Definition: crouzeix_raviart.h:112
Lagrange Shape Function Set without virtual functions and fixed order.
Definition: crouzeix_raviart.h:92
Definition: local_dof_set.h:58
reference element for a hexahedron
Definition: reference_element.h:648
static const int dim
dimension of reference element
Definition: reference_element.h:654
Definition: reference_element.h:595
Definition: reference_element.h:546
Definition: reference_element.h:445
Definition: reference_element.h:494
static const int dim
dimension of reference element
Definition: reference_element.h:500
Definition: reference_element.h:394
static const int dim
#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
int ReferenceElementDimension(ReferenceObjectID roid)
returns the reference element dimension at run-time
Definition: reference_element_util.h:40