$extrastylesheet
00001 // The libMesh Finite Element Library. 00002 // Copyright (C) 2002-2014 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner 00003 00004 // This library is free software; you can redistribute it and/or 00005 // modify it under the terms of the GNU Lesser General Public 00006 // License as published by the Free Software Foundation; either 00007 // version 2.1 of the License, or (at your option) any later version. 00008 00009 // This library is distributed in the hope that it will be useful, 00010 // but WITHOUT ANY WARRANTY; without even the implied warranty of 00011 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU 00012 // Lesser General Public License for more details. 00013 00014 // You should have received a copy of the GNU Lesser General Public 00015 // License along with this library; if not, write to the Free Software 00016 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA 00017 00018 00019 // C++ inlcludes 00020 00021 // Local includes 00022 #include "libmesh/fe.h" 00023 #include "libmesh/elem.h" 00024 00025 namespace libMesh 00026 { 00027 00028 00029 00030 00031 template <> 00032 Real FE<2,L2_LAGRANGE>::shape(const ElemType type, 00033 const Order order, 00034 const unsigned int i, 00035 const Point& p) 00036 { 00037 #if LIBMESH_DIM > 1 00038 00039 switch (order) 00040 { 00041 // linear Lagrange shape functions 00042 case FIRST: 00043 { 00044 switch (type) 00045 { 00046 case QUAD4: 00047 case QUAD8: 00048 case QUAD9: 00049 { 00050 // Compute quad shape functions as a tensor-product 00051 const Real xi = p(0); 00052 const Real eta = p(1); 00053 00054 libmesh_assert_less (i, 4); 00055 00056 // 0 1 2 3 00057 static const unsigned int i0[] = {0, 1, 1, 0}; 00058 static const unsigned int i1[] = {0, 0, 1, 1}; 00059 00060 return (FE<1,L2_LAGRANGE>::shape(EDGE2, FIRST, i0[i], xi)* 00061 FE<1,L2_LAGRANGE>::shape(EDGE2, FIRST, i1[i], eta)); 00062 } 00063 00064 case TRI3: 00065 case TRI6: 00066 { 00067 const Real zeta1 = p(0); 00068 const Real zeta2 = p(1); 00069 const Real zeta0 = 1. - zeta1 - zeta2; 00070 00071 libmesh_assert_less (i, 3); 00072 00073 switch(i) 00074 { 00075 case 0: 00076 return zeta0; 00077 00078 case 1: 00079 return zeta1; 00080 00081 case 2: 00082 return zeta2; 00083 00084 default: 00085 libmesh_error_msg("Invalid shape function index i = " << i); 00086 } 00087 } 00088 00089 default: 00090 libmesh_error_msg("ERROR: Unsupported 2D element type: " << type); 00091 } 00092 } 00093 00094 00095 // quadratic Lagrange shape functions 00096 case SECOND: 00097 { 00098 switch (type) 00099 { 00100 case QUAD8: 00101 { 00102 const Real xi = p(0); 00103 const Real eta = p(1); 00104 00105 libmesh_assert_less (i, 8); 00106 00107 switch (i) 00108 { 00109 case 0: 00110 return .25*(1. - xi)*(1. - eta)*(-1. - xi - eta); 00111 00112 case 1: 00113 return .25*(1. + xi)*(1. - eta)*(-1. + xi - eta); 00114 00115 case 2: 00116 return .25*(1. + xi)*(1. + eta)*(-1. + xi + eta); 00117 00118 case 3: 00119 return .25*(1. - xi)*(1. + eta)*(-1. - xi + eta); 00120 00121 case 4: 00122 return .5*(1. - xi*xi)*(1. - eta); 00123 00124 case 5: 00125 return .5*(1. + xi)*(1. - eta*eta); 00126 00127 case 6: 00128 return .5*(1. - xi*xi)*(1. + eta); 00129 00130 case 7: 00131 return .5*(1. - xi)*(1. - eta*eta); 00132 00133 default: 00134 libmesh_error_msg("Invalid shape function index i = " << i); 00135 } 00136 } 00137 00138 case QUAD9: 00139 { 00140 // Compute quad shape functions as a tensor-product 00141 const Real xi = p(0); 00142 const Real eta = p(1); 00143 00144 libmesh_assert_less (i, 9); 00145 00146 // 0 1 2 3 4 5 6 7 8 00147 static const unsigned int i0[] = {0, 1, 1, 0, 2, 1, 2, 0, 2}; 00148 static const unsigned int i1[] = {0, 0, 1, 1, 0, 2, 1, 2, 2}; 00149 00150 return (FE<1,L2_LAGRANGE>::shape(EDGE3, SECOND, i0[i], xi)* 00151 FE<1,L2_LAGRANGE>::shape(EDGE3, SECOND, i1[i], eta)); 00152 } 00153 00154 case TRI6: 00155 { 00156 const Real zeta1 = p(0); 00157 const Real zeta2 = p(1); 00158 const Real zeta0 = 1. - zeta1 - zeta2; 00159 00160 libmesh_assert_less (i, 6); 00161 00162 switch(i) 00163 { 00164 case 0: 00165 return 2.*zeta0*(zeta0-0.5); 00166 00167 case 1: 00168 return 2.*zeta1*(zeta1-0.5); 00169 00170 case 2: 00171 return 2.*zeta2*(zeta2-0.5); 00172 00173 case 3: 00174 return 4.*zeta0*zeta1; 00175 00176 case 4: 00177 return 4.*zeta1*zeta2; 00178 00179 case 5: 00180 return 4.*zeta2*zeta0; 00181 00182 default: 00183 libmesh_error_msg("Invalid shape function index i = " << i); 00184 } 00185 } 00186 00187 default: 00188 libmesh_error_msg("ERROR: Unsupported 2D element type: " << type); 00189 } 00190 } 00191 00192 00193 00194 // unsupported order 00195 default: 00196 libmesh_error_msg("ERROR: Unsupported 2D FE order: " << order); 00197 } 00198 00199 libmesh_error_msg("We'll never get here!"); 00200 return 0.; 00201 00202 #endif 00203 } 00204 00205 00206 00207 template <> 00208 Real FE<2,L2_LAGRANGE>::shape(const Elem* elem, 00209 const Order order, 00210 const unsigned int i, 00211 const Point& p) 00212 { 00213 libmesh_assert(elem); 00214 00215 // call the orientation-independent shape functions 00216 return FE<2,L2_LAGRANGE>::shape(elem->type(), static_cast<Order>(order + elem->p_level()), i, p); 00217 } 00218 00219 00220 00221 template <> 00222 Real FE<2,L2_LAGRANGE>::shape_deriv(const ElemType type, 00223 const Order order, 00224 const unsigned int i, 00225 const unsigned int j, 00226 const Point& p) 00227 { 00228 #if LIBMESH_DIM > 1 00229 00230 00231 libmesh_assert_less (j, 2); 00232 00233 switch (order) 00234 { 00235 // linear Lagrange shape functions 00236 case FIRST: 00237 { 00238 switch (type) 00239 { 00240 case QUAD4: 00241 case QUAD8: 00242 case QUAD9: 00243 { 00244 // Compute quad shape functions as a tensor-product 00245 const Real xi = p(0); 00246 const Real eta = p(1); 00247 00248 libmesh_assert_less (i, 4); 00249 00250 // 0 1 2 3 00251 static const unsigned int i0[] = {0, 1, 1, 0}; 00252 static const unsigned int i1[] = {0, 0, 1, 1}; 00253 00254 switch (j) 00255 { 00256 // d()/dxi 00257 case 0: 00258 return (FE<1,L2_LAGRANGE>::shape_deriv(EDGE2, FIRST, i0[i], 0, xi)* 00259 FE<1,L2_LAGRANGE>::shape (EDGE2, FIRST, i1[i], eta)); 00260 00261 // d()/deta 00262 case 1: 00263 return (FE<1,L2_LAGRANGE>::shape (EDGE2, FIRST, i0[i], xi)* 00264 FE<1,L2_LAGRANGE>::shape_deriv(EDGE2, FIRST, i1[i], 0, eta)); 00265 00266 default: 00267 libmesh_error_msg("ERROR: Invalid derivative index j = " << j); 00268 } 00269 } 00270 00271 case TRI3: 00272 case TRI6: 00273 { 00274 libmesh_assert_less (i, 3); 00275 00276 const Real dzeta0dxi = -1.; 00277 const Real dzeta1dxi = 1.; 00278 const Real dzeta2dxi = 0.; 00279 00280 const Real dzeta0deta = -1.; 00281 const Real dzeta1deta = 0.; 00282 const Real dzeta2deta = 1.; 00283 00284 switch (j) 00285 { 00286 // d()/dxi 00287 case 0: 00288 { 00289 switch(i) 00290 { 00291 case 0: 00292 return dzeta0dxi; 00293 00294 case 1: 00295 return dzeta1dxi; 00296 00297 case 2: 00298 return dzeta2dxi; 00299 00300 default: 00301 libmesh_error_msg("Invalid shape function index i = " << i); 00302 } 00303 } 00304 // d()/deta 00305 case 1: 00306 { 00307 switch(i) 00308 { 00309 case 0: 00310 return dzeta0deta; 00311 00312 case 1: 00313 return dzeta1deta; 00314 00315 case 2: 00316 return dzeta2deta; 00317 00318 default: 00319 libmesh_error_msg("Invalid shape function index i = " << i); 00320 } 00321 } 00322 default: 00323 libmesh_error_msg("ERROR: Invalid derivative index j = " << j); 00324 } 00325 } 00326 00327 default: 00328 libmesh_error_msg("ERROR: Unsupported 2D element type: " << type); 00329 } 00330 } 00331 00332 00333 // quadratic Lagrange shape functions 00334 case SECOND: 00335 { 00336 switch (type) 00337 { 00338 case QUAD8: 00339 { 00340 const Real xi = p(0); 00341 const Real eta = p(1); 00342 00343 libmesh_assert_less (i, 8); 00344 00345 switch (j) 00346 { 00347 // d/dxi 00348 case 0: 00349 switch (i) 00350 { 00351 case 0: 00352 return .25*(1. - eta)*((1. - xi)*(-1.) + 00353 (-1.)*(-1. - xi - eta)); 00354 00355 case 1: 00356 return .25*(1. - eta)*((1. + xi)*(1.) + 00357 (1.)*(-1. + xi - eta)); 00358 00359 case 2: 00360 return .25*(1. + eta)*((1. + xi)*(1.) + 00361 (1.)*(-1. + xi + eta)); 00362 00363 case 3: 00364 return .25*(1. + eta)*((1. - xi)*(-1.) + 00365 (-1.)*(-1. - xi + eta)); 00366 00367 case 4: 00368 return .5*(-2.*xi)*(1. - eta); 00369 00370 case 5: 00371 return .5*(1.)*(1. - eta*eta); 00372 00373 case 6: 00374 return .5*(-2.*xi)*(1. + eta); 00375 00376 case 7: 00377 return .5*(-1.)*(1. - eta*eta); 00378 00379 default: 00380 libmesh_error_msg("Invalid shape function index i = " << i); 00381 } 00382 00383 // d/deta 00384 case 1: 00385 switch (i) 00386 { 00387 case 0: 00388 return .25*(1. - xi)*((1. - eta)*(-1.) + 00389 (-1.)*(-1. - xi - eta)); 00390 00391 case 1: 00392 return .25*(1. + xi)*((1. - eta)*(-1.) + 00393 (-1.)*(-1. + xi - eta)); 00394 00395 case 2: 00396 return .25*(1. + xi)*((1. + eta)*(1.) + 00397 (1.)*(-1. + xi + eta)); 00398 00399 case 3: 00400 return .25*(1. - xi)*((1. + eta)*(1.) + 00401 (1.)*(-1. - xi + eta)); 00402 00403 case 4: 00404 return .5*(1. - xi*xi)*(-1.); 00405 00406 case 5: 00407 return .5*(1. + xi)*(-2.*eta); 00408 00409 case 6: 00410 return .5*(1. - xi*xi)*(1.); 00411 00412 case 7: 00413 return .5*(1. - xi)*(-2.*eta); 00414 00415 default: 00416 libmesh_error_msg("Invalid shape function index i = " << i); 00417 } 00418 00419 default: 00420 libmesh_error_msg("ERROR: Invalid derivative index j = " << j); 00421 } 00422 } 00423 00424 case QUAD9: 00425 { 00426 // Compute quad shape functions as a tensor-product 00427 const Real xi = p(0); 00428 const Real eta = p(1); 00429 00430 libmesh_assert_less (i, 9); 00431 00432 // 0 1 2 3 4 5 6 7 8 00433 static const unsigned int i0[] = {0, 1, 1, 0, 2, 1, 2, 0, 2}; 00434 static const unsigned int i1[] = {0, 0, 1, 1, 0, 2, 1, 2, 2}; 00435 00436 switch (j) 00437 { 00438 // d()/dxi 00439 case 0: 00440 return (FE<1,L2_LAGRANGE>::shape_deriv(EDGE3, SECOND, i0[i], 0, xi)* 00441 FE<1,L2_LAGRANGE>::shape (EDGE3, SECOND, i1[i], eta)); 00442 00443 // d()/deta 00444 case 1: 00445 return (FE<1,L2_LAGRANGE>::shape (EDGE3, SECOND, i0[i], xi)* 00446 FE<1,L2_LAGRANGE>::shape_deriv(EDGE3, SECOND, i1[i], 0, eta)); 00447 00448 default: 00449 libmesh_error_msg("ERROR: Invalid derivative index j = " << j); 00450 } 00451 } 00452 00453 case TRI6: 00454 { 00455 libmesh_assert_less (i, 6); 00456 00457 const Real zeta1 = p(0); 00458 const Real zeta2 = p(1); 00459 const Real zeta0 = 1. - zeta1 - zeta2; 00460 00461 const Real dzeta0dxi = -1.; 00462 const Real dzeta1dxi = 1.; 00463 const Real dzeta2dxi = 0.; 00464 00465 const Real dzeta0deta = -1.; 00466 const Real dzeta1deta = 0.; 00467 const Real dzeta2deta = 1.; 00468 00469 switch(j) 00470 { 00471 case 0: 00472 { 00473 switch(i) 00474 { 00475 case 0: 00476 return (4.*zeta0-1.)*dzeta0dxi; 00477 00478 case 1: 00479 return (4.*zeta1-1.)*dzeta1dxi; 00480 00481 case 2: 00482 return (4.*zeta2-1.)*dzeta2dxi; 00483 00484 case 3: 00485 return 4.*zeta1*dzeta0dxi + 4.*zeta0*dzeta1dxi; 00486 00487 case 4: 00488 return 4.*zeta2*dzeta1dxi + 4.*zeta1*dzeta2dxi; 00489 00490 case 5: 00491 return 4.*zeta2*dzeta0dxi + 4*zeta0*dzeta2dxi; 00492 00493 default: 00494 libmesh_error_msg("Invalid shape function index i = " << i); 00495 } 00496 } 00497 00498 case 1: 00499 { 00500 switch(i) 00501 { 00502 case 0: 00503 return (4.*zeta0-1.)*dzeta0deta; 00504 00505 case 1: 00506 return (4.*zeta1-1.)*dzeta1deta; 00507 00508 case 2: 00509 return (4.*zeta2-1.)*dzeta2deta; 00510 00511 case 3: 00512 return 4.*zeta1*dzeta0deta + 4.*zeta0*dzeta1deta; 00513 00514 case 4: 00515 return 4.*zeta2*dzeta1deta + 4.*zeta1*dzeta2deta; 00516 00517 case 5: 00518 return 4.*zeta2*dzeta0deta + 4*zeta0*dzeta2deta; 00519 00520 default: 00521 libmesh_error_msg("Invalid shape function index i = " << i); 00522 } 00523 } 00524 default: 00525 libmesh_error_msg("ERROR: Invalid derivative index j = " << j); 00526 } 00527 } 00528 00529 default: 00530 libmesh_error_msg("ERROR: Unsupported 2D element type: " << type); 00531 } 00532 } 00533 00534 00535 00536 // unsupported order 00537 default: 00538 libmesh_error_msg("ERROR: Unsupported 2D FE order: " << order); 00539 } 00540 00541 00542 libmesh_error_msg("We'll never get here!"); 00543 return 0.; 00544 00545 #endif 00546 } 00547 00548 00549 00550 template <> 00551 Real FE<2,L2_LAGRANGE>::shape_deriv(const Elem* elem, 00552 const Order order, 00553 const unsigned int i, 00554 const unsigned int j, 00555 const Point& p) 00556 { 00557 libmesh_assert(elem); 00558 00559 00560 // call the orientation-independent shape functions 00561 return FE<2,L2_LAGRANGE>::shape_deriv(elem->type(), static_cast<Order>(order + elem->p_level()), i, j, p); 00562 } 00563 00564 00565 00566 00567 template <> 00568 Real FE<2,L2_LAGRANGE>::shape_second_deriv(const ElemType type, 00569 const Order order, 00570 const unsigned int i, 00571 const unsigned int j, 00572 const Point& p) 00573 { 00574 #if LIBMESH_DIM > 1 00575 00576 // j = 0 ==> d^2 phi / dxi^2 00577 // j = 1 ==> d^2 phi / dxi deta 00578 // j = 2 ==> d^2 phi / deta^2 00579 libmesh_assert_less (j, 3); 00580 00581 switch (order) 00582 { 00583 // linear Lagrange shape functions 00584 case FIRST: 00585 { 00586 switch (type) 00587 { 00588 case QUAD4: 00589 case QUAD8: 00590 case QUAD9: 00591 { 00592 // Compute quad shape functions as a tensor-product 00593 const Real xi = p(0); 00594 const Real eta = p(1); 00595 00596 libmesh_assert_less (i, 4); 00597 00598 // 0 1 2 3 00599 static const unsigned int i0[] = {0, 1, 1, 0}; 00600 static const unsigned int i1[] = {0, 0, 1, 1}; 00601 00602 switch (j) 00603 { 00604 // d^2() / dxi^2 00605 case 0: 00606 return 0.; 00607 00608 // d^2() / dxi deta 00609 case 1: 00610 return (FE<1,L2_LAGRANGE>::shape_deriv(EDGE2, FIRST, i0[i], 0, xi)* 00611 FE<1,L2_LAGRANGE>::shape_deriv(EDGE2, FIRST, i1[i], 0, eta)); 00612 00613 // d^2() / deta^2 00614 case 2: 00615 return 0.; 00616 00617 default: 00618 libmesh_error_msg("ERROR: Invalid derivative index j = " << j); 00619 } 00620 } 00621 00622 case TRI3: 00623 case TRI6: 00624 { 00625 // All second derivatives for linear triangles are zero. 00626 return 0.; 00627 } 00628 00629 default: 00630 libmesh_error_msg("ERROR: Unsupported 2D element type: " << type); 00631 00632 } // end switch (type) 00633 } // end case FIRST 00634 00635 00636 // quadratic Lagrange shape functions 00637 case SECOND: 00638 { 00639 switch (type) 00640 { 00641 case QUAD8: 00642 { 00643 const Real xi = p(0); 00644 const Real eta = p(1); 00645 00646 libmesh_assert_less (j, 3); 00647 00648 switch (j) 00649 { 00650 // d^2() / dxi^2 00651 case 0: 00652 { 00653 switch (i) 00654 { 00655 case 0: 00656 case 1: 00657 return 0.5*(1.-eta); 00658 00659 case 2: 00660 case 3: 00661 return 0.5*(1.+eta); 00662 00663 case 4: 00664 return eta - 1.; 00665 00666 case 5: 00667 case 7: 00668 return 0.0; 00669 00670 case 6: 00671 return -1. - eta; 00672 00673 default: 00674 libmesh_error_msg("Invalid shape function index i = " << i); 00675 } 00676 } 00677 00678 // d^2() / dxi deta 00679 case 1: 00680 { 00681 switch (i) 00682 { 00683 case 0: 00684 return 0.25*( 1. - 2.*xi - 2.*eta); 00685 00686 case 1: 00687 return 0.25*(-1. - 2.*xi + 2.*eta); 00688 00689 case 2: 00690 return 0.25*( 1. + 2.*xi + 2.*eta); 00691 00692 case 3: 00693 return 0.25*(-1. + 2.*xi - 2.*eta); 00694 00695 case 4: 00696 return xi; 00697 00698 case 5: 00699 return -eta; 00700 00701 case 6: 00702 return -xi; 00703 00704 case 7: 00705 return eta; 00706 00707 default: 00708 libmesh_error_msg("Invalid shape function index i = " << i); 00709 } 00710 } 00711 00712 // d^2() / deta^2 00713 case 2: 00714 { 00715 switch (i) 00716 { 00717 case 0: 00718 case 3: 00719 return 0.5*(1.-xi); 00720 00721 case 1: 00722 case 2: 00723 return 0.5*(1.+xi); 00724 00725 case 4: 00726 case 6: 00727 return 0.0; 00728 00729 case 5: 00730 return -1.0 - xi; 00731 00732 case 7: 00733 return xi - 1.0; 00734 00735 default: 00736 libmesh_error_msg("Invalid shape function index i = " << i); 00737 } 00738 } 00739 00740 00741 default: 00742 libmesh_error_msg("ERROR: Invalid derivative index j = " << j); 00743 } // end switch (j) 00744 } // end case QUAD8 00745 00746 case QUAD9: 00747 { 00748 // Compute QUAD9 second derivatives as tensor product 00749 const Real xi = p(0); 00750 const Real eta = p(1); 00751 00752 libmesh_assert_less (i, 9); 00753 00754 // 0 1 2 3 4 5 6 7 8 00755 static const unsigned int i0[] = {0, 1, 1, 0, 2, 1, 2, 0, 2}; 00756 static const unsigned int i1[] = {0, 0, 1, 1, 0, 2, 1, 2, 2}; 00757 00758 switch (j) 00759 { 00760 // d^2() / dxi^2 00761 case 0: 00762 return (FE<1,L2_LAGRANGE>::shape_second_deriv(EDGE3, SECOND, i0[i], 0, xi)* 00763 FE<1,L2_LAGRANGE>::shape (EDGE3, SECOND, i1[i], eta)); 00764 00765 // d^2() / dxi deta 00766 case 1: 00767 return (FE<1,L2_LAGRANGE>::shape_deriv(EDGE3, SECOND, i0[i], 0, xi)* 00768 FE<1,L2_LAGRANGE>::shape_deriv(EDGE3, SECOND, i1[i], 0, eta)); 00769 00770 // d^2() / deta^2 00771 case 2: 00772 return (FE<1,L2_LAGRANGE>::shape (EDGE3, SECOND, i0[i], xi)* 00773 FE<1,L2_LAGRANGE>::shape_second_deriv(EDGE3, SECOND, i1[i], 0, eta)); 00774 00775 default: 00776 libmesh_error_msg("ERROR: Invalid derivative index j = " << j); 00777 } // end switch (j) 00778 } // end case QUAD9 00779 00780 case TRI6: 00781 { 00782 const Real dzeta0dxi = -1.; 00783 const Real dzeta1dxi = 1.; 00784 const Real dzeta2dxi = 0.; 00785 00786 const Real dzeta0deta = -1.; 00787 const Real dzeta1deta = 0.; 00788 const Real dzeta2deta = 1.; 00789 00790 libmesh_assert_less (j, 3); 00791 00792 switch (j) 00793 { 00794 // d^2() / dxi^2 00795 case 0: 00796 { 00797 switch (i) 00798 { 00799 case 0: 00800 return 4.*dzeta0dxi*dzeta0dxi; 00801 00802 case 1: 00803 return 4.*dzeta1dxi*dzeta1dxi; 00804 00805 case 2: 00806 return 4.*dzeta2dxi*dzeta2dxi; 00807 00808 case 3: 00809 return 8.*dzeta0dxi*dzeta1dxi; 00810 00811 case 4: 00812 return 8.*dzeta1dxi*dzeta2dxi; 00813 00814 case 5: 00815 return 8.*dzeta0dxi*dzeta2dxi; 00816 00817 default: 00818 libmesh_error_msg("Invalid shape function index i = " << i); 00819 } 00820 } 00821 00822 // d^2() / dxi deta 00823 case 1: 00824 { 00825 switch (i) 00826 { 00827 case 0: 00828 return 4.*dzeta0dxi*dzeta0deta; 00829 00830 case 1: 00831 return 4.*dzeta1dxi*dzeta1deta; 00832 00833 case 2: 00834 return 4.*dzeta2dxi*dzeta2deta; 00835 00836 case 3: 00837 return 4.*dzeta1deta*dzeta0dxi + 4.*dzeta0deta*dzeta1dxi; 00838 00839 case 4: 00840 return 4.*dzeta2deta*dzeta1dxi + 4.*dzeta1deta*dzeta2dxi; 00841 00842 case 5: 00843 return 4.*dzeta2deta*dzeta0dxi + 4.*dzeta0deta*dzeta2dxi; 00844 00845 default: 00846 libmesh_error_msg("Invalid shape function index i = " << i); 00847 } 00848 } 00849 00850 // d^2() / deta^2 00851 case 2: 00852 { 00853 switch (i) 00854 { 00855 case 0: 00856 return 4.*dzeta0deta*dzeta0deta; 00857 00858 case 1: 00859 return 4.*dzeta1deta*dzeta1deta; 00860 00861 case 2: 00862 return 4.*dzeta2deta*dzeta2deta; 00863 00864 case 3: 00865 return 8.*dzeta0deta*dzeta1deta; 00866 00867 case 4: 00868 return 8.*dzeta1deta*dzeta2deta; 00869 00870 case 5: 00871 return 8.*dzeta0deta*dzeta2deta; 00872 00873 default: 00874 libmesh_error_msg("Invalid shape function index i = " << i); 00875 } 00876 } 00877 00878 default: 00879 libmesh_error_msg("ERROR: Invalid derivative index j = " << j); 00880 } // end switch (j) 00881 } // end case TRI6 00882 00883 default: 00884 libmesh_error_msg("ERROR: Unsupported 2D element type: " << type); 00885 } 00886 } // end case SECOND 00887 00888 00889 00890 // unsupported order 00891 default: 00892 libmesh_error_msg("ERROR: Unsupported 2D FE order: " << order); 00893 00894 } // end switch (order) 00895 00896 00897 libmesh_error_msg("We'll never get here!"); 00898 return 0.; 00899 #endif // LIBMESH_DIM > 1 00900 } 00901 00902 00903 00904 template <> 00905 Real FE<2,L2_LAGRANGE>::shape_second_deriv(const Elem* elem, 00906 const Order order, 00907 const unsigned int i, 00908 const unsigned int j, 00909 const Point& p) 00910 { 00911 libmesh_assert(elem); 00912 00913 // call the orientation-independent shape functions 00914 return FE<2,L2_LAGRANGE>::shape_second_deriv(elem->type(), static_cast<Order>(order + elem->p_level()), i, j, p); 00915 } 00916 00917 } // namespace libMesh