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