fe_lagrange_shape_2D.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2018 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 
29 
30 
31 template <>
33  const Order order,
34  const unsigned int i,
35  const Point & p)
36 {
37 #if LIBMESH_DIM > 1
38 
39  switch (order)
40  {
41  // linear Lagrange shape functions
42  case FIRST:
43  {
44  switch (type)
45  {
46  case QUAD4:
47  case QUADSHELL4:
48  case QUAD8:
49  case QUADSHELL8:
50  case QUAD9:
51  {
52  // Compute quad shape functions as a tensor-product
53  const Real xi = p(0);
54  const Real eta = p(1);
55 
56  libmesh_assert_less (i, 4);
57 
58  // 0 1 2 3
59  static const unsigned int i0[] = {0, 1, 1, 0};
60  static const unsigned int i1[] = {0, 0, 1, 1};
61 
62  return (FE<1,LAGRANGE>::shape(EDGE2, FIRST, i0[i], xi)*
63  FE<1,LAGRANGE>::shape(EDGE2, FIRST, i1[i], eta));
64  }
65 
66  case TRI3:
67  case TRISHELL3:
68  case TRI6:
69  {
70  const Real zeta1 = p(0);
71  const Real zeta2 = p(1);
72  const Real zeta0 = 1. - zeta1 - zeta2;
73 
74  libmesh_assert_less (i, 3);
75 
76  switch(i)
77  {
78  case 0:
79  return zeta0;
80 
81  case 1:
82  return zeta1;
83 
84  case 2:
85  return zeta2;
86 
87  default:
88  libmesh_error_msg("Invalid shape function index i = " << i);
89  }
90  }
91 
92  default:
93  libmesh_error_msg("ERROR: Unsupported 2D element type: " << type);
94  }
95  }
96 
97 
98  // quadratic Lagrange shape functions
99  case SECOND:
100  {
101  switch (type)
102  {
103  case QUAD8:
104  case QUADSHELL8:
105  {
106  const Real xi = p(0);
107  const Real eta = p(1);
108 
109  libmesh_assert_less (i, 8);
110 
111  switch (i)
112  {
113  case 0:
114  return .25*(1. - xi)*(1. - eta)*(-1. - xi - eta);
115 
116  case 1:
117  return .25*(1. + xi)*(1. - eta)*(-1. + xi - eta);
118 
119  case 2:
120  return .25*(1. + xi)*(1. + eta)*(-1. + xi + eta);
121 
122  case 3:
123  return .25*(1. - xi)*(1. + eta)*(-1. - xi + eta);
124 
125  case 4:
126  return .5*(1. - xi*xi)*(1. - eta);
127 
128  case 5:
129  return .5*(1. + xi)*(1. - eta*eta);
130 
131  case 6:
132  return .5*(1. - xi*xi)*(1. + eta);
133 
134  case 7:
135  return .5*(1. - xi)*(1. - eta*eta);
136 
137  default:
138  libmesh_error_msg("Invalid shape function index i = " << i);
139  }
140  }
141 
142  case QUAD9:
143  {
144  // Compute quad shape functions as a tensor-product
145  const Real xi = p(0);
146  const Real eta = p(1);
147 
148  libmesh_assert_less (i, 9);
149 
150  // 0 1 2 3 4 5 6 7 8
151  static const unsigned int i0[] = {0, 1, 1, 0, 2, 1, 2, 0, 2};
152  static const unsigned int i1[] = {0, 0, 1, 1, 0, 2, 1, 2, 2};
153 
154  return (FE<1,LAGRANGE>::shape(EDGE3, SECOND, i0[i], xi)*
155  FE<1,LAGRANGE>::shape(EDGE3, SECOND, i1[i], eta));
156  }
157 
158  case TRI6:
159  {
160  const Real zeta1 = p(0);
161  const Real zeta2 = p(1);
162  const Real zeta0 = 1. - zeta1 - zeta2;
163 
164  libmesh_assert_less (i, 6);
165 
166  switch(i)
167  {
168  case 0:
169  return 2.*zeta0*(zeta0-0.5);
170 
171  case 1:
172  return 2.*zeta1*(zeta1-0.5);
173 
174  case 2:
175  return 2.*zeta2*(zeta2-0.5);
176 
177  case 3:
178  return 4.*zeta0*zeta1;
179 
180  case 4:
181  return 4.*zeta1*zeta2;
182 
183  case 5:
184  return 4.*zeta2*zeta0;
185 
186  default:
187  libmesh_error_msg("Invalid shape function index i = " << i);
188  }
189  }
190 
191  default:
192  libmesh_error_msg("ERROR: Unsupported 2D element type: " << type);
193  }
194  }
195 
196 
197 
198  // unsupported order
199  default:
200  libmesh_error_msg("ERROR: Unsupported 2D FE order: " << order);
201  }
202 #else // LIBMESH_DIM > 1
203  return 0.;
204 #endif
205 }
206 
207 
208 
209 template <>
211  const Order order,
212  const unsigned int i,
213  const Point & p)
214 {
215  libmesh_assert(elem);
216 
217  // call the orientation-independent shape functions
218  return FE<2,LAGRANGE>::shape(elem->type(), static_cast<Order>(order + elem->p_level()), i, p);
219 }
220 
221 
222 
223 template <>
225  const Order order,
226  const unsigned int i,
227  const unsigned int j,
228  const Point & p)
229 {
230 #if LIBMESH_DIM > 1
231 
232 
233  libmesh_assert_less (j, 2);
234 
235  switch (order)
236  {
237  // linear Lagrange shape functions
238  case FIRST:
239  {
240  switch (type)
241  {
242  case QUAD4:
243  case QUADSHELL4:
244  case QUAD8:
245  case QUADSHELL8:
246  case QUAD9:
247  {
248  // Compute quad shape functions as a tensor-product
249  const Real xi = p(0);
250  const Real eta = p(1);
251 
252  libmesh_assert_less (i, 4);
253 
254  // 0 1 2 3
255  static const unsigned int i0[] = {0, 1, 1, 0};
256  static const unsigned int i1[] = {0, 0, 1, 1};
257 
258  switch (j)
259  {
260  // d()/dxi
261  case 0:
262  return (FE<1,LAGRANGE>::shape_deriv(EDGE2, FIRST, i0[i], 0, xi)*
263  FE<1,LAGRANGE>::shape (EDGE2, FIRST, i1[i], eta));
264 
265  // d()/deta
266  case 1:
267  return (FE<1,LAGRANGE>::shape (EDGE2, FIRST, i0[i], xi)*
268  FE<1,LAGRANGE>::shape_deriv(EDGE2, FIRST, i1[i], 0, eta));
269 
270  default:
271  libmesh_error_msg("ERROR: Invalid derivative index j = " << j);
272  }
273  }
274 
275  case TRI3:
276  case TRISHELL3:
277  case TRI6:
278  {
279  libmesh_assert_less (i, 3);
280 
281  const Real dzeta0dxi = -1.;
282  const Real dzeta1dxi = 1.;
283  const Real dzeta2dxi = 0.;
284 
285  const Real dzeta0deta = -1.;
286  const Real dzeta1deta = 0.;
287  const Real dzeta2deta = 1.;
288 
289  switch (j)
290  {
291  // d()/dxi
292  case 0:
293  {
294  switch(i)
295  {
296  case 0:
297  return dzeta0dxi;
298 
299  case 1:
300  return dzeta1dxi;
301 
302  case 2:
303  return dzeta2dxi;
304 
305  default:
306  libmesh_error_msg("Invalid shape function index i = " << i);
307  }
308  }
309  // d()/deta
310  case 1:
311  {
312  switch(i)
313  {
314  case 0:
315  return dzeta0deta;
316 
317  case 1:
318  return dzeta1deta;
319 
320  case 2:
321  return dzeta2deta;
322 
323  default:
324  libmesh_error_msg("Invalid shape function index i = " << i);
325  }
326  }
327  default:
328  libmesh_error_msg("ERROR: Invalid derivative index j = " << j);
329  }
330  }
331 
332  default:
333  libmesh_error_msg("ERROR: Unsupported 2D element type: " << type);
334  }
335  }
336 
337 
338  // quadratic Lagrange shape functions
339  case SECOND:
340  {
341  switch (type)
342  {
343  case QUAD8:
344  case QUADSHELL8:
345  {
346  const Real xi = p(0);
347  const Real eta = p(1);
348 
349  libmesh_assert_less (i, 8);
350 
351  switch (j)
352  {
353  // d/dxi
354  case 0:
355  switch (i)
356  {
357  case 0:
358  return .25*(1. - eta)*((1. - xi)*(-1.) +
359  (-1.)*(-1. - xi - eta));
360 
361  case 1:
362  return .25*(1. - eta)*((1. + xi)*(1.) +
363  (1.)*(-1. + xi - eta));
364 
365  case 2:
366  return .25*(1. + eta)*((1. + xi)*(1.) +
367  (1.)*(-1. + xi + eta));
368 
369  case 3:
370  return .25*(1. + eta)*((1. - xi)*(-1.) +
371  (-1.)*(-1. - xi + eta));
372 
373  case 4:
374  return .5*(-2.*xi)*(1. - eta);
375 
376  case 5:
377  return .5*(1.)*(1. - eta*eta);
378 
379  case 6:
380  return .5*(-2.*xi)*(1. + eta);
381 
382  case 7:
383  return .5*(-1.)*(1. - eta*eta);
384 
385  default:
386  libmesh_error_msg("Invalid shape function index i = " << i);
387  }
388 
389  // d/deta
390  case 1:
391  switch (i)
392  {
393  case 0:
394  return .25*(1. - xi)*((1. - eta)*(-1.) +
395  (-1.)*(-1. - xi - eta));
396 
397  case 1:
398  return .25*(1. + xi)*((1. - eta)*(-1.) +
399  (-1.)*(-1. + xi - eta));
400 
401  case 2:
402  return .25*(1. + xi)*((1. + eta)*(1.) +
403  (1.)*(-1. + xi + eta));
404 
405  case 3:
406  return .25*(1. - xi)*((1. + eta)*(1.) +
407  (1.)*(-1. - xi + eta));
408 
409  case 4:
410  return .5*(1. - xi*xi)*(-1.);
411 
412  case 5:
413  return .5*(1. + xi)*(-2.*eta);
414 
415  case 6:
416  return .5*(1. - xi*xi)*(1.);
417 
418  case 7:
419  return .5*(1. - xi)*(-2.*eta);
420 
421  default:
422  libmesh_error_msg("Invalid shape function index i = " << i);
423  }
424 
425  default:
426  libmesh_error_msg("ERROR: Invalid derivative index j = " << j);
427  }
428  }
429 
430  case QUAD9:
431  {
432  // Compute quad shape functions as a tensor-product
433  const Real xi = p(0);
434  const Real eta = p(1);
435 
436  libmesh_assert_less (i, 9);
437 
438  // 0 1 2 3 4 5 6 7 8
439  static const unsigned int i0[] = {0, 1, 1, 0, 2, 1, 2, 0, 2};
440  static const unsigned int i1[] = {0, 0, 1, 1, 0, 2, 1, 2, 2};
441 
442  switch (j)
443  {
444  // d()/dxi
445  case 0:
446  return (FE<1,LAGRANGE>::shape_deriv(EDGE3, SECOND, i0[i], 0, xi)*
447  FE<1,LAGRANGE>::shape (EDGE3, SECOND, i1[i], eta));
448 
449  // d()/deta
450  case 1:
451  return (FE<1,LAGRANGE>::shape (EDGE3, SECOND, i0[i], xi)*
452  FE<1,LAGRANGE>::shape_deriv(EDGE3, SECOND, i1[i], 0, eta));
453 
454  default:
455  libmesh_error_msg("ERROR: Invalid derivative index j = " << j);
456  }
457  }
458 
459  case TRI6:
460  {
461  libmesh_assert_less (i, 6);
462 
463  const Real zeta1 = p(0);
464  const Real zeta2 = p(1);
465  const Real zeta0 = 1. - zeta1 - zeta2;
466 
467  const Real dzeta0dxi = -1.;
468  const Real dzeta1dxi = 1.;
469  const Real dzeta2dxi = 0.;
470 
471  const Real dzeta0deta = -1.;
472  const Real dzeta1deta = 0.;
473  const Real dzeta2deta = 1.;
474 
475  switch(j)
476  {
477  case 0:
478  {
479  switch(i)
480  {
481  case 0:
482  return (4.*zeta0-1.)*dzeta0dxi;
483 
484  case 1:
485  return (4.*zeta1-1.)*dzeta1dxi;
486 
487  case 2:
488  return (4.*zeta2-1.)*dzeta2dxi;
489 
490  case 3:
491  return 4.*zeta1*dzeta0dxi + 4.*zeta0*dzeta1dxi;
492 
493  case 4:
494  return 4.*zeta2*dzeta1dxi + 4.*zeta1*dzeta2dxi;
495 
496  case 5:
497  return 4.*zeta2*dzeta0dxi + 4*zeta0*dzeta2dxi;
498 
499  default:
500  libmesh_error_msg("Invalid shape function index i = " << i);
501  }
502  }
503 
504  case 1:
505  {
506  switch(i)
507  {
508  case 0:
509  return (4.*zeta0-1.)*dzeta0deta;
510 
511  case 1:
512  return (4.*zeta1-1.)*dzeta1deta;
513 
514  case 2:
515  return (4.*zeta2-1.)*dzeta2deta;
516 
517  case 3:
518  return 4.*zeta1*dzeta0deta + 4.*zeta0*dzeta1deta;
519 
520  case 4:
521  return 4.*zeta2*dzeta1deta + 4.*zeta1*dzeta2deta;
522 
523  case 5:
524  return 4.*zeta2*dzeta0deta + 4*zeta0*dzeta2deta;
525 
526  default:
527  libmesh_error_msg("Invalid shape function index i = " << i);
528  }
529  }
530  default:
531  libmesh_error_msg("ERROR: Invalid derivative index j = " << j);
532  }
533  }
534 
535  default:
536  libmesh_error_msg("ERROR: Unsupported 2D element type: " << type);
537  }
538  }
539 
540  // unsupported order
541  default:
542  libmesh_error_msg("ERROR: Unsupported 2D FE order: " << order);
543  }
544 #else // LIBMESH_DIM > 1
545  return 0.;
546 #endif
547 }
548 
549 
550 
551 template <>
553  const Order order,
554  const unsigned int i,
555  const unsigned int j,
556  const Point & p)
557 {
558  libmesh_assert(elem);
559 
560 
561  // call the orientation-independent shape functions
562  return FE<2,LAGRANGE>::shape_deriv(elem->type(), static_cast<Order>(order + elem->p_level()), i, j, p);
563 }
564 
565 
566 
567 
568 template <>
570  const Order order,
571  const unsigned int i,
572  const unsigned int j,
573  const Point & p)
574 {
575 #if LIBMESH_DIM > 1
576 
577  // j = 0 ==> d^2 phi / dxi^2
578  // j = 1 ==> d^2 phi / dxi deta
579  // j = 2 ==> d^2 phi / deta^2
580  libmesh_assert_less (j, 3);
581 
582  switch (order)
583  {
584  // linear Lagrange shape functions
585  case FIRST:
586  {
587  switch (type)
588  {
589  case QUAD4:
590  case QUADSHELL4:
591  case QUAD8:
592  case QUADSHELL8:
593  case QUAD9:
594  {
595  // Compute quad shape functions as a tensor-product
596  const Real xi = p(0);
597  const Real eta = p(1);
598 
599  libmesh_assert_less (i, 4);
600 
601  // 0 1 2 3
602  static const unsigned int i0[] = {0, 1, 1, 0};
603  static const unsigned int i1[] = {0, 0, 1, 1};
604 
605  switch (j)
606  {
607  // d^2() / dxi^2
608  case 0:
609  return 0.;
610 
611  // d^2() / dxi deta
612  case 1:
613  return (FE<1,LAGRANGE>::shape_deriv(EDGE2, FIRST, i0[i], 0, xi)*
614  FE<1,LAGRANGE>::shape_deriv(EDGE2, FIRST, i1[i], 0, eta));
615 
616  // d^2() / deta^2
617  case 2:
618  return 0.;
619 
620  default:
621  libmesh_error_msg("ERROR: Invalid derivative index j = " << j);
622  }
623  }
624 
625  case TRI3:
626  case TRISHELL3:
627  case TRI6:
628  {
629  // All second derivatives for linear triangles are zero.
630  return 0.;
631  }
632 
633  default:
634  libmesh_error_msg("ERROR: Unsupported 2D element type: " << type);
635 
636  } // end switch (type)
637  } // end case FIRST
638 
639 
640  // quadratic Lagrange shape functions
641  case SECOND:
642  {
643  switch (type)
644  {
645  case QUAD8:
646  case QUADSHELL8:
647  {
648  const Real xi = p(0);
649  const Real eta = p(1);
650 
651  libmesh_assert_less (j, 3);
652 
653  switch (j)
654  {
655  // d^2() / dxi^2
656  case 0:
657  {
658  switch (i)
659  {
660  case 0:
661  case 1:
662  return 0.5*(1.-eta);
663 
664  case 2:
665  case 3:
666  return 0.5*(1.+eta);
667 
668  case 4:
669  return eta - 1.;
670 
671  case 5:
672  case 7:
673  return 0.0;
674 
675  case 6:
676  return -1. - eta;
677 
678  default:
679  libmesh_error_msg("Invalid shape function index i = " << i);
680  }
681  }
682 
683  // d^2() / dxi deta
684  case 1:
685  {
686  switch (i)
687  {
688  case 0:
689  return 0.25*( 1. - 2.*xi - 2.*eta);
690 
691  case 1:
692  return 0.25*(-1. - 2.*xi + 2.*eta);
693 
694  case 2:
695  return 0.25*( 1. + 2.*xi + 2.*eta);
696 
697  case 3:
698  return 0.25*(-1. + 2.*xi - 2.*eta);
699 
700  case 4:
701  return xi;
702 
703  case 5:
704  return -eta;
705 
706  case 6:
707  return -xi;
708 
709  case 7:
710  return eta;
711 
712  default:
713  libmesh_error_msg("Invalid shape function index i = " << i);
714  }
715  }
716 
717  // d^2() / deta^2
718  case 2:
719  {
720  switch (i)
721  {
722  case 0:
723  case 3:
724  return 0.5*(1.-xi);
725 
726  case 1:
727  case 2:
728  return 0.5*(1.+xi);
729 
730  case 4:
731  case 6:
732  return 0.0;
733 
734  case 5:
735  return -1.0 - xi;
736 
737  case 7:
738  return xi - 1.0;
739 
740  default:
741  libmesh_error_msg("Invalid shape function index i = " << i);
742  }
743  }
744 
745  default:
746  libmesh_error_msg("ERROR: Invalid derivative index j = " << j);
747  } // end switch (j)
748  } // end case QUAD8
749 
750  case QUAD9:
751  {
752  // Compute QUAD9 second derivatives as tensor product
753  const Real xi = p(0);
754  const Real eta = p(1);
755 
756  libmesh_assert_less (i, 9);
757 
758  // 0 1 2 3 4 5 6 7 8
759  static const unsigned int i0[] = {0, 1, 1, 0, 2, 1, 2, 0, 2};
760  static const unsigned int i1[] = {0, 0, 1, 1, 0, 2, 1, 2, 2};
761 
762  switch (j)
763  {
764  // d^2() / dxi^2
765  case 0:
766  return (FE<1,LAGRANGE>::shape_second_deriv(EDGE3, SECOND, i0[i], 0, xi)*
767  FE<1,LAGRANGE>::shape (EDGE3, SECOND, i1[i], eta));
768 
769  // d^2() / dxi deta
770  case 1:
771  return (FE<1,LAGRANGE>::shape_deriv(EDGE3, SECOND, i0[i], 0, xi)*
772  FE<1,LAGRANGE>::shape_deriv(EDGE3, SECOND, i1[i], 0, eta));
773 
774  // d^2() / deta^2
775  case 2:
776  return (FE<1,LAGRANGE>::shape (EDGE3, SECOND, i0[i], xi)*
778 
779  default:
780  libmesh_error_msg("ERROR: Invalid derivative index j = " << j);
781  } // end switch (j)
782  } // end case QUAD9
783 
784  case TRI6:
785  {
786  const Real dzeta0dxi = -1.;
787  const Real dzeta1dxi = 1.;
788  const Real dzeta2dxi = 0.;
789 
790  const Real dzeta0deta = -1.;
791  const Real dzeta1deta = 0.;
792  const Real dzeta2deta = 1.;
793 
794  libmesh_assert_less (j, 3);
795 
796  switch (j)
797  {
798  // d^2() / dxi^2
799  case 0:
800  {
801  switch (i)
802  {
803  case 0:
804  return 4.*dzeta0dxi*dzeta0dxi;
805 
806  case 1:
807  return 4.*dzeta1dxi*dzeta1dxi;
808 
809  case 2:
810  return 4.*dzeta2dxi*dzeta2dxi;
811 
812  case 3:
813  return 8.*dzeta0dxi*dzeta1dxi;
814 
815  case 4:
816  return 8.*dzeta1dxi*dzeta2dxi;
817 
818  case 5:
819  return 8.*dzeta0dxi*dzeta2dxi;
820 
821  default:
822  libmesh_error_msg("Invalid shape function index i = " << i);
823  }
824  }
825 
826  // d^2() / dxi deta
827  case 1:
828  {
829  switch (i)
830  {
831  case 0:
832  return 4.*dzeta0dxi*dzeta0deta;
833 
834  case 1:
835  return 4.*dzeta1dxi*dzeta1deta;
836 
837  case 2:
838  return 4.*dzeta2dxi*dzeta2deta;
839 
840  case 3:
841  return 4.*dzeta1deta*dzeta0dxi + 4.*dzeta0deta*dzeta1dxi;
842 
843  case 4:
844  return 4.*dzeta2deta*dzeta1dxi + 4.*dzeta1deta*dzeta2dxi;
845 
846  case 5:
847  return 4.*dzeta2deta*dzeta0dxi + 4.*dzeta0deta*dzeta2dxi;
848 
849  default:
850  libmesh_error_msg("Invalid shape function index i = " << i);
851  }
852  }
853 
854  // d^2() / deta^2
855  case 2:
856  {
857  switch (i)
858  {
859  case 0:
860  return 4.*dzeta0deta*dzeta0deta;
861 
862  case 1:
863  return 4.*dzeta1deta*dzeta1deta;
864 
865  case 2:
866  return 4.*dzeta2deta*dzeta2deta;
867 
868  case 3:
869  return 8.*dzeta0deta*dzeta1deta;
870 
871  case 4:
872  return 8.*dzeta1deta*dzeta2deta;
873 
874  case 5:
875  return 8.*dzeta0deta*dzeta2deta;
876 
877  default:
878  libmesh_error_msg("Invalid shape function index i = " << i);
879  }
880  }
881 
882  default:
883  libmesh_error_msg("ERROR: Invalid derivative index j = " << j);
884  } // end switch (j)
885  } // end case TRI6
886 
887  default:
888  libmesh_error_msg("ERROR: Unsupported 2D element type: " << type);
889  }
890  } // end case SECOND
891 
892 
893 
894  // unsupported order
895  default:
896  libmesh_error_msg("ERROR: Unsupported 2D FE order: " << order);
897 
898  } // end switch (order)
899 
900 #else // LIBMESH_DIM > 1
901  return 0.;
902 #endif
903 }
904 
905 
906 
907 template <>
909  const Order order,
910  const unsigned int i,
911  const unsigned int j,
912  const Point & p)
913 {
914  libmesh_assert(elem);
915 
916  // call the orientation-independent shape functions
917  return FE<2,LAGRANGE>::shape_second_deriv(elem->type(), static_cast<Order>(order + elem->p_level()), i, j, p);
918 }
919 
920 } // namespace libMesh
static OutputShape shape(const ElemType t, const Order o, const unsigned int i, const Point &p)
The base class for all geometric element types.
Definition: elem.h:100
static OutputShape shape_deriv(const ElemType t, const Order o, const unsigned int i, const unsigned int j, const Point &p)
unsigned int p_level() const
Definition: elem.h:2555
Template class which generates the different FE families and orders.
Definition: fe.h:89
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual ElemType type() const =0
A geometric point in (x,y,z) space.
Definition: point.h:38
static OutputShape shape_second_deriv(const ElemType t, const Order o, const unsigned int i, const unsigned int j, const Point &p)