$extrastylesheet
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 #include "libmesh/utility.h" 00025 00026 namespace 00027 { 00028 using namespace libMesh; 00029 00030 // Compute the static coefficients for an element 00031 void hermite_compute_coefs(const Elem* elem, Real & d1xd1x, Real & d2xd2x) 00032 { 00033 const Order mapping_order (elem->default_order()); 00034 const ElemType mapping_elem_type (elem->type()); 00035 const int n_mapping_shape_functions = 00036 FE<1,LAGRANGE>::n_shape_functions(mapping_elem_type, 00037 mapping_order); 00038 00039 // Degrees of freedom are at vertices and edge midpoints 00040 std::vector<Point> dofpt; 00041 dofpt.push_back(Point(-1)); 00042 dofpt.push_back(Point(1)); 00043 00044 // Mapping functions - first derivatives at each dofpt 00045 std::vector<Real> dxdxi(2); 00046 std::vector<Real> dxidx(2); 00047 00048 for (int p = 0; p != 2; ++p) 00049 { 00050 dxdxi[p] = 0; 00051 for (int i = 0; i != n_mapping_shape_functions; ++i) 00052 { 00053 const Real ddxi = FE<1,LAGRANGE>::shape_deriv 00054 (mapping_elem_type, mapping_order, i, 0, dofpt[p]); 00055 dxdxi[p] += elem->point(i)(0) * ddxi; 00056 } 00057 } 00058 00059 // Calculate derivative scaling factors 00060 00061 d1xd1x = dxdxi[0]; 00062 d2xd2x = dxdxi[1]; 00063 } 00064 00065 00066 } // end anonymous namespace 00067 00068 00069 namespace libMesh 00070 { 00071 00072 00073 template<> 00074 Real FEHermite<1>::hermite_raw_shape_second_deriv 00075 (const unsigned int i, const Real xi) 00076 { 00077 using Utility::pow; 00078 00079 switch (i) 00080 { 00081 case 0: 00082 return 1.5 * xi; 00083 case 1: 00084 return -1.5 * xi; 00085 case 2: 00086 return 0.5 * (-1. + 3.*xi); 00087 case 3: 00088 return 0.5 * (1. + 3.*xi); 00089 case 4: 00090 return (8.*xi*xi + 4.*(xi*xi-1.))/24.; 00091 case 5: 00092 return (8.*xi*xi*xi + 12.*xi*(xi*xi-1.))/120.; 00093 // case 6: 00094 // return (8.*pow<4>(xi) + 20.*xi*xi*(xi*xi-1.) + 00095 // 2.*(xi*xi-1)*(xi*xi-1))/720.; 00096 default: 00097 Real denominator = 720., xipower = 1.; 00098 for (unsigned n=6; n != i; ++n) 00099 { 00100 xipower *= xi; 00101 denominator *= (n+1); 00102 } 00103 return (8.*pow<4>(xi)*xipower + 00104 (8.*(i-4)+4.)*xi*xi*xipower*(xi*xi-1.) + 00105 (i-4)*(i-5)*xipower*(xi*xi-1.)*(xi*xi-1.))/denominator; 00106 } 00107 00108 libmesh_error_msg("We'll never get here!"); 00109 return 0.; 00110 } 00111 00112 00113 00114 template<> 00115 Real FEHermite<1>::hermite_raw_shape_deriv 00116 (const unsigned int i, const Real xi) 00117 { 00118 switch (i) 00119 { 00120 case 0: 00121 return 0.75 * (-1. + xi*xi); 00122 case 1: 00123 return 0.75 * (1. - xi*xi); 00124 case 2: 00125 return 0.25 * (-1. - 2.*xi + 3.*xi*xi); 00126 case 3: 00127 return 0.25 * (-1. + 2.*xi + 3.*xi*xi); 00128 case 4: 00129 return 4.*xi * (xi*xi-1.)/24.; 00130 case 5: 00131 return (4*xi*xi*(xi*xi-1.) + (xi*xi-1.)*(xi*xi-1.))/120.; 00132 // case 6: 00133 // return (4*xi*xi*xi*(xi*xi-1.) + 2*xi*(xi*xi-1.)*(xi*xi-1.))/720.; 00134 default: 00135 Real denominator = 720., xipower = 1.; 00136 for (unsigned n=6; n != i; ++n) 00137 { 00138 xipower *= xi; 00139 denominator *= (n+1); 00140 } 00141 return (4*xi*xi*xi*xipower*(xi*xi-1.) + 00142 (i-4)*xi*xipower*(xi*xi-1.)*(xi*xi-1.))/denominator; 00143 } 00144 00145 libmesh_error_msg("We'll never get here!"); 00146 return 0.; 00147 } 00148 00149 template<> 00150 Real FEHermite<1>::hermite_raw_shape 00151 (const unsigned int i, const Real xi) 00152 { 00153 switch (i) 00154 { 00155 case 0: 00156 return 0.25 * (2. - 3.*xi + xi*xi*xi); 00157 case 1: 00158 return 0.25 * (2. + 3.*xi - xi*xi*xi); 00159 case 2: 00160 return 0.25 * (1. - xi - xi*xi + xi*xi*xi); 00161 case 3: 00162 return 0.25 * (-1. - xi + xi*xi + xi*xi*xi); 00163 // All high order terms have the form x^(p-4)(x^2-1)^2/p! 00164 case 4: 00165 return (xi*xi-1.) * (xi*xi-1.)/24.; 00166 case 5: 00167 return xi * (xi*xi-1.) * (xi*xi-1.)/120.; 00168 // case 6: 00169 // return xi*xi * (xi*xi-1.) * (xi*xi-1.)/720.; 00170 default: 00171 Real denominator = 720., xipower = 1.; 00172 for (unsigned n=6; n != i; ++n) 00173 { 00174 xipower *= xi; 00175 denominator *= (n+1); 00176 } 00177 return (xi*xi*xipower*(xi*xi-1.)*(xi*xi-1.))/denominator; 00178 00179 } 00180 00181 libmesh_error_msg("We'll never get here!"); 00182 return 0.; 00183 } 00184 00185 00186 template <> 00187 Real FE<1,HERMITE>::shape(const ElemType, 00188 const Order, 00189 const unsigned int, 00190 const Point&) 00191 { 00192 libmesh_error_msg("Hermite elements require the real element \nto construct gradient-based degrees of freedom."); 00193 return 0.; 00194 } 00195 00196 00197 00198 template <> 00199 Real FE<1,HERMITE>::shape(const Elem* elem, 00200 const Order order, 00201 const unsigned int i, 00202 const Point& p) 00203 { 00204 libmesh_assert(elem); 00205 00206 // Coefficient naming: d(1)d(2n) is the coefficient of the 00207 // global shape function corresponding to value 1 in terms of the 00208 // local shape function corresponding to normal derivative 2 00209 Real d1xd1x, d2xd2x; 00210 00211 hermite_compute_coefs(elem, d1xd1x, d2xd2x); 00212 00213 const ElemType type = elem->type(); 00214 00215 const Order totalorder = static_cast<Order>(order + elem->p_level()); 00216 00217 switch (totalorder) 00218 { 00219 // Hermite cubic shape functions 00220 case THIRD: 00221 { 00222 switch (type) 00223 { 00224 // C1 functions on the C1 cubic edge 00225 case EDGE2: 00226 case EDGE3: 00227 { 00228 libmesh_assert_less (i, 4); 00229 00230 switch (i) 00231 { 00232 case 0: 00233 return FEHermite<1>::hermite_raw_shape(0, p(0)); 00234 case 1: 00235 return d1xd1x * FEHermite<1>::hermite_raw_shape(2, p(0)); 00236 case 2: 00237 return FEHermite<1>::hermite_raw_shape(1, p(0)); 00238 case 3: 00239 return d2xd2x * FEHermite<1>::hermite_raw_shape(3, p(0)); 00240 default: 00241 return FEHermite<1>::hermite_raw_shape(i, p(0)); 00242 } 00243 } 00244 default: 00245 libmesh_error_msg("ERROR: Unsupported element type = " << type); 00246 } 00247 } 00248 // by default throw an error 00249 default: 00250 libmesh_error_msg("ERROR: Unsupported polynomial order = " << totalorder); 00251 } 00252 00253 libmesh_error_msg("We'll never get here!"); 00254 return 0.; 00255 } 00256 00257 00258 00259 template <> 00260 Real FE<1,HERMITE>::shape_deriv(const ElemType, 00261 const Order, 00262 const unsigned int, 00263 const unsigned int, 00264 const Point&) 00265 { 00266 libmesh_error_msg("Hermite elements require the real element \nto construct gradient-based degrees of freedom."); 00267 return 0.; 00268 } 00269 00270 00271 00272 template <> 00273 Real FE<1,HERMITE>::shape_deriv(const Elem* elem, 00274 const Order order, 00275 const unsigned int i, 00276 const unsigned int, 00277 const Point& p) 00278 { 00279 libmesh_assert(elem); 00280 00281 // Coefficient naming: d(1)d(2n) is the coefficient of the 00282 // global shape function corresponding to value 1 in terms of the 00283 // local shape function corresponding to normal derivative 2 00284 Real d1xd1x, d2xd2x; 00285 00286 hermite_compute_coefs(elem, d1xd1x, d2xd2x); 00287 00288 const ElemType type = elem->type(); 00289 00290 const Order totalorder = static_cast<Order>(order + elem->p_level()); 00291 00292 switch (totalorder) 00293 { 00294 // Hermite cubic shape functions 00295 case THIRD: 00296 { 00297 switch (type) 00298 { 00299 // C1 functions on the C1 cubic edge 00300 case EDGE2: 00301 case EDGE3: 00302 { 00303 switch (i) 00304 { 00305 case 0: 00306 return FEHermite<1>::hermite_raw_shape_deriv(0, p(0)); 00307 case 1: 00308 return d1xd1x * FEHermite<1>::hermite_raw_shape_deriv(2, p(0)); 00309 case 2: 00310 return FEHermite<1>::hermite_raw_shape_deriv(1, p(0)); 00311 case 3: 00312 return d2xd2x * FEHermite<1>::hermite_raw_shape_deriv(3, p(0)); 00313 default: 00314 return FEHermite<1>::hermite_raw_shape_deriv(i, p(0)); 00315 } 00316 } 00317 default: 00318 libmesh_error_msg("ERROR: Unsupported element type = " << type); 00319 } 00320 } 00321 // by default throw an error 00322 default: 00323 libmesh_error_msg("ERROR: Unsupported polynomial order = " << totalorder); 00324 } 00325 00326 libmesh_error_msg("We'll never get here!"); 00327 return 0.; 00328 } 00329 00330 00331 00332 template <> 00333 Real FE<1,HERMITE>::shape_second_deriv(const Elem* elem, 00334 const Order order, 00335 const unsigned int i, 00336 const unsigned int, 00337 const Point& p) 00338 { 00339 libmesh_assert(elem); 00340 00341 // Coefficient naming: d(1)d(2n) is the coefficient of the 00342 // global shape function corresponding to value 1 in terms of the 00343 // local shape function corresponding to normal derivative 2 00344 Real d1xd1x, d2xd2x; 00345 00346 hermite_compute_coefs(elem, d1xd1x, d2xd2x); 00347 00348 const ElemType type = elem->type(); 00349 00350 const Order totalorder = static_cast<Order>(order + elem->p_level()); 00351 00352 switch (totalorder) 00353 { 00354 // Hermite cubic shape functions 00355 case THIRD: 00356 { 00357 switch (type) 00358 { 00359 // C1 functions on the C1 cubic edge 00360 case EDGE2: 00361 case EDGE3: 00362 { 00363 switch (i) 00364 { 00365 case 0: 00366 return FEHermite<1>::hermite_raw_shape_second_deriv(0, p(0)); 00367 case 1: 00368 return d1xd1x * FEHermite<1>::hermite_raw_shape_second_deriv(2, p(0)); 00369 case 2: 00370 return FEHermite<1>::hermite_raw_shape_second_deriv(1, p(0)); 00371 case 3: 00372 return d2xd2x * FEHermite<1>::hermite_raw_shape_second_deriv(3, p(0)); 00373 default: 00374 return FEHermite<1>::hermite_raw_shape_second_deriv(i, p(0)); 00375 } 00376 } 00377 default: 00378 libmesh_error_msg("ERROR: Unsupported element type = " << type); 00379 } 00380 } 00381 // by default throw an error 00382 default: 00383 libmesh_error_msg("ERROR: Unsupported polynomial order = " << totalorder); 00384 } 00385 00386 libmesh_error_msg("We'll never get here!"); 00387 return 0.; 00388 } 00389 00390 } // namespace libMesh