libMesh
fe_nedelec_one_shape_3D.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2019 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 // C++ includes
20 
21 // Local includes
22 #include "libmesh/fe.h"
23 #include "libmesh/elem.h"
24 
25 namespace libMesh
26 {
27 
28 template <>
30  const Order,
31  const unsigned int,
32  const Point &)
33 {
34  libmesh_error_msg("Nedelec elements require the element type \nbecause edge orientation is needed.");
35  return RealGradient();
36 }
37 
38 
39 
40 template <>
42  const Order order,
43  const unsigned int i,
44  const Point & p,
45  const bool add_p_level)
46 {
47 #if LIBMESH_DIM == 3
48  libmesh_assert(elem);
49 
50  const Order totalorder = static_cast<Order>(order + add_p_level * elem->p_level());
51 
52  switch (totalorder)
53  {
54  // linear Lagrange shape functions
55  case FIRST:
56  {
57  switch (elem->type())
58  {
59  case HEX20:
60  case HEX27:
61  {
62  libmesh_assert_less (i, 12);
63 
64  const Real xi = p(0);
65  const Real eta = p(1);
66  const Real zeta = p(2);
67 
68  // Even with a loose inverse_map tolerance we ought to
69  // be nearly on the element interior in master
70  // coordinates
71  libmesh_assert_less_equal ( std::fabs(xi), 1.0+10*TOLERANCE );
72  libmesh_assert_less_equal ( std::fabs(eta), 1.0+10*TOLERANCE );
73  libmesh_assert_less_equal ( std::fabs(zeta), 1.0+10*TOLERANCE );
74 
75  switch(i)
76  {
77  case 0:
78  {
79  if (elem->point(0) > elem->point(1))
80  return RealGradient( -0.125*(1.0-eta-zeta+eta*zeta), 0.0, 0.0 );
81  else
82  return RealGradient( 0.125*(1.0-eta-zeta+eta*zeta), 0.0, 0.0 );
83  }
84  case 1:
85  {
86  if (elem->point(1) > elem->point(2))
87  return RealGradient( 0.0, -0.125*(1.0+xi-zeta-xi*zeta), 0.0 );
88  else
89  return RealGradient( 0.0, 0.125*(1.0+xi-zeta-xi*zeta), 0.0 );
90  }
91  case 2:
92  {
93  if (elem->point(2) > elem->point(3))
94  return RealGradient( 0.125*(1.0+eta-zeta-eta*zeta), 0.0, 0.0 );
95  else
96  return RealGradient( -0.125*(1.0+eta-zeta-eta*zeta), 0.0, 0.0 );
97  }
98  case 3:
99  {
100  if (elem->point(3) > elem->point(0))
101  return RealGradient( 0.0, 0.125*(1.0-xi-zeta+xi*zeta), 0.0 );
102  else
103  return RealGradient( 0.0, -0.125*(1.0-xi-zeta+xi*zeta), 0.0 );
104  }
105  case 4:
106  {
107  if (elem->point(0) > elem->point(4))
108  return RealGradient( 0.0, 0.0, -0.125*(1.0-xi-eta+xi*eta) );
109  else
110  return RealGradient( 0.0, 0.0, 0.125*(1.0-xi-eta+xi*eta) );
111  }
112  case 5:
113  {
114  if (elem->point(1) > elem->point(5))
115  return RealGradient( 0.0, 0.0, -0.125*(1.0+xi-eta-xi*eta) );
116  else
117  return RealGradient( 0.0, 0.0, 0.125*(1.0+xi-eta-xi*eta) );
118  }
119  case 6:
120  {
121  if (elem->point(2) > elem->point(6))
122  return RealGradient( 0.0, 0.0, -0.125*(1.0+xi+eta+xi*eta) );
123  else
124  return RealGradient( 0.0, 0.0, 0.125*(1.0+xi+eta+xi*eta) );
125  }
126  case 7:
127  {
128  if (elem->point(3) > elem->point(7))
129  return RealGradient( 0.0, 0.0, -0.125*(1.0-xi+eta-xi*eta) );
130  else
131  return RealGradient( 0.0, 0.0, 0.125*(1.0-xi+eta-xi*eta) );
132  }
133  case 8:
134  {
135  if (elem->point(4) > elem->point(5))
136  return RealGradient( -0.125*(1.0-eta+zeta-eta*zeta), 0.0, 0.0 );
137  else
138  return RealGradient( 0.125*(1.0-eta+zeta-eta*zeta), 0.0, 0.0 );
139  }
140  case 9:
141  {
142  if (elem->point(5) > elem->point(6))
143  return RealGradient( 0.0, -0.125*(1.0+xi+zeta+xi*zeta), 0.0 );
144  else
145  return RealGradient( 0.0, 0.125*(1.0+xi+zeta+xi*zeta), 0.0 );
146  }
147  case 10:
148  {
149  if (elem->point(7) > elem->point(6))
150  return RealGradient( -0.125*(1.0+eta+zeta+eta*zeta), 0.0, 0.0 );
151  else
152  return RealGradient( 0.125*(1.0+eta+zeta+eta*zeta), 0.0, 0.0 );
153  }
154  case 11:
155  {
156  if (elem->point(4) > elem->point(7))
157  return RealGradient( 0.0, -0.125*(1.0-xi+zeta-xi*zeta), 0.0 );
158  else
159  return RealGradient( 0.0, 0.125*(1.0-xi+zeta-xi*zeta), 0.0 );
160  }
161 
162  default:
163  libmesh_error_msg("Invalid i = " << i);
164  }
165 
166  return RealGradient();
167  }
168 
169  case TET10:
170  {
171  libmesh_assert_less (i, 6);
172 
173  libmesh_not_implemented();
174  return RealGradient();
175  }
176 
177  default:
178  libmesh_error_msg("ERROR: Unsupported 3D element type!: " << elem->type());
179  }
180  }
181 
182  // unsupported order
183  default:
184  libmesh_error_msg("ERROR: Unsupported 3D FE order!: " << totalorder);
185  }
186 #else // LIBMESH_DIM != 3
187  libmesh_ignore(elem, order, i, p, add_p_level);
188  libmesh_not_implemented();
189 #endif
190 }
191 
192 
193 
194 
195 template <>
197  const Order,
198  const unsigned int,
199  const unsigned int,
200  const Point &)
201 {
202  libmesh_error_msg("Nedelec elements require the element type \nbecause edge orientation is needed.");
203  return RealGradient();
204 }
205 
206 template <>
208  const Order order,
209  const unsigned int i,
210  const unsigned int j,
211  const Point & p,
212  const bool add_p_level)
213 {
214 #if LIBMESH_DIM == 3
215  libmesh_assert(elem);
216  libmesh_assert_less (j, 3);
217 
218  const Order totalorder = static_cast<Order>(order + add_p_level * elem->p_level());
219 
220  switch (totalorder)
221  {
222  case FIRST:
223  {
224  switch (elem->type())
225  {
226  case HEX20:
227  case HEX27:
228  {
229  libmesh_assert_less (i, 12);
230 
231  const Real xi = p(0);
232  const Real eta = p(1);
233  const Real zeta = p(2);
234 
235  // Even with a loose inverse_map tolerance we ought to
236  // be nearly on the element interior in master
237  // coordinates
238  libmesh_assert_less_equal ( std::fabs(xi), 1.0+TOLERANCE );
239  libmesh_assert_less_equal ( std::fabs(eta), 1.0+TOLERANCE );
240  libmesh_assert_less_equal ( std::fabs(zeta), 1.0+TOLERANCE );
241 
242  switch (j)
243  {
244  // d()/dxi
245  case 0:
246  {
247  switch(i)
248  {
249  case 0:
250  case 2:
251  case 8:
252  case 10:
253  return RealGradient();
254  case 1:
255  {
256  if (elem->point(1) > elem->point(2))
257  return RealGradient( 0.0, -0.125*(1.0-zeta) );
258  else
259  return RealGradient( 0.0, 0.125*(1.0-zeta) );
260  }
261  case 3:
262  {
263  if (elem->point(3) > elem->point(0))
264  return RealGradient( 0.0, 0.125*(-1.0+zeta) );
265  else
266  return RealGradient( 0.0, -0.125*(-1.0+zeta) );
267  }
268  case 4:
269  {
270  if (elem->point(0) > elem->point(4))
271  return RealGradient( 0.0, 0.0, -0.125*(-1.0+eta) );
272  else
273  return RealGradient( 0.0, 0.0, 0.125*(-1.0+eta) );
274  }
275  case 5:
276  {
277  if (elem->point(1) > elem->point(5))
278  return RealGradient( 0.0, 0.0, -0.125*(1.0-eta) );
279  else
280  return RealGradient( 0.0, 0.0, 0.125*(1.0-eta) );
281  }
282  case 6:
283  {
284  if (elem->point(2) > elem->point(6))
285  return RealGradient( 0.0, 0.0, -0.125*(1.0+eta) );
286  else
287  return RealGradient( 0.0, 0.0, 0.125*(1.0+eta) );
288  }
289  case 7:
290  {
291  if (elem->point(3) > elem->point(7))
292  return RealGradient( 0.0, 0.0, -0.125*(-1.0-eta) );
293  else
294  return RealGradient( 0.0, 0.0, 0.125*(-1.0-eta) );
295  }
296  case 9:
297  {
298  if (elem->point(5) > elem->point(6))
299  return RealGradient( 0.0, -0.125*(1.0+zeta), 0.0 );
300  else
301  return RealGradient( 0.0, 0.125*(1.0+zeta), 0.0 );
302  }
303  case 11:
304  {
305  if (elem->point(4) > elem->point(7))
306  return RealGradient( 0.0, -0.125*(-1.0-zeta), 0.0 );
307  else
308  return RealGradient( 0.0, 0.125*(-1.0-zeta), 0.0 );
309  }
310  default:
311  libmesh_error_msg("Invalid i = " << i);
312  } // switch(i)
313 
314  } // j=0
315 
316  // d()/deta
317  case 1:
318  {
319  switch(i)
320  {
321  case 1:
322  case 3:
323  case 9:
324  case 11:
325  return RealGradient();
326  case 0:
327  {
328  if (elem->point(0) > elem->point(1))
329  return RealGradient( -0.125*(-1.0+zeta), 0.0, 0.0 );
330  else
331  return RealGradient( 0.125*(-1.0+zeta), 0.0, 0.0 );
332  }
333  case 2:
334  {
335  if (elem->point(2) > elem->point(3))
336  return RealGradient( 0.125*(1.0-zeta), 0.0, 0.0 );
337  else
338  return RealGradient( -0.125*(1.0-zeta), 0.0, 0.0 );
339  }
340  case 4:
341  {
342  if (elem->point(0) > elem->point(4))
343  return RealGradient( 0.0, 0.0, -0.125*(-1.0+xi) );
344  else
345  return RealGradient( 0.0, 0.0, 0.125*(-1.0+xi) );
346  }
347  case 5:
348  {
349  if (elem->point(1) > elem->point(5))
350  return RealGradient( 0.0, 0.0, -0.125*(-1.0-xi) );
351  else
352  return RealGradient( 0.0, 0.0, 0.125*(-1.0-xi) );
353  }
354  case 6:
355  {
356  if (elem->point(2) > elem->point(6))
357  return RealGradient( 0.0, 0.0, -0.125*(1.0+xi) );
358  else
359  return RealGradient( 0.0, 0.0, 0.125*(1.0+xi) );
360  }
361  case 7:
362  {
363  if (elem->point(3) > elem->point(7))
364  return RealGradient( 0.0, 0.0, -0.125*(1.0-xi) );
365  else
366  return RealGradient( 0.0, 0.0, 0.125*(1.0-xi) );
367  }
368  case 8:
369  {
370  if (elem->point(4) > elem->point(5))
371  return RealGradient( -0.125*(-1.0-zeta), 0.0, 0.0 );
372  else
373  return RealGradient( 0.125*(-1.0-zeta), 0.0, 0.0 );
374  }
375  case 10:
376  {
377  if (elem->point(7) > elem->point(6))
378  return RealGradient( -0.125*(1.0+zeta), 0.0, 0.0 );
379  else
380  return RealGradient( 0.125*(1.0+zeta), 0.0, 0.0 );
381  }
382  default:
383  libmesh_error_msg("Invalid i = " << i);
384  } // switch(i)
385 
386  } // j=1
387 
388  // d()/dzeta
389  case 2:
390  {
391  switch(i)
392  {
393  case 4:
394  case 5:
395  case 6:
396  case 7:
397  return RealGradient();
398 
399  case 0:
400  {
401  if (elem->point(0) > elem->point(1))
402  return RealGradient( -0.125*(-1.0+eta), 0.0, 0.0 );
403  else
404  return RealGradient( 0.125*(-1.0+eta), 0.0, 0.0 );
405  }
406  case 1:
407  {
408  if (elem->point(1) > elem->point(2))
409  return RealGradient( 0.0, -0.125*(-1.0-xi), 0.0 );
410  else
411  return RealGradient( 0.0, 0.125*(-1.0-xi), 0.0 );
412  }
413  case 2:
414  {
415  if (elem->point(2) > elem->point(3))
416  return RealGradient( 0.125*(-1.0-eta), 0.0, 0.0 );
417  else
418  return RealGradient( -0.125*(-1.0-eta), 0.0, 0.0 );
419  }
420  case 3:
421  {
422  if (elem->point(3) > elem->point(0))
423  return RealGradient( 0.0, 0.125*(-1.0+xi), 0.0 );
424  else
425  return RealGradient( 0.0, -0.125*(-1.0+xi), 0.0 );
426  }
427  case 8:
428  {
429  if (elem->point(4) > elem->point(5))
430  return RealGradient( -0.125*(1.0-eta), 0.0, 0.0 );
431  else
432  return RealGradient( 0.125*(1.0-eta), 0.0, 0.0 );
433  }
434  case 9:
435  {
436  if (elem->point(5) > elem->point(6))
437  return RealGradient( 0.0, -0.125*(1.0+xi), 0.0 );
438  else
439  return RealGradient( 0.0, 0.125*(1.0+xi), 0.0 );
440  }
441  case 10:
442  {
443  if (elem->point(7) > elem->point(6))
444  return RealGradient( -0.125*(1.0+eta), 0.0, 0.0 );
445  else
446  return RealGradient( 0.125*(1.0+eta), 0.0, 0.0 );
447  }
448  case 11:
449  {
450  if (elem->point(4) > elem->point(7))
451  return RealGradient( 0.0, -0.125*(1.0-xi), 0.0 );
452  else
453  return RealGradient( 0.0, 0.125*(1.0-xi), 0.0 );
454  }
455  default:
456  libmesh_error_msg("Invalid i = " << i);
457  } // switch(i)
458 
459  } // j = 2
460 
461  default:
462  libmesh_error_msg("Invalid j = " << j);
463  }
464 
465  return RealGradient();
466  }
467 
468  case TET10:
469  {
470  libmesh_assert_less (i, 6);
471 
472  libmesh_not_implemented();
473  return RealGradient();
474  }
475 
476  default:
477  libmesh_error_msg("ERROR: Unsupported 3D element type!: " << elem->type());
478  }
479  }
480 
481  // unsupported order
482  default:
483  libmesh_error_msg("ERROR: Unsupported 3D FE order!: " << totalorder);
484  }
485 
486 #else // LIBMESH_DIM != 3
487  libmesh_ignore(elem, order, i, j, p, add_p_level);
488  libmesh_not_implemented();
489 #endif
490 }
491 
492 
493 
494 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
495 template <>
497  const Order,
498  const unsigned int,
499  const unsigned int,
500  const Point &)
501 {
502  libmesh_error_msg("Nedelec elements require the element type \nbecause edge orientation is needed.");
503  return RealGradient();
504 }
505 
506 
507 
508 template <>
510  const Order order,
511  const unsigned int i,
512  const unsigned int j,
513  const Point & libmesh_dbg_var(p),
514  const bool add_p_level)
515 {
516 #if LIBMESH_DIM == 3
517 
518  libmesh_assert(elem);
519 
520  // j = 0 ==> d^2 phi / dxi^2
521  // j = 1 ==> d^2 phi / dxi deta
522  // j = 2 ==> d^2 phi / deta^2
523  // j = 3 ==> d^2 phi / dxi dzeta
524  // j = 4 ==> d^2 phi / deta dzeta
525  // j = 5 ==> d^2 phi / dzeta^2
526  libmesh_assert_less (j, 6);
527 
528  const Order totalorder = static_cast<Order>(order + add_p_level * elem->p_level());
529 
530  switch (totalorder)
531  {
532  // linear Lagrange shape functions
533  case FIRST:
534  {
535  switch (elem->type())
536  {
537  case HEX20:
538  case HEX27:
539  {
540  libmesh_assert_less (i, 12);
541 
542 #ifndef NDEBUG
543  const Real xi = p(0);
544  const Real eta = p(1);
545  const Real zeta = p(2);
546 #endif
547 
548  libmesh_assert_less_equal ( std::fabs(xi), 1.0+TOLERANCE );
549  libmesh_assert_less_equal ( std::fabs(eta), 1.0+TOLERANCE );
550  libmesh_assert_less_equal ( std::fabs(zeta), 1.0+TOLERANCE );
551 
552  switch (j)
553  {
554  // d^2()/dxi^2
555  case 0:
556  {
557  // All d^2()/dxi^2 derivatives for linear hexes are zero.
558  return RealGradient();
559  } // j=0
560 
561  // d^2()/dxideta
562  case 1:
563  {
564  switch(i)
565  {
566  case 0:
567  case 1:
568  case 2:
569  case 3:
570  case 8:
571  case 9:
572  case 10:
573  case 11:
574  return RealGradient();
575  case 4:
576  {
577  if (elem->point(0) > elem->point(4))
578  return RealGradient( 0.0, 0.0, -0.125 );
579  else
580  return RealGradient( 0.0, 0.0, 0.125 );
581  }
582  case 5:
583  {
584  if (elem->point(1) > elem->point(5))
585  return RealGradient( 0.0, 0.0, 0.125 );
586  else
587  return RealGradient( 0.0, 0.0, -0.125 );
588  }
589  case 6:
590  {
591  if (elem->point(2) > elem->point(6))
592  return RealGradient( 0.0, 0.0, -0.125 );
593  else
594  return RealGradient( 0.0, 0.0, 0.125 );
595  }
596  case 7:
597  {
598  if (elem->point(3) > elem->point(7))
599  return RealGradient( 0.0, 0.0, 0.125 );
600  else
601  return RealGradient( 0.0, 0.0, -0.125 );
602  }
603  default:
604  libmesh_error_msg("Invalid i = " << i);
605  } // switch(i)
606 
607  } // j=1
608 
609  // d^2()/deta^2
610  case 2:
611  {
612  // All d^2()/deta^2 derivatives for linear hexes are zero.
613  return RealGradient();
614  } // j = 2
615 
616  // d^2()/dxidzeta
617  case 3:
618  {
619  switch(i)
620  {
621  case 0:
622  case 2:
623  case 4:
624  case 5:
625  case 6:
626  case 7:
627  case 8:
628  case 10:
629  return RealGradient();
630 
631  case 1:
632  {
633  if (elem->point(1) > elem->point(2))
634  return RealGradient( 0.0, 0.125 );
635  else
636  return RealGradient( 0.0, -0.125 );
637  }
638  case 3:
639  {
640  if (elem->point(3) > elem->point(0))
641  return RealGradient( 0.0, -0.125 );
642  else
643  return RealGradient( 0.0, 0.125 );
644  }
645  case 9:
646  {
647  if (elem->point(5) > elem->point(6))
648  return RealGradient( 0.0, -0.125, 0.0 );
649  else
650  return RealGradient( 0.0, 0.125, 0.0 );
651  }
652  case 11:
653  {
654  if (elem->point(4) > elem->point(7))
655  return RealGradient( 0.0, 0.125, 0.0 );
656  else
657  return RealGradient( 0.0, -0.125, 0.0 );
658  }
659  default:
660  libmesh_error_msg("Invalid i = " << i);
661  } // switch(i)
662 
663  } // j = 3
664 
665  // d^2()/detadzeta
666  case 4:
667  {
668  switch(i)
669  {
670  case 1:
671  case 3:
672  case 4:
673  case 5:
674  case 6:
675  case 7:
676  case 9:
677  case 11:
678  return RealGradient();
679 
680  case 0:
681  {
682  if (elem->point(0) > elem->point(1))
683  return RealGradient( -0.125, 0.0, 0.0 );
684  else
685  return RealGradient( 0.125, 0.0, 0.0 );
686  }
687  case 2:
688  {
689  if (elem->point(2) > elem->point(3))
690  return RealGradient( 0.125, 0.0, 0.0 );
691  else
692  return RealGradient( -0.125, 0.0, 0.0 );
693  }
694  case 8:
695  {
696  if (elem->point(4) > elem->point(5))
697  return RealGradient( 0.125, 0.0, 0.0 );
698  else
699  return RealGradient( -0.125, 0.0, 0.0 );
700  }
701  case 10:
702  {
703  if (elem->point(7) > elem->point(6))
704  return RealGradient( -0.125, 0.0, 0.0 );
705  else
706  return RealGradient( 0.125, 0.0, 0.0 );
707  }
708  default:
709  libmesh_error_msg("Invalid i = " << i);
710  } // switch(i)
711 
712  } // j = 4
713 
714  // d^2()/dzeta^2
715  case 5:
716  {
717  // All d^2()/dzeta^2 derivatives for linear hexes are zero.
718  return RealGradient();
719  } // j = 5
720 
721  default:
722  libmesh_error_msg("Invalid j = " << j);
723  }
724 
725  return RealGradient();
726  }
727 
728  case TET10:
729  {
730  libmesh_assert_less (i, 6);
731 
732  libmesh_not_implemented();
733  return RealGradient();
734  }
735 
736  default:
737  libmesh_error_msg("ERROR: Unsupported 3D element type!: " << elem->type());
738 
739  } //switch(type)
740 
741  } // case FIRST:
742  // unsupported order
743  default:
744  libmesh_error_msg("ERROR: Unsupported 3D FE order!: " << totalorder);
745  }
746 
747 #else // LIBMESH_DIM != 3
748  libmesh_assert(true || p(0));
749  libmesh_ignore(elem, order, i, j, add_p_level);
750  libmesh_not_implemented();
751 #endif
752 }
753 
754 #endif
755 
756 } // namespace libMesh
libMesh::FE::shape_second_deriv
static OutputShape shape_second_deriv(const ElemType t, const Order o, const unsigned int i, const unsigned int j, const Point &p)
libMesh::HEX20
Definition: enum_elem_type.h:48
libMesh::RealGradient
RealVectorValue RealGradient
Definition: hp_coarsentest.h:49
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::TET10
Definition: enum_elem_type.h:46
libMesh::Order
Order
Definition: enum_order.h:40
libMesh::TOLERANCE
static const Real TOLERANCE
Definition: libmesh_common.h:128
libMesh::Elem::p_level
unsigned int p_level() const
Definition: elem.h:2512
libMesh::Elem::point
const Point & point(const unsigned int i) const
Definition: elem.h:1955
libMesh::VectorValue< Real >
libMesh::libmesh_assert
libmesh_assert(ctx)
libMesh::FE::shape
static OutputShape shape(const ElemType t, const Order o, const unsigned int i, const Point &p)
libMesh::HEX27
Definition: enum_elem_type.h:49
libMesh::FE::shape_deriv
static OutputShape shape_deriv(const ElemType t, const Order o, const unsigned int i, const unsigned int j, const Point &p)
libMesh::libmesh_ignore
void libmesh_ignore(const Args &...)
Definition: libmesh_common.h:526
libMesh::Point
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:38
libMesh::Elem
This is the base class from which all geometric element types are derived.
Definition: elem.h:100
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::Elem::type
virtual ElemType type() const =0
libMesh::FIRST
Definition: enum_order.h:42
libMesh::ElemType
ElemType
Defines an enum for geometric element types.
Definition: enum_elem_type.h:33