fe_monomial_shape_2D.C
Go to the documentation of this file.00001 // The libMesh Finite Element Library. 00002 // Copyright (C) 2002-2012 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(); 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 00294 libmesh_error(); 00295 return 0.; 00296 00297 #endif 00298 } 00299 00300 00301 00302 template <> 00303 Real FE<2,MONOMIAL>::shape_deriv(const Elem* elem, 00304 const Order order, 00305 const unsigned int i, 00306 const unsigned int j, 00307 const Point& p) 00308 { 00309 libmesh_assert(elem); 00310 00311 // by default call the orientation-independent shape functions 00312 return FE<2,MONOMIAL>::shape_deriv(elem->type(), static_cast<Order>(order + elem->p_level()), i, j, p); 00313 } 00314 00315 00316 00317 template <> 00318 Real FE<2,MONOMIAL>::shape_second_deriv(const ElemType, 00319 const Order libmesh_dbg_var(order), 00320 const unsigned int i, 00321 const unsigned int j, 00322 const Point& p) 00323 { 00324 #if LIBMESH_DIM > 1 00325 00326 00327 libmesh_assert_less_equal (j, 2); 00328 00329 libmesh_assert_less (i, (static_cast<unsigned int>(order)+1)* 00330 (static_cast<unsigned int>(order)+2)/2); 00331 00332 const Real xi = p(0); 00333 const Real eta = p(1); 00334 00335 // monomials. since they are hierarchic we only need one case block. 00336 00337 switch (j) 00338 { 00339 // d^2()/dxi^2 00340 case 0: 00341 { 00342 switch (i) 00343 { 00344 // constants 00345 case 0: 00346 // linears 00347 case 1: 00348 case 2: 00349 return 0.; 00350 00351 // quadratics 00352 case 3: 00353 return 2.; 00354 00355 case 4: 00356 case 5: 00357 return 0.; 00358 00359 // cubics 00360 case 6: 00361 return 6.*xi; 00362 00363 case 7: 00364 return 2.*eta; 00365 00366 case 8: 00367 case 9: 00368 return 0.; 00369 00370 // quartics 00371 case 10: 00372 return 12.*xi*xi; 00373 00374 case 11: 00375 return 6.*xi*eta; 00376 00377 case 12: 00378 return 2.*eta*eta; 00379 00380 case 13: 00381 case 14: 00382 return 0.; 00383 00384 default: 00385 unsigned int o = 0; 00386 for (; i >= (o+1)*(o+2)/2; o++) { } 00387 unsigned int ny = i - (o*(o+1)/2); 00388 unsigned int nx = o - ny; 00389 Real val = nx * (nx - 1); 00390 for (unsigned int index=2; index < nx; index++) 00391 val *= xi; 00392 for (unsigned int index=0; index != ny; index++) 00393 val *= eta; 00394 return val; 00395 } 00396 } 00397 00398 // d^2()/dxideta 00399 case 1: 00400 { 00401 switch (i) 00402 { 00403 // constants 00404 case 0: 00405 00406 // linears 00407 case 1: 00408 case 2: 00409 return 0.; 00410 00411 // quadratics 00412 case 3: 00413 return 0.; 00414 00415 case 4: 00416 return 1.; 00417 00418 case 5: 00419 return 0.; 00420 00421 // cubics 00422 case 6: 00423 return 0.; 00424 case 7: 00425 return 2.*xi; 00426 00427 case 8: 00428 return 2.*eta; 00429 00430 case 9: 00431 return 0.; 00432 00433 // quartics 00434 case 10: 00435 return 0.; 00436 00437 case 11: 00438 return 3.*xi*xi; 00439 00440 case 12: 00441 return 4.*xi*eta; 00442 00443 case 13: 00444 return 3.*eta*eta; 00445 00446 case 14: 00447 return 0.; 00448 00449 default: 00450 unsigned int o = 0; 00451 for (; i >= (o+1)*(o+2)/2; o++) { } 00452 unsigned int ny = i - (o*(o+1)/2); 00453 unsigned int nx = o - ny; 00454 Real val = nx * ny; 00455 for (unsigned int index=1; index < nx; index++) 00456 val *= xi; 00457 for (unsigned int index=1; index < ny; index++) 00458 val *= eta; 00459 return val; 00460 } 00461 } 00462 00463 // d^2()/deta^2 00464 case 2: 00465 { 00466 switch (i) 00467 { 00468 // constants 00469 case 0: 00470 00471 // linears 00472 case 1: 00473 case 2: 00474 return 0.; 00475 00476 // quadratics 00477 case 3: 00478 case 4: 00479 return 0.; 00480 00481 case 5: 00482 return 2.; 00483 00484 // cubics 00485 case 6: 00486 return 0.; 00487 00488 case 7: 00489 return 0.; 00490 00491 case 8: 00492 return 2.*xi; 00493 00494 case 9: 00495 return 6.*eta; 00496 00497 // quartics 00498 case 10: 00499 case 11: 00500 return 0.; 00501 00502 case 12: 00503 return 2.*xi*xi; 00504 00505 case 13: 00506 return 6.*xi*eta; 00507 00508 case 14: 00509 return 12.*eta*eta; 00510 00511 default: 00512 unsigned int o = 0; 00513 for (; i >= (o+1)*(o+2)/2; o++) { } 00514 unsigned int ny = i - (o*(o+1)/2); 00515 unsigned int nx = o - ny; 00516 Real val = ny * (ny - 1); 00517 for (unsigned int index=0; index != nx; index++) 00518 val *= xi; 00519 for (unsigned int index=2; index < ny; index++) 00520 val *= eta; 00521 return val; 00522 } 00523 } 00524 } 00525 00526 libmesh_error(); 00527 return 0.; 00528 00529 #endif 00530 } 00531 00532 00533 00534 template <> 00535 Real FE<2,MONOMIAL>::shape_second_deriv(const Elem* elem, 00536 const Order order, 00537 const unsigned int i, 00538 const unsigned int j, 00539 const Point& p) 00540 { 00541 libmesh_assert(elem); 00542 00543 // by default call the orientation-independent shape functions 00544 return FE<2,MONOMIAL>::shape_second_deriv(elem->type(), static_cast<Order>(order + elem->p_level()), i, j, p); 00545 } 00546 00547 00548 } // namespace libMesh
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:46 UTC
Hosted By: