fe_l2_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,L2_LAGRANGE>::shape(EDGE2, FIRST, i0[i], xi)*
63  FE<1,L2_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,L2_LAGRANGE>::shape(EDGE3, SECOND, i0[i], xi)*
155  FE<1,L2_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 #endif
203 }
204 
205 
206 
207 template <>
209  const Order order,
210  const unsigned int i,
211  const Point & p)
212 {
213  libmesh_assert(elem);
214 
215  // call the orientation-independent shape functions
216  return FE<2,L2_LAGRANGE>::shape(elem->type(), static_cast<Order>(order + elem->p_level()), i, p);
217 }
218 
219 
220 
221 template <>
223  const Order order,
224  const unsigned int i,
225  const unsigned int j,
226  const Point & p)
227 {
228 #if LIBMESH_DIM > 1
229 
230 
231  libmesh_assert_less (j, 2);
232 
233  switch (order)
234  {
235  // linear Lagrange shape functions
236  case FIRST:
237  {
238  switch (type)
239  {
240  case QUAD4:
241  case QUADSHELL4:
242  case QUAD8:
243  case QUADSHELL8:
244  case QUAD9:
245  {
246  // Compute quad shape functions as a tensor-product
247  const Real xi = p(0);
248  const Real eta = p(1);
249 
250  libmesh_assert_less (i, 4);
251 
252  // 0 1 2 3
253  static const unsigned int i0[] = {0, 1, 1, 0};
254  static const unsigned int i1[] = {0, 0, 1, 1};
255 
256  switch (j)
257  {
258  // d()/dxi
259  case 0:
260  return (FE<1,L2_LAGRANGE>::shape_deriv(EDGE2, FIRST, i0[i], 0, xi)*
261  FE<1,L2_LAGRANGE>::shape (EDGE2, FIRST, i1[i], eta));
262 
263  // d()/deta
264  case 1:
265  return (FE<1,L2_LAGRANGE>::shape (EDGE2, FIRST, i0[i], xi)*
266  FE<1,L2_LAGRANGE>::shape_deriv(EDGE2, FIRST, i1[i], 0, eta));
267 
268  default:
269  libmesh_error_msg("ERROR: Invalid derivative index j = " << j);
270  }
271  }
272 
273  case TRI3:
274  case TRISHELL3:
275  case TRI6:
276  {
277  libmesh_assert_less (i, 3);
278 
279  const Real dzeta0dxi = -1.;
280  const Real dzeta1dxi = 1.;
281  const Real dzeta2dxi = 0.;
282 
283  const Real dzeta0deta = -1.;
284  const Real dzeta1deta = 0.;
285  const Real dzeta2deta = 1.;
286 
287  switch (j)
288  {
289  // d()/dxi
290  case 0:
291  {
292  switch(i)
293  {
294  case 0:
295  return dzeta0dxi;
296 
297  case 1:
298  return dzeta1dxi;
299 
300  case 2:
301  return dzeta2dxi;
302 
303  default:
304  libmesh_error_msg("Invalid shape function index i = " << i);
305  }
306  }
307  // d()/deta
308  case 1:
309  {
310  switch(i)
311  {
312  case 0:
313  return dzeta0deta;
314 
315  case 1:
316  return dzeta1deta;
317 
318  case 2:
319  return dzeta2deta;
320 
321  default:
322  libmesh_error_msg("Invalid shape function index i = " << i);
323  }
324  }
325  default:
326  libmesh_error_msg("ERROR: Invalid derivative index j = " << j);
327  }
328  }
329 
330  default:
331  libmesh_error_msg("ERROR: Unsupported 2D element type: " << type);
332  }
333  }
334 
335 
336  // quadratic Lagrange shape functions
337  case SECOND:
338  {
339  switch (type)
340  {
341  case QUAD8:
342  case QUADSHELL8:
343  {
344  const Real xi = p(0);
345  const Real eta = p(1);
346 
347  libmesh_assert_less (i, 8);
348 
349  switch (j)
350  {
351  // d/dxi
352  case 0:
353  switch (i)
354  {
355  case 0:
356  return .25*(1. - eta)*((1. - xi)*(-1.) +
357  (-1.)*(-1. - xi - eta));
358 
359  case 1:
360  return .25*(1. - eta)*((1. + xi)*(1.) +
361  (1.)*(-1. + xi - eta));
362 
363  case 2:
364  return .25*(1. + eta)*((1. + xi)*(1.) +
365  (1.)*(-1. + xi + eta));
366 
367  case 3:
368  return .25*(1. + eta)*((1. - xi)*(-1.) +
369  (-1.)*(-1. - xi + eta));
370 
371  case 4:
372  return .5*(-2.*xi)*(1. - eta);
373 
374  case 5:
375  return .5*(1.)*(1. - eta*eta);
376 
377  case 6:
378  return .5*(-2.*xi)*(1. + eta);
379 
380  case 7:
381  return .5*(-1.)*(1. - eta*eta);
382 
383  default:
384  libmesh_error_msg("Invalid shape function index i = " << i);
385  }
386 
387  // d/deta
388  case 1:
389  switch (i)
390  {
391  case 0:
392  return .25*(1. - xi)*((1. - eta)*(-1.) +
393  (-1.)*(-1. - xi - eta));
394 
395  case 1:
396  return .25*(1. + xi)*((1. - eta)*(-1.) +
397  (-1.)*(-1. + xi - eta));
398 
399  case 2:
400  return .25*(1. + xi)*((1. + eta)*(1.) +
401  (1.)*(-1. + xi + eta));
402 
403  case 3:
404  return .25*(1. - xi)*((1. + eta)*(1.) +
405  (1.)*(-1. - xi + eta));
406 
407  case 4:
408  return .5*(1. - xi*xi)*(-1.);
409 
410  case 5:
411  return .5*(1. + xi)*(-2.*eta);
412 
413  case 6:
414  return .5*(1. - xi*xi)*(1.);
415 
416  case 7:
417  return .5*(1. - xi)*(-2.*eta);
418 
419  default:
420  libmesh_error_msg("Invalid shape function index i = " << i);
421  }
422 
423  default:
424  libmesh_error_msg("ERROR: Invalid derivative index j = " << j);
425  }
426  }
427 
428  case QUAD9:
429  {
430  // Compute quad shape functions as a tensor-product
431  const Real xi = p(0);
432  const Real eta = p(1);
433 
434  libmesh_assert_less (i, 9);
435 
436  // 0 1 2 3 4 5 6 7 8
437  static const unsigned int i0[] = {0, 1, 1, 0, 2, 1, 2, 0, 2};
438  static const unsigned int i1[] = {0, 0, 1, 1, 0, 2, 1, 2, 2};
439 
440  switch (j)
441  {
442  // d()/dxi
443  case 0:
444  return (FE<1,L2_LAGRANGE>::shape_deriv(EDGE3, SECOND, i0[i], 0, xi)*
445  FE<1,L2_LAGRANGE>::shape (EDGE3, SECOND, i1[i], eta));
446 
447  // d()/deta
448  case 1:
449  return (FE<1,L2_LAGRANGE>::shape (EDGE3, SECOND, i0[i], xi)*
450  FE<1,L2_LAGRANGE>::shape_deriv(EDGE3, SECOND, i1[i], 0, eta));
451 
452  default:
453  libmesh_error_msg("ERROR: Invalid derivative index j = " << j);
454  }
455  }
456 
457  case TRI6:
458  {
459  libmesh_assert_less (i, 6);
460 
461  const Real zeta1 = p(0);
462  const Real zeta2 = p(1);
463  const Real zeta0 = 1. - zeta1 - zeta2;
464 
465  const Real dzeta0dxi = -1.;
466  const Real dzeta1dxi = 1.;
467  const Real dzeta2dxi = 0.;
468 
469  const Real dzeta0deta = -1.;
470  const Real dzeta1deta = 0.;
471  const Real dzeta2deta = 1.;
472 
473  switch(j)
474  {
475  case 0:
476  {
477  switch(i)
478  {
479  case 0:
480  return (4.*zeta0-1.)*dzeta0dxi;
481 
482  case 1:
483  return (4.*zeta1-1.)*dzeta1dxi;
484 
485  case 2:
486  return (4.*zeta2-1.)*dzeta2dxi;
487 
488  case 3:
489  return 4.*zeta1*dzeta0dxi + 4.*zeta0*dzeta1dxi;
490 
491  case 4:
492  return 4.*zeta2*dzeta1dxi + 4.*zeta1*dzeta2dxi;
493 
494  case 5:
495  return 4.*zeta2*dzeta0dxi + 4*zeta0*dzeta2dxi;
496 
497  default:
498  libmesh_error_msg("Invalid shape function index i = " << i);
499  }
500  }
501 
502  case 1:
503  {
504  switch(i)
505  {
506  case 0:
507  return (4.*zeta0-1.)*dzeta0deta;
508 
509  case 1:
510  return (4.*zeta1-1.)*dzeta1deta;
511 
512  case 2:
513  return (4.*zeta2-1.)*dzeta2deta;
514 
515  case 3:
516  return 4.*zeta1*dzeta0deta + 4.*zeta0*dzeta1deta;
517 
518  case 4:
519  return 4.*zeta2*dzeta1deta + 4.*zeta1*dzeta2deta;
520 
521  case 5:
522  return 4.*zeta2*dzeta0deta + 4*zeta0*dzeta2deta;
523 
524  default:
525  libmesh_error_msg("Invalid shape function index i = " << i);
526  }
527  }
528  default:
529  libmesh_error_msg("ERROR: Invalid derivative index j = " << j);
530  }
531  }
532 
533  default:
534  libmesh_error_msg("ERROR: Unsupported 2D element type: " << type);
535  }
536  }
537 
538 
539 
540  // unsupported order
541  default:
542  libmesh_error_msg("ERROR: Unsupported 2D FE order: " << order);
543  }
544 #endif
545 }
546 
547 
548 
549 template <>
551  const Order order,
552  const unsigned int i,
553  const unsigned int j,
554  const Point & p)
555 {
556  libmesh_assert(elem);
557 
558 
559  // call the orientation-independent shape functions
560  return FE<2,L2_LAGRANGE>::shape_deriv(elem->type(), static_cast<Order>(order + elem->p_level()), i, j, p);
561 }
562 
563 
564 
565 
566 template <>
568  const Order order,
569  const unsigned int i,
570  const unsigned int j,
571  const Point & p)
572 {
573 #if LIBMESH_DIM > 1
574 
575  // j = 0 ==> d^2 phi / dxi^2
576  // j = 1 ==> d^2 phi / dxi deta
577  // j = 2 ==> d^2 phi / deta^2
578  libmesh_assert_less (j, 3);
579 
580  switch (order)
581  {
582  // linear Lagrange shape functions
583  case FIRST:
584  {
585  switch (type)
586  {
587  case QUAD4:
588  case QUADSHELL4:
589  case QUAD8:
590  case QUADSHELL8:
591  case QUAD9:
592  {
593  // Compute quad shape functions as a tensor-product
594  const Real xi = p(0);
595  const Real eta = p(1);
596 
597  libmesh_assert_less (i, 4);
598 
599  // 0 1 2 3
600  static const unsigned int i0[] = {0, 1, 1, 0};
601  static const unsigned int i1[] = {0, 0, 1, 1};
602 
603  switch (j)
604  {
605  // d^2() / dxi^2
606  case 0:
607  return 0.;
608 
609  // d^2() / dxi deta
610  case 1:
611  return (FE<1,L2_LAGRANGE>::shape_deriv(EDGE2, FIRST, i0[i], 0, xi)*
612  FE<1,L2_LAGRANGE>::shape_deriv(EDGE2, FIRST, i1[i], 0, eta));
613 
614  // d^2() / deta^2
615  case 2:
616  return 0.;
617 
618  default:
619  libmesh_error_msg("ERROR: Invalid derivative index j = " << j);
620  }
621  }
622 
623  case TRI3:
624  case TRISHELL3:
625  case TRI6:
626  {
627  // All second derivatives for linear triangles are zero.
628  return 0.;
629  }
630 
631  default:
632  libmesh_error_msg("ERROR: Unsupported 2D element type: " << type);
633 
634  } // end switch (type)
635  } // end case FIRST
636 
637 
638  // quadratic Lagrange shape functions
639  case SECOND:
640  {
641  switch (type)
642  {
643  case QUAD8:
644  {
645  const Real xi = p(0);
646  const Real eta = p(1);
647 
648  libmesh_assert_less (j, 3);
649 
650  switch (j)
651  {
652  // d^2() / dxi^2
653  case 0:
654  {
655  switch (i)
656  {
657  case 0:
658  case 1:
659  return 0.5*(1.-eta);
660 
661  case 2:
662  case 3:
663  return 0.5*(1.+eta);
664 
665  case 4:
666  return eta - 1.;
667 
668  case 5:
669  case 7:
670  return 0.0;
671 
672  case 6:
673  return -1. - eta;
674 
675  default:
676  libmesh_error_msg("Invalid shape function index i = " << i);
677  }
678  }
679 
680  // d^2() / dxi deta
681  case 1:
682  {
683  switch (i)
684  {
685  case 0:
686  return 0.25*( 1. - 2.*xi - 2.*eta);
687 
688  case 1:
689  return 0.25*(-1. - 2.*xi + 2.*eta);
690 
691  case 2:
692  return 0.25*( 1. + 2.*xi + 2.*eta);
693 
694  case 3:
695  return 0.25*(-1. + 2.*xi - 2.*eta);
696 
697  case 4:
698  return xi;
699 
700  case 5:
701  return -eta;
702 
703  case 6:
704  return -xi;
705 
706  case 7:
707  return eta;
708 
709  default:
710  libmesh_error_msg("Invalid shape function index i = " << i);
711  }
712  }
713 
714  // d^2() / deta^2
715  case 2:
716  {
717  switch (i)
718  {
719  case 0:
720  case 3:
721  return 0.5*(1.-xi);
722 
723  case 1:
724  case 2:
725  return 0.5*(1.+xi);
726 
727  case 4:
728  case 6:
729  return 0.0;
730 
731  case 5:
732  return -1.0 - xi;
733 
734  case 7:
735  return xi - 1.0;
736 
737  default:
738  libmesh_error_msg("Invalid shape function index i = " << i);
739  }
740  }
741 
742 
743  default:
744  libmesh_error_msg("ERROR: Invalid derivative index j = " << j);
745  } // end switch (j)
746  } // end case QUAD8
747 
748  case QUAD9:
749  {
750  // Compute QUAD9 second derivatives as tensor product
751  const Real xi = p(0);
752  const Real eta = p(1);
753 
754  libmesh_assert_less (i, 9);
755 
756  // 0 1 2 3 4 5 6 7 8
757  static const unsigned int i0[] = {0, 1, 1, 0, 2, 1, 2, 0, 2};
758  static const unsigned int i1[] = {0, 0, 1, 1, 0, 2, 1, 2, 2};
759 
760  switch (j)
761  {
762  // d^2() / dxi^2
763  case 0:
764  return (FE<1,L2_LAGRANGE>::shape_second_deriv(EDGE3, SECOND, i0[i], 0, xi)*
765  FE<1,L2_LAGRANGE>::shape (EDGE3, SECOND, i1[i], eta));
766 
767  // d^2() / dxi deta
768  case 1:
769  return (FE<1,L2_LAGRANGE>::shape_deriv(EDGE3, SECOND, i0[i], 0, xi)*
770  FE<1,L2_LAGRANGE>::shape_deriv(EDGE3, SECOND, i1[i], 0, eta));
771 
772  // d^2() / deta^2
773  case 2:
774  return (FE<1,L2_LAGRANGE>::shape (EDGE3, SECOND, i0[i], xi)*
776 
777  default:
778  libmesh_error_msg("ERROR: Invalid derivative index j = " << j);
779  } // end switch (j)
780  } // end case QUAD9
781 
782  case TRI6:
783  {
784  const Real dzeta0dxi = -1.;
785  const Real dzeta1dxi = 1.;
786  const Real dzeta2dxi = 0.;
787 
788  const Real dzeta0deta = -1.;
789  const Real dzeta1deta = 0.;
790  const Real dzeta2deta = 1.;
791 
792  libmesh_assert_less (j, 3);
793 
794  switch (j)
795  {
796  // d^2() / dxi^2
797  case 0:
798  {
799  switch (i)
800  {
801  case 0:
802  return 4.*dzeta0dxi*dzeta0dxi;
803 
804  case 1:
805  return 4.*dzeta1dxi*dzeta1dxi;
806 
807  case 2:
808  return 4.*dzeta2dxi*dzeta2dxi;
809 
810  case 3:
811  return 8.*dzeta0dxi*dzeta1dxi;
812 
813  case 4:
814  return 8.*dzeta1dxi*dzeta2dxi;
815 
816  case 5:
817  return 8.*dzeta0dxi*dzeta2dxi;
818 
819  default:
820  libmesh_error_msg("Invalid shape function index i = " << i);
821  }
822  }
823 
824  // d^2() / dxi deta
825  case 1:
826  {
827  switch (i)
828  {
829  case 0:
830  return 4.*dzeta0dxi*dzeta0deta;
831 
832  case 1:
833  return 4.*dzeta1dxi*dzeta1deta;
834 
835  case 2:
836  return 4.*dzeta2dxi*dzeta2deta;
837 
838  case 3:
839  return 4.*dzeta1deta*dzeta0dxi + 4.*dzeta0deta*dzeta1dxi;
840 
841  case 4:
842  return 4.*dzeta2deta*dzeta1dxi + 4.*dzeta1deta*dzeta2dxi;
843 
844  case 5:
845  return 4.*dzeta2deta*dzeta0dxi + 4.*dzeta0deta*dzeta2dxi;
846 
847  default:
848  libmesh_error_msg("Invalid shape function index i = " << i);
849  }
850  }
851 
852  // d^2() / deta^2
853  case 2:
854  {
855  switch (i)
856  {
857  case 0:
858  return 4.*dzeta0deta*dzeta0deta;
859 
860  case 1:
861  return 4.*dzeta1deta*dzeta1deta;
862 
863  case 2:
864  return 4.*dzeta2deta*dzeta2deta;
865 
866  case 3:
867  return 8.*dzeta0deta*dzeta1deta;
868 
869  case 4:
870  return 8.*dzeta1deta*dzeta2deta;
871 
872  case 5:
873  return 8.*dzeta0deta*dzeta2deta;
874 
875  default:
876  libmesh_error_msg("Invalid shape function index i = " << i);
877  }
878  }
879 
880  default:
881  libmesh_error_msg("ERROR: Invalid derivative index j = " << j);
882  } // end switch (j)
883  } // end case TRI6
884 
885  default:
886  libmesh_error_msg("ERROR: Unsupported 2D element type: " << type);
887  }
888  } // end case SECOND
889 
890 
891 
892  // unsupported order
893  default:
894  libmesh_error_msg("ERROR: Unsupported 2D FE order: " << order);
895 
896  } // end switch (order)
897 
898 #endif // LIBMESH_DIM > 1
899 }
900 
901 
902 
903 template <>
905  const Order order,
906  const unsigned int i,
907  const unsigned int j,
908  const Point & p)
909 {
910  libmesh_assert(elem);
911 
912  // call the orientation-independent shape functions
913  return FE<2,L2_LAGRANGE>::shape_second_deriv(elem->type(), static_cast<Order>(order + elem->p_level()), i, j, p);
914 }
915 
916 } // 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)