fe_hermite_shape_1D.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 #include "libmesh/utility.h" 00025 00026 #ifdef LIBMESH_HAVE_TBB_API 00027 #include "tbb/enumerable_thread_specific.h" 00028 #endif 00029 00030 00031 // Anonymous namespace for persistant variables. 00032 // This allows us to cache the global-to-local mapping transformation 00033 // This caching is turned off when TBB is enabled... 00034 namespace 00035 { 00036 using namespace libMesh; 00037 00038 #ifndef LIBMESH_HAVE_TBB_API 00039 static dof_id_type old_elem_id = libMesh::invalid_uint; 00040 // Coefficient naming: d(1)d(2n) is the coefficient of the 00041 // global shape function corresponding to value 1 in terms of the 00042 // local shape function corresponding to normal derivative 2 00043 static Real d1xd1x, d2xd2x; 00044 #else //LIBMESH_HAVE_TBB_API 00045 static tbb::enumerable_thread_specific< dof_id_type > old_elem_id_tls; 00046 static tbb::enumerable_thread_specific< Real > d1xd1x_tls; 00047 static tbb::enumerable_thread_specific< Real > d2xd2x_tls; 00048 #endif //LIBMESH_HAVE_TBB_API 00049 00050 // Compute the static coefficients for an element 00051 void hermite_compute_coefs(const Elem* elem) 00052 { 00053 #ifndef LIBMESH_HAVE_TBB_API 00054 // Coefficients are cached from old elements 00055 // ... except that we can't be sure that a renumbering didn't 00056 // cause the first element id in a new assembly to match the 00057 // last element id in an old assembly ... 00058 // if (elem->id() == old_elem_id) 00059 // return; 00060 00061 old_elem_id = elem->id(); 00062 #else 00063 bool old_elem_id_exists = false; 00064 dof_id_type & old_elem_id = old_elem_id_tls.local(old_elem_id_exists); 00065 00066 if(!old_elem_id_exists) 00067 old_elem_id = libMesh::invalid_uint; 00068 00069 // Coefficients are cached from old elements 00070 // ... except that we can't be sure that a renumbering didn't 00071 // cause the first element id in a new assembly to match the 00072 // last element id in an old assembly ... 00073 // if (elem->id() == old_elem_id) 00074 // return; 00075 00076 old_elem_id = elem->id(); 00077 00078 Real & d1xd1x = d1xd1x_tls.local(); 00079 Real & d2xd2x = d2xd2x_tls.local(); 00080 #endif //LIBMESH_HAVE_TBB_API 00081 00082 const Order mapping_order (elem->default_order()); 00083 const ElemType mapping_elem_type (elem->type()); 00084 const int n_mapping_shape_functions = 00085 FE<1,LAGRANGE>::n_shape_functions(mapping_elem_type, 00086 mapping_order); 00087 00088 // Degrees of freedom are at vertices and edge midpoints 00089 std::vector<Point> dofpt; 00090 dofpt.push_back(Point(-1)); 00091 dofpt.push_back(Point(1)); 00092 00093 // Mapping functions - first derivatives at each dofpt 00094 std::vector<Real> dxdxi(2); 00095 std::vector<Real> dxidx(2); 00096 00097 for (int p = 0; p != 2; ++p) 00098 { 00099 dxdxi[p] = 0; 00100 for (int i = 0; i != n_mapping_shape_functions; ++i) 00101 { 00102 const Real ddxi = FE<1,LAGRANGE>::shape_deriv 00103 (mapping_elem_type, mapping_order, i, 0, dofpt[p]); 00104 dxdxi[p] += elem->point(i)(0) * ddxi; 00105 } 00106 } 00107 00108 // Calculate derivative scaling factors 00109 00110 d1xd1x = dxdxi[0]; 00111 d2xd2x = dxdxi[1]; 00112 } 00113 00114 00115 } // end anonymous namespace 00116 00117 00118 namespace libMesh 00119 { 00120 00121 00122 template<> 00123 Real FEHermite<1>::hermite_raw_shape_second_deriv 00124 (const unsigned int i, const Real xi) 00125 { 00126 using Utility::pow; 00127 00128 switch (i) 00129 { 00130 case 0: 00131 return 1.5 * xi; 00132 case 1: 00133 return -1.5 * xi; 00134 case 2: 00135 return 0.5 * (-1. + 3.*xi); 00136 case 3: 00137 return 0.5 * (1. + 3.*xi); 00138 case 4: 00139 return (8.*xi*xi + 4.*(xi*xi-1.))/24.; 00140 case 5: 00141 return (8.*xi*xi*xi + 12.*xi*(xi*xi-1.))/120.; 00142 // case 6: 00143 // return (8.*pow<4>(xi) + 20.*xi*xi*(xi*xi-1.) + 00144 // 2.*(xi*xi-1)*(xi*xi-1))/720.; 00145 default: 00146 Real denominator = 720., xipower = 1.; 00147 for (unsigned n=6; n != i; ++n) 00148 { 00149 xipower *= xi; 00150 denominator *= (n+1); 00151 } 00152 return (8.*pow<4>(xi)*xipower + 00153 (8.*(i-4)+4.)*xi*xi*xipower*(xi*xi-1.) + 00154 (i-4)*(i-5)*xipower*(xi*xi-1.)*(xi*xi-1.))/denominator; 00155 } 00156 00157 libmesh_error(); 00158 return 0.; 00159 } 00160 00161 00162 00163 template<> 00164 Real FEHermite<1>::hermite_raw_shape_deriv 00165 (const unsigned int i, const Real xi) 00166 { 00167 switch (i) 00168 { 00169 case 0: 00170 return 0.75 * (-1. + xi*xi); 00171 case 1: 00172 return 0.75 * (1. - xi*xi); 00173 case 2: 00174 return 0.25 * (-1. - 2.*xi + 3.*xi*xi); 00175 case 3: 00176 return 0.25 * (-1. + 2.*xi + 3.*xi*xi); 00177 case 4: 00178 return 4.*xi * (xi*xi-1.)/24.; 00179 case 5: 00180 return (4*xi*xi*(xi*xi-1.) + (xi*xi-1.)*(xi*xi-1.))/120.; 00181 // case 6: 00182 // return (4*xi*xi*xi*(xi*xi-1.) + 2*xi*(xi*xi-1.)*(xi*xi-1.))/720.; 00183 default: 00184 Real denominator = 720., xipower = 1.; 00185 for (unsigned n=6; n != i; ++n) 00186 { 00187 xipower *= xi; 00188 denominator *= (n+1); 00189 } 00190 return (4*xi*xi*xi*xipower*(xi*xi-1.) + 00191 (i-4)*xi*xipower*(xi*xi-1.)*(xi*xi-1.))/denominator; 00192 } 00193 00194 libmesh_error(); 00195 return 0.; 00196 } 00197 00198 template<> 00199 Real FEHermite<1>::hermite_raw_shape 00200 (const unsigned int i, const Real xi) 00201 { 00202 switch (i) 00203 { 00204 case 0: 00205 return 0.25 * (2. - 3.*xi + xi*xi*xi); 00206 case 1: 00207 return 0.25 * (2. + 3.*xi - xi*xi*xi); 00208 case 2: 00209 return 0.25 * (1. - xi - xi*xi + xi*xi*xi); 00210 case 3: 00211 return 0.25 * (-1. - xi + xi*xi + xi*xi*xi); 00212 // All high order terms have the form x^(p-4)(x^2-1)^2/p! 00213 case 4: 00214 return (xi*xi-1.) * (xi*xi-1.)/24.; 00215 case 5: 00216 return xi * (xi*xi-1.) * (xi*xi-1.)/120.; 00217 // case 6: 00218 // return xi*xi * (xi*xi-1.) * (xi*xi-1.)/720.; 00219 default: 00220 Real denominator = 720., xipower = 1.; 00221 for (unsigned n=6; n != i; ++n) 00222 { 00223 xipower *= xi; 00224 denominator *= (n+1); 00225 } 00226 return (xi*xi*xipower*(xi*xi-1.)*(xi*xi-1.))/denominator; 00227 00228 } 00229 00230 libmesh_error(); 00231 return 0.; 00232 } 00233 00234 00235 template <> 00236 Real FE<1,HERMITE>::shape(const ElemType, 00237 const Order, 00238 const unsigned int, 00239 const Point&) 00240 { 00241 libMesh::err << "Hermite elements require the real element\n" 00242 << "to construct gradient-based degrees of freedom." 00243 << std::endl; 00244 00245 libmesh_error(); 00246 return 0.; 00247 } 00248 00249 00250 00251 template <> 00252 Real FE<1,HERMITE>::shape(const Elem* elem, 00253 const Order order, 00254 const unsigned int i, 00255 const Point& p) 00256 { 00257 libmesh_assert(elem); 00258 00259 hermite_compute_coefs(elem); 00260 00261 #ifdef LIBMESH_HAVE_TBB_API 00262 Real & d1xd1x = d1xd1x_tls.local(); 00263 Real & d2xd2x = d2xd2x_tls.local(); 00264 #endif //LIBMESH_HAVE_TBB_API 00265 00266 const ElemType type = elem->type(); 00267 00268 const Order totalorder = static_cast<Order>(order + elem->p_level()); 00269 00270 switch (totalorder) 00271 { 00272 // Hermite cubic shape functions 00273 case THIRD: 00274 { 00275 switch (type) 00276 { 00277 // C1 functions on the C1 cubic edge 00278 case EDGE2: 00279 case EDGE3: 00280 { 00281 libmesh_assert_less (i, 4); 00282 00283 switch (i) 00284 { 00285 case 0: 00286 return FEHermite<1>::hermite_raw_shape(0, p(0)); 00287 case 1: 00288 return d1xd1x * FEHermite<1>::hermite_raw_shape(2, p(0)); 00289 case 2: 00290 return FEHermite<1>::hermite_raw_shape(1, p(0)); 00291 case 3: 00292 return d2xd2x * FEHermite<1>::hermite_raw_shape(3, p(0)); 00293 default: 00294 return FEHermite<1>::hermite_raw_shape(i, p(0)); 00295 } 00296 } 00297 default: 00298 libMesh::err << "ERROR: Unsupported element type!" << std::endl; 00299 libmesh_error(); 00300 } 00301 } 00302 // by default throw an error 00303 default: 00304 libMesh::err << "ERROR: Unsupported polynomial order!" << std::endl; 00305 libmesh_error(); 00306 } 00307 00308 libmesh_error(); 00309 return 0.; 00310 } 00311 00312 00313 00314 template <> 00315 Real FE<1,HERMITE>::shape_deriv(const ElemType, 00316 const Order, 00317 const unsigned int, 00318 const unsigned int, 00319 const Point&) 00320 { 00321 libMesh::err << "Hermite elements require the real element\n" 00322 << "to construct gradient-based degrees of freedom." 00323 << std::endl; 00324 00325 libmesh_error(); 00326 return 0.; 00327 } 00328 00329 00330 00331 template <> 00332 Real FE<1,HERMITE>::shape_deriv(const Elem* elem, 00333 const Order order, 00334 const unsigned int i, 00335 const unsigned int, 00336 const Point& p) 00337 { 00338 libmesh_assert(elem); 00339 00340 hermite_compute_coefs(elem); 00341 00342 #ifdef LIBMESH_HAVE_TBB_API 00343 Real & d1xd1x = d1xd1x_tls.local(); 00344 Real & d2xd2x = d2xd2x_tls.local(); 00345 #endif //LIBMESH_HAVE_TBB_API 00346 00347 const ElemType type = elem->type(); 00348 00349 const Order totalorder = static_cast<Order>(order + elem->p_level()); 00350 00351 switch (totalorder) 00352 { 00353 // Hermite cubic shape functions 00354 case THIRD: 00355 { 00356 switch (type) 00357 { 00358 // C1 functions on the C1 cubic edge 00359 case EDGE2: 00360 case EDGE3: 00361 { 00362 switch (i) 00363 { 00364 case 0: 00365 return FEHermite<1>::hermite_raw_shape_deriv(0, p(0)); 00366 case 1: 00367 return d1xd1x * FEHermite<1>::hermite_raw_shape_deriv(2, p(0)); 00368 case 2: 00369 return FEHermite<1>::hermite_raw_shape_deriv(1, p(0)); 00370 case 3: 00371 return d2xd2x * FEHermite<1>::hermite_raw_shape_deriv(3, p(0)); 00372 default: 00373 return FEHermite<1>::hermite_raw_shape_deriv(i, p(0)); 00374 } 00375 } 00376 default: 00377 libMesh::err << "ERROR: Unsupported element type!" << std::endl; 00378 libmesh_error(); 00379 } 00380 } 00381 // by default throw an error 00382 default: 00383 libMesh::err << "ERROR: Unsupported polynomial order!" << std::endl; 00384 libmesh_error(); 00385 } 00386 00387 libmesh_error(); 00388 return 0.; 00389 } 00390 00391 00392 00393 template <> 00394 Real FE<1,HERMITE>::shape_second_deriv(const Elem* elem, 00395 const Order order, 00396 const unsigned int i, 00397 const unsigned int, 00398 const Point& p) 00399 { 00400 libmesh_assert(elem); 00401 00402 hermite_compute_coefs(elem); 00403 00404 #ifdef LIBMESH_HAVE_TBB_API 00405 Real & d1xd1x = d1xd1x_tls.local(); 00406 Real & d2xd2x = d2xd2x_tls.local(); 00407 #endif //LIBMESH_HAVE_TBB_API 00408 00409 const ElemType type = elem->type(); 00410 00411 const Order totalorder = static_cast<Order>(order + elem->p_level()); 00412 00413 switch (totalorder) 00414 { 00415 // Hermite cubic shape functions 00416 case THIRD: 00417 { 00418 switch (type) 00419 { 00420 // C1 functions on the C1 cubic edge 00421 case EDGE2: 00422 case EDGE3: 00423 { 00424 switch (i) 00425 { 00426 case 0: 00427 return FEHermite<1>::hermite_raw_shape_second_deriv(0, p(0)); 00428 case 1: 00429 return d1xd1x * FEHermite<1>::hermite_raw_shape_second_deriv(2, p(0)); 00430 case 2: 00431 return FEHermite<1>::hermite_raw_shape_second_deriv(1, p(0)); 00432 case 3: 00433 return d2xd2x * FEHermite<1>::hermite_raw_shape_second_deriv(3, p(0)); 00434 default: 00435 return FEHermite<1>::hermite_raw_shape_second_deriv(i, p(0)); 00436 } 00437 } 00438 default: 00439 libMesh::err << "ERROR: Unsupported element type!" << std::endl; 00440 libmesh_error(); 00441 } 00442 } 00443 // by default throw an error 00444 default: 00445 libMesh::err << "ERROR: Unsupported polynomial order!" << std::endl; 00446 libmesh_error(); 00447 } 00448 00449 libmesh_error(); 00450 return 0.; 00451 } 00452 00453 } // namespace libMesh
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:46 UTC
Hosted By: