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