$extrastylesheet
fe_monomial_shape_2D.C
Go to the documentation of this file.
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,MONOMIAL>::shape(const ElemType,
00033                            const Order libmesh_dbg_var(order),
00034                            const unsigned int i,
00035                            const Point& p)
00036 {
00037 #if LIBMESH_DIM > 1
00038 
00039   libmesh_assert_less (i, (static_cast<unsigned int>(order)+1)*
00040                        (static_cast<unsigned int>(order)+2)/2);
00041 
00042   const Real xi  = p(0);
00043   const Real eta = p(1);
00044 
00045   switch (i)
00046     {
00047       // constant
00048     case 0:
00049       return 1.;
00050 
00051       // linear
00052     case 1:
00053       return xi;
00054 
00055     case 2:
00056       return eta;
00057 
00058       // quadratics
00059     case 3:
00060       return xi*xi;
00061 
00062     case 4:
00063       return xi*eta;
00064 
00065     case 5:
00066       return eta*eta;
00067 
00068       // cubics
00069     case 6:
00070       return xi*xi*xi;
00071 
00072     case 7:
00073       return xi*xi*eta;
00074 
00075     case 8:
00076       return xi*eta*eta;
00077 
00078     case 9:
00079       return eta*eta*eta;
00080 
00081       // quartics
00082     case 10:
00083       return xi*xi*xi*xi;
00084 
00085     case 11:
00086       return xi*xi*xi*eta;
00087 
00088     case 12:
00089       return xi*xi*eta*eta;
00090 
00091     case 13:
00092       return xi*eta*eta*eta;
00093 
00094     case 14:
00095       return eta*eta*eta*eta;
00096 
00097     default:
00098       unsigned int o = 0;
00099       for (; i >= (o+1)*(o+2)/2; o++) { }
00100       unsigned int ny = i - (o*(o+1)/2);
00101       unsigned int nx = o - ny;
00102       Real val = 1.;
00103       for (unsigned int index=0; index != nx; index++)
00104         val *= xi;
00105       for (unsigned int index=0; index != ny; index++)
00106         val *= eta;
00107       return val;
00108     }
00109 
00110   libmesh_error_msg("We'll never get here!");
00111   return 0.;
00112 
00113 #endif
00114 }
00115 
00116 
00117 
00118 template <>
00119 Real FE<2,MONOMIAL>::shape(const Elem* elem,
00120                            const Order order,
00121                            const unsigned int i,
00122                            const Point& p)
00123 {
00124   libmesh_assert(elem);
00125 
00126   // by default call the orientation-independent shape functions
00127   return FE<2,MONOMIAL>::shape(elem->type(), static_cast<Order>(order + elem->p_level()), i, p);
00128 }
00129 
00130 
00131 
00132 template <>
00133 Real FE<2,MONOMIAL>::shape_deriv(const ElemType,
00134                                  const Order libmesh_dbg_var(order),
00135                                  const unsigned int i,
00136                                  const unsigned int j,
00137                                  const Point& p)
00138 {
00139 #if LIBMESH_DIM > 1
00140 
00141 
00142   libmesh_assert_less (j, 2);
00143 
00144   libmesh_assert_less (i, (static_cast<unsigned int>(order)+1)*
00145                        (static_cast<unsigned int>(order)+2)/2);
00146 
00147   const Real xi  = p(0);
00148   const Real eta = p(1);
00149 
00150   // monomials. since they are hierarchic we only need one case block.
00151 
00152   switch (j)
00153     {
00154       // d()/dxi
00155     case 0:
00156       {
00157         switch (i)
00158           {
00159             // constants
00160           case 0:
00161             return 0.;
00162 
00163             // linears
00164           case 1:
00165             return 1.;
00166 
00167           case 2:
00168             return 0.;
00169 
00170             // quadratics
00171           case 3:
00172             return 2.*xi;
00173 
00174           case 4:
00175             return eta;
00176 
00177           case 5:
00178             return 0.;
00179 
00180             // cubics
00181           case 6:
00182             return 3.*xi*xi;
00183 
00184           case 7:
00185             return 2.*xi*eta;
00186 
00187           case 8:
00188             return eta*eta;
00189 
00190           case 9:
00191             return 0.;
00192 
00193             // quartics
00194           case 10:
00195             return 4.*xi*xi*xi;
00196 
00197           case 11:
00198             return 3.*xi*xi*eta;
00199 
00200           case 12:
00201             return 2.*xi*eta*eta;
00202 
00203           case 13:
00204             return eta*eta*eta;
00205 
00206           case 14:
00207             return 0.;
00208 
00209           default:
00210             unsigned int o = 0;
00211             for (; i >= (o+1)*(o+2)/2; o++) { }
00212             unsigned int ny = i - (o*(o+1)/2);
00213             unsigned int nx = o - ny;
00214             Real val = nx;
00215             for (unsigned int index=1; index < nx; index++)
00216               val *= xi;
00217             for (unsigned int index=0; index != ny; index++)
00218               val *= eta;
00219             return val;
00220           }
00221       }
00222 
00223 
00224       // d()/deta
00225     case 1:
00226       {
00227         switch (i)
00228           {
00229             // constants
00230           case 0:
00231             return 0.;
00232 
00233             // linears
00234           case 1:
00235             return 0.;
00236 
00237           case 2:
00238             return 1.;
00239 
00240             // quadratics
00241           case 3:
00242             return 0.;
00243 
00244           case 4:
00245             return xi;
00246 
00247           case 5:
00248             return 2.*eta;
00249 
00250             // cubics
00251           case 6:
00252             return 0.;
00253 
00254           case 7:
00255             return xi*xi;
00256 
00257           case 8:
00258             return 2.*xi*eta;
00259 
00260           case 9:
00261             return 3.*eta*eta;
00262 
00263             // quartics
00264           case 10:
00265             return 0.;
00266 
00267           case 11:
00268             return xi*xi*xi;
00269 
00270           case 12:
00271             return 2.*xi*xi*eta;
00272 
00273           case 13:
00274             return 3.*xi*eta*eta;
00275 
00276           case 14:
00277             return 4.*eta*eta*eta;
00278 
00279           default:
00280             unsigned int o = 0;
00281             for (; i >= (o+1)*(o+2)/2; o++) { }
00282             unsigned int ny = i - (o*(o+1)/2);
00283             unsigned int nx = o - ny;
00284             Real val = ny;
00285             for (unsigned int index=0; index != nx; index++)
00286               val *= xi;
00287             for (unsigned int index=1; index < ny; index++)
00288               val *= eta;
00289             return val;
00290           }
00291       }
00292 
00293     default:
00294       libmesh_error_msg("Invalid shape function derivative j = " << j);
00295     }
00296 
00297   libmesh_error_msg("We'll never get here!");
00298   return 0.;
00299 
00300 #endif
00301 }
00302 
00303 
00304 
00305 template <>
00306 Real FE<2,MONOMIAL>::shape_deriv(const Elem* elem,
00307                                  const Order order,
00308                                  const unsigned int i,
00309                                  const unsigned int j,
00310                                  const Point& p)
00311 {
00312   libmesh_assert(elem);
00313 
00314   // by default call the orientation-independent shape functions
00315   return FE<2,MONOMIAL>::shape_deriv(elem->type(), static_cast<Order>(order + elem->p_level()), i, j, p);
00316 }
00317 
00318 
00319 
00320 template <>
00321 Real FE<2,MONOMIAL>::shape_second_deriv(const ElemType,
00322                                         const Order libmesh_dbg_var(order),
00323                                         const unsigned int i,
00324                                         const unsigned int j,
00325                                         const Point& p)
00326 {
00327 #if LIBMESH_DIM > 1
00328 
00329 
00330   libmesh_assert_less_equal (j, 2);
00331 
00332   libmesh_assert_less (i, (static_cast<unsigned int>(order)+1)*
00333                        (static_cast<unsigned int>(order)+2)/2);
00334 
00335   const Real xi  = p(0);
00336   const Real eta = p(1);
00337 
00338   // monomials. since they are hierarchic we only need one case block.
00339 
00340   switch (j)
00341     {
00342       // d^2()/dxi^2
00343     case 0:
00344       {
00345         switch (i)
00346           {
00347             // constants
00348           case 0:
00349             // linears
00350           case 1:
00351           case 2:
00352             return 0.;
00353 
00354             // quadratics
00355           case 3:
00356             return 2.;
00357 
00358           case 4:
00359           case 5:
00360             return 0.;
00361 
00362             // cubics
00363           case 6:
00364             return 6.*xi;
00365 
00366           case 7:
00367             return 2.*eta;
00368 
00369           case 8:
00370           case 9:
00371             return 0.;
00372 
00373             // quartics
00374           case 10:
00375             return 12.*xi*xi;
00376 
00377           case 11:
00378             return 6.*xi*eta;
00379 
00380           case 12:
00381             return 2.*eta*eta;
00382 
00383           case 13:
00384           case 14:
00385             return 0.;
00386 
00387           default:
00388             unsigned int o = 0;
00389             for (; i >= (o+1)*(o+2)/2; o++) { }
00390             unsigned int ny = i - (o*(o+1)/2);
00391             unsigned int nx = o - ny;
00392             Real val = nx * (nx - 1);
00393             for (unsigned int index=2; index < nx; index++)
00394               val *= xi;
00395             for (unsigned int index=0; index != ny; index++)
00396               val *= eta;
00397             return val;
00398           }
00399       }
00400 
00401       // d^2()/dxideta
00402     case 1:
00403       {
00404         switch (i)
00405           {
00406             // constants
00407           case 0:
00408 
00409             // linears
00410           case 1:
00411           case 2:
00412             return 0.;
00413 
00414             // quadratics
00415           case 3:
00416             return 0.;
00417 
00418           case 4:
00419             return 1.;
00420 
00421           case 5:
00422             return 0.;
00423 
00424             // cubics
00425           case 6:
00426             return 0.;
00427           case 7:
00428             return 2.*xi;
00429 
00430           case 8:
00431             return 2.*eta;
00432 
00433           case 9:
00434             return 0.;
00435 
00436             // quartics
00437           case 10:
00438             return 0.;
00439 
00440           case 11:
00441             return 3.*xi*xi;
00442 
00443           case 12:
00444             return 4.*xi*eta;
00445 
00446           case 13:
00447             return 3.*eta*eta;
00448 
00449           case 14:
00450             return 0.;
00451 
00452           default:
00453             unsigned int o = 0;
00454             for (; i >= (o+1)*(o+2)/2; o++) { }
00455             unsigned int ny = i - (o*(o+1)/2);
00456             unsigned int nx = o - ny;
00457             Real val = nx * ny;
00458             for (unsigned int index=1; index < nx; index++)
00459               val *= xi;
00460             for (unsigned int index=1; index < ny; index++)
00461               val *= eta;
00462             return val;
00463           }
00464       }
00465 
00466       // d^2()/deta^2
00467     case 2:
00468       {
00469         switch (i)
00470           {
00471             // constants
00472           case 0:
00473 
00474             // linears
00475           case 1:
00476           case 2:
00477             return 0.;
00478 
00479             // quadratics
00480           case 3:
00481           case 4:
00482             return 0.;
00483 
00484           case 5:
00485             return 2.;
00486 
00487             // cubics
00488           case 6:
00489             return 0.;
00490 
00491           case 7:
00492             return 0.;
00493 
00494           case 8:
00495             return 2.*xi;
00496 
00497           case 9:
00498             return 6.*eta;
00499 
00500             // quartics
00501           case 10:
00502           case 11:
00503             return 0.;
00504 
00505           case 12:
00506             return 2.*xi*xi;
00507 
00508           case 13:
00509             return 6.*xi*eta;
00510 
00511           case 14:
00512             return 12.*eta*eta;
00513 
00514           default:
00515             unsigned int o = 0;
00516             for (; i >= (o+1)*(o+2)/2; o++) { }
00517             unsigned int ny = i - (o*(o+1)/2);
00518             unsigned int nx = o - ny;
00519             Real val = ny * (ny - 1);
00520             for (unsigned int index=0; index != nx; index++)
00521               val *= xi;
00522             for (unsigned int index=2; index < ny; index++)
00523               val *= eta;
00524             return val;
00525           }
00526       }
00527 
00528     default:
00529       libmesh_error_msg("Invalid shape function derivative j = " << j);
00530     }
00531 
00532   libmesh_error_msg("We'll never get here!");
00533   return 0.;
00534 
00535 #endif
00536 }
00537 
00538 
00539 
00540 template <>
00541 Real FE<2,MONOMIAL>::shape_second_deriv(const Elem* elem,
00542                                         const Order order,
00543                                         const unsigned int i,
00544                                         const unsigned int j,
00545                                         const Point& p)
00546 {
00547   libmesh_assert(elem);
00548 
00549   // by default call the orientation-independent shape functions
00550   return FE<2,MONOMIAL>::shape_second_deriv(elem->type(), static_cast<Order>(order + elem->p_level()), i, j, p);
00551 }
00552 
00553 
00554 } // namespace libMesh