fe_xyz_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 00026 00027 00028 // Anonymous namespace for persistant variables. 00029 // This allows us to determine when the centroid needs 00030 // to be recalculated. 00031 namespace 00032 { 00033 using namespace libMesh; 00034 00035 static dof_id_type old_elem_id = DofObject::invalid_id; 00036 static Point centroid; 00037 static Point max_distance; 00038 } 00039 00040 00041 namespace libMesh 00042 { 00043 00044 00045 template <> 00046 Real FE<2,XYZ>::shape(const ElemType, 00047 const Order, 00048 const unsigned int, 00049 const Point&) 00050 { 00051 libMesh::err << "XYZ polynomials require the element\n" 00052 << "because the centroid is needed." 00053 << std::endl; 00054 00055 libmesh_error(); 00056 return 0.; 00057 } 00058 00059 00060 00061 template <> 00062 Real FE<2,XYZ>::shape(const Elem* elem, 00063 const Order libmesh_dbg_var(order), 00064 const unsigned int i, 00065 const Point& point_in) 00066 { 00067 #if LIBMESH_DIM > 1 00068 00069 libmesh_assert(elem); 00070 00071 // Only recompute the centroid if the element 00072 // has changed from the last one we computed. 00073 // This avoids repeated centroid calculations 00074 // when called in succession with the same element. 00075 if (elem->id() != old_elem_id) 00076 { 00077 centroid = elem->centroid(); 00078 old_elem_id = elem->id(); 00079 max_distance = Point(0.,0.,0.); 00080 for (unsigned int p = 0; p < elem->n_nodes(); p++) 00081 for (unsigned int d = 0; d < 2; d++) 00082 { 00083 const Real distance = std::abs(centroid(d) - elem->point(p)(d)); 00084 max_distance(d) = std::max(distance, max_distance(d)); 00085 } 00086 } 00087 00088 // Using static globals for old_elem_id, etc. will fail 00089 // horribly with more than one thread. 00090 libmesh_assert_equal_to (libMesh::n_threads(), 1); 00091 00092 const Real x = point_in(0); 00093 const Real y = point_in(1); 00094 const Real xc = centroid(0); 00095 const Real yc = centroid(1); 00096 const Real distx = max_distance(0); 00097 const Real disty = max_distance(1); 00098 const Real dx = (x - xc)/distx; 00099 const Real dy = (y - yc)/disty; 00100 00101 #ifndef NDEBUG 00102 // totalorder is only used in the assertion below, so 00103 // we avoid declaring it when asserts are not active. 00104 const unsigned int totalorder = order + elem->p_level(); 00105 #endif 00106 libmesh_assert_less (i, (totalorder+1)*(totalorder+2)/2); 00107 00108 00109 // monomials. since they are hierarchic we only need one case block. 00110 switch (i) 00111 { 00112 // constant 00113 case 0: 00114 return 1.; 00115 00116 // linear 00117 case 1: 00118 return dx; 00119 00120 case 2: 00121 return dy; 00122 00123 // quadratics 00124 case 3: 00125 return dx*dx; 00126 00127 case 4: 00128 return dx*dy; 00129 00130 case 5: 00131 return dy*dy; 00132 00133 // cubics 00134 case 6: 00135 return dx*dx*dx; 00136 00137 case 7: 00138 return dx*dx*dy; 00139 00140 case 8: 00141 return dx*dy*dy; 00142 00143 case 9: 00144 return dy*dy*dy; 00145 00146 // quartics 00147 case 10: 00148 return dx*dx*dx*dx; 00149 00150 case 11: 00151 return dx*dx*dx*dy; 00152 00153 case 12: 00154 return dx*dx*dy*dy; 00155 00156 case 13: 00157 return dx*dy*dy*dy; 00158 00159 case 14: 00160 return dy*dy*dy*dy; 00161 00162 default: 00163 unsigned int o = 0; 00164 for (; i >= (o+1)*(o+2)/2; o++) { } 00165 unsigned int i2 = i - (o*(o+1)/2); 00166 Real val = 1.; 00167 for (unsigned int index=i2; index != o; index++) 00168 val *= dx; 00169 for (unsigned int index=0; index != i2; index++) 00170 val *= dy; 00171 return val; 00172 } 00173 00174 libmesh_error(); 00175 return 0.; 00176 00177 #endif 00178 } 00179 00180 00181 00182 template <> 00183 Real FE<2,XYZ>::shape_deriv(const ElemType, 00184 const Order, 00185 const unsigned int, 00186 const unsigned int, 00187 const Point&) 00188 { 00189 libMesh::err << "XYZ polynomials require the element\n" 00190 << "because the centroid is needed." 00191 << std::endl; 00192 00193 libmesh_error(); 00194 return 0.; 00195 } 00196 00197 00198 00199 template <> 00200 Real FE<2,XYZ>::shape_deriv(const Elem* elem, 00201 const Order libmesh_dbg_var(order), 00202 const unsigned int i, 00203 const unsigned int j, 00204 const Point& point_in) 00205 { 00206 #if LIBMESH_DIM > 1 00207 00208 00209 libmesh_assert_less (j, 2); 00210 libmesh_assert(elem); 00211 00212 // Only recompute the centroid if the element 00213 // has changed from the last one we computed. 00214 // This avoids repeated centroid calculations 00215 // when called in succession with the same element. 00216 if (elem->id() != old_elem_id) 00217 { 00218 centroid = elem->centroid(); 00219 old_elem_id = elem->id(); 00220 max_distance = Point(0.,0.,0.); 00221 for (unsigned int p = 0; p < elem->n_nodes(); p++) 00222 for (unsigned int d = 0; d < 2; d++) 00223 { 00224 const Real distance = std::abs(centroid(d) - elem->point(p)(d)); 00225 max_distance(d) = std::max(distance, max_distance(d)); 00226 } 00227 } 00228 00229 // Using static globals for old_elem_id, etc. will fail 00230 // horribly with more than one thread. 00231 libmesh_assert_equal_to (libMesh::n_threads(), 1); 00232 00233 const Real x = point_in(0); 00234 const Real y = point_in(1); 00235 const Real xc = centroid(0); 00236 const Real yc = centroid(1); 00237 const Real distx = max_distance(0); 00238 const Real disty = max_distance(1); 00239 const Real dx = (x - xc)/distx; 00240 const Real dy = (y - yc)/disty; 00241 00242 #ifndef NDEBUG 00243 // totalorder is only used in the assertion below, so 00244 // we avoid declaring it when asserts are not active. 00245 const unsigned int totalorder = order + elem->p_level(); 00246 #endif 00247 libmesh_assert_less (i, (totalorder+1)*(totalorder+2)/2); 00248 00249 // monomials. since they are hierarchic we only need one case block. 00250 00251 switch (j) 00252 { 00253 // d()/dx 00254 case 0: 00255 { 00256 switch (i) 00257 { 00258 // constants 00259 case 0: 00260 return 0.; 00261 00262 // linears 00263 case 1: 00264 return 1./distx; 00265 00266 case 2: 00267 return 0.; 00268 00269 // quadratics 00270 case 3: 00271 return 2.*dx/distx; 00272 00273 case 4: 00274 return dy/distx; 00275 00276 case 5: 00277 return 0.; 00278 00279 // cubics 00280 case 6: 00281 return 3.*dx*dx/distx; 00282 00283 case 7: 00284 return 2.*dx*dy/distx; 00285 00286 case 8: 00287 return dy*dy/distx; 00288 00289 case 9: 00290 return 0.; 00291 00292 // quartics 00293 case 10: 00294 return 4.*dx*dx*dx/distx; 00295 00296 case 11: 00297 return 3.*dx*dx*dy/distx; 00298 00299 case 12: 00300 return 2.*dx*dy*dy/distx; 00301 00302 case 13: 00303 return dy*dy*dy/distx; 00304 00305 case 14: 00306 return 0.; 00307 00308 default: 00309 unsigned int o = 0; 00310 for (; i >= (o+1)*(o+2)/2; o++) { } 00311 unsigned int i2 = i - (o*(o+1)/2); 00312 Real val = o - i2; 00313 for (unsigned int index=i2+1; index < o; index++) 00314 val *= dx; 00315 for (unsigned int index=0; index != i2; index++) 00316 val *= dy; 00317 return val/distx; 00318 } 00319 } 00320 00321 00322 // d()/dy 00323 case 1: 00324 { 00325 switch (i) 00326 { 00327 // constants 00328 case 0: 00329 return 0.; 00330 00331 // linears 00332 case 1: 00333 return 0.; 00334 00335 case 2: 00336 return 1./disty; 00337 00338 // quadratics 00339 case 3: 00340 return 0.; 00341 00342 case 4: 00343 return dx/disty; 00344 00345 case 5: 00346 return 2.*dy/disty; 00347 00348 // cubics 00349 case 6: 00350 return 0.; 00351 00352 case 7: 00353 return dx*dx/disty; 00354 00355 case 8: 00356 return 2.*dx*dy/disty; 00357 00358 case 9: 00359 return 3.*dy*dy/disty; 00360 00361 // quartics 00362 case 10: 00363 return 0.; 00364 00365 case 11: 00366 return dx*dx*dx/disty; 00367 00368 case 12: 00369 return 2.*dx*dx*dy/disty; 00370 00371 case 13: 00372 return 3.*dx*dy*dy/disty; 00373 00374 case 14: 00375 return 4.*dy*dy*dy/disty; 00376 00377 default: 00378 unsigned int o = 0; 00379 for (; i >= (o+1)*(o+2)/2; o++) { } 00380 unsigned int i2 = i - (o*(o+1)/2); 00381 Real val = i2; 00382 for (unsigned int index=i2; index != o; index++) 00383 val *= dx; 00384 for (unsigned int index=1; index <= i2; index++) 00385 val *= dy; 00386 return val/disty; 00387 } 00388 } 00389 00390 00391 default: 00392 libmesh_error(); 00393 } 00394 00395 libmesh_error(); 00396 return 0.; 00397 00398 #endif 00399 } 00400 00401 00402 00403 template <> 00404 Real FE<2,XYZ>::shape_second_deriv(const ElemType, 00405 const Order, 00406 const unsigned int, 00407 const unsigned int, 00408 const Point&) 00409 { 00410 libMesh::err << "XYZ polynomials require the element\n" 00411 << "because the centroid is needed." 00412 << std::endl; 00413 00414 libmesh_error(); 00415 return 0.; 00416 } 00417 00418 00419 00420 template <> 00421 Real FE<2,XYZ>::shape_second_deriv(const Elem* elem, 00422 const Order libmesh_dbg_var(order), 00423 const unsigned int i, 00424 const unsigned int j, 00425 const Point& point_in) 00426 { 00427 #if LIBMESH_DIM > 1 00428 00429 libmesh_assert_less_equal (j, 2); 00430 libmesh_assert(elem); 00431 00432 // Only recompute the centroid if the element 00433 // has changed from the last one we computed. 00434 // This avoids repeated centroid calculations 00435 // when called in succession with the same element. 00436 if (elem->id() != old_elem_id) 00437 { 00438 centroid = elem->centroid(); 00439 old_elem_id = elem->id(); 00440 max_distance = Point(0.,0.,0.); 00441 for (unsigned int p = 0; p < elem->n_nodes(); p++) 00442 for (unsigned int d = 0; d < 2; d++) 00443 { 00444 const Real distance = std::abs(centroid(d) - elem->point(p)(d)); 00445 max_distance(d) = std::max(distance, max_distance(d)); 00446 } 00447 } 00448 00449 // Using static globals for old_elem_id, etc. will fail 00450 // horribly with more than one thread. 00451 libmesh_assert_equal_to (libMesh::n_threads(), 1); 00452 00453 const Real x = point_in(0); 00454 const Real y = point_in(1); 00455 const Real xc = centroid(0); 00456 const Real yc = centroid(1); 00457 const Real distx = max_distance(0); 00458 const Real disty = max_distance(1); 00459 const Real dx = (x - xc)/distx; 00460 const Real dy = (y - yc)/disty; 00461 const Real dist2x = pow(distx,2.); 00462 const Real dist2y = pow(disty,2.); 00463 const Real distxy = distx * disty; 00464 00465 #ifndef NDEBUG 00466 // totalorder is only used in the assertion below, so 00467 // we avoid declaring it when asserts are not active. 00468 const unsigned int totalorder = order + elem->p_level(); 00469 #endif 00470 libmesh_assert_less (i, (totalorder+1)*(totalorder+2)/2); 00471 00472 // monomials. since they are hierarchic we only need one case block. 00473 00474 switch (j) 00475 { 00476 // d^2()/dx^2 00477 case 0: 00478 { 00479 switch (i) 00480 { 00481 // constants 00482 case 0: 00483 // linears 00484 case 1: 00485 case 2: 00486 return 0.; 00487 00488 // quadratics 00489 case 3: 00490 return 2./dist2x; 00491 00492 case 4: 00493 case 5: 00494 return 0.; 00495 00496 // cubics 00497 case 6: 00498 return 6.*dx/dist2x; 00499 00500 case 7: 00501 return 2.*dy/dist2x; 00502 00503 case 8: 00504 case 9: 00505 return 0.; 00506 00507 // quartics 00508 case 10: 00509 return 12.*dx*dx/dist2x; 00510 00511 case 11: 00512 return 6.*dx*dy/dist2x; 00513 00514 case 12: 00515 return 2.*dy*dy/dist2x; 00516 00517 case 13: 00518 case 14: 00519 return 0.; 00520 00521 default: 00522 unsigned int o = 0; 00523 for (; i >= (o+1)*(o+2)/2; o++) { } 00524 unsigned int i2 = i - (o*(o+1)/2); 00525 Real val = (o - i2) * (o - i2 - 1); 00526 for (unsigned int index=i2+2; index < o; index++) 00527 val *= dx; 00528 for (unsigned int index=0; index != i2; index++) 00529 val *= dy; 00530 return val/dist2x; 00531 } 00532 } 00533 00534 // d^2()/dxdy 00535 case 1: 00536 { 00537 switch (i) 00538 { 00539 // constants 00540 case 0: 00541 00542 // linears 00543 case 1: 00544 case 2: 00545 return 0.; 00546 00547 // quadratics 00548 case 3: 00549 return 0.; 00550 00551 case 4: 00552 return 1./distxy; 00553 00554 case 5: 00555 return 0.; 00556 00557 // cubics 00558 case 6: 00559 return 0.; 00560 case 7: 00561 return 2.*dx/distxy; 00562 00563 case 8: 00564 return 2.*dy/distxy; 00565 00566 case 9: 00567 return 0.; 00568 00569 // quartics 00570 case 10: 00571 return 0.; 00572 00573 case 11: 00574 return 3.*dx*dx/distxy; 00575 00576 case 12: 00577 return 4.*dx*dy/distxy; 00578 00579 case 13: 00580 return 3.*dy*dy/distxy; 00581 00582 case 14: 00583 return 0.; 00584 00585 default: 00586 unsigned int o = 0; 00587 for (; i >= (o+1)*(o+2)/2; o++) { } 00588 unsigned int i2 = i - (o*(o+1)/2); 00589 Real val = (o - i2) * i2; 00590 for (unsigned int index=i2+1; index < o; index++) 00591 val *= dx; 00592 for (unsigned int index=1; index < i2; index++) 00593 val *= dy; 00594 return val/distxy; 00595 } 00596 } 00597 00598 // d^2()/dy^2 00599 case 2: 00600 { 00601 switch (i) 00602 { 00603 // constants 00604 case 0: 00605 00606 // linears 00607 case 1: 00608 case 2: 00609 return 0.; 00610 00611 // quadratics 00612 case 3: 00613 case 4: 00614 return 0.; 00615 00616 case 5: 00617 return 2./dist2y; 00618 00619 // cubics 00620 case 6: 00621 return 0.; 00622 00623 case 7: 00624 return 0.; 00625 00626 case 8: 00627 return 2.*dx/dist2y; 00628 00629 case 9: 00630 return 6.*dy/dist2y; 00631 00632 // quartics 00633 case 10: 00634 case 11: 00635 return 0.; 00636 00637 case 12: 00638 return 2.*dx*dx/dist2y; 00639 00640 case 13: 00641 return 6.*dx*dy/dist2y; 00642 00643 case 14: 00644 return 12.*dy*dy/dist2y; 00645 00646 default: 00647 unsigned int o = 0; 00648 for (; i >= (o+1)*(o+2)/2; o++) { } 00649 unsigned int i2 = i - (o*(o+1)/2); 00650 Real val = i2 * (i2 - 1); 00651 for (unsigned int index=i2; index != o; index++) 00652 val *= dx; 00653 for (unsigned int index=2; index < i2; index++) 00654 val *= dy; 00655 return val/dist2y; 00656 } 00657 } 00658 } 00659 00660 libmesh_error(); 00661 return 0.; 00662 00663 #endif 00664 } 00665 00666 } // namespace libMesh
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:47 UTC
Hosted By: