fe_xyz_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 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 libMesh::Point centroid; 00037 static Real max_distance; 00038 } 00039 00040 00041 namespace libMesh 00042 { 00043 00044 00045 template <> 00046 Real FE<1,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<1,XYZ>::shape(const Elem* elem, 00063 const Order libmesh_dbg_var(order), 00064 const unsigned int i, 00065 const Point& point_in) 00066 { 00067 libmesh_assert(elem); 00068 libmesh_assert_less_equal (i, order + elem->p_level()); 00069 00070 // Only recompute the centroid if the element 00071 // has changed from the last one we computed. 00072 // This avoids repeated centroid calculations 00073 // when called in succession with the same element. 00074 if (elem->id() != old_elem_id) 00075 { 00076 centroid = elem->centroid(); 00077 old_elem_id = elem->id(); 00078 max_distance = 0.; 00079 for (unsigned int p = 0; p < elem->n_nodes(); p++) 00080 { 00081 const Real distance = std::abs(centroid(0) - elem->point(p)(0)); 00082 max_distance = std::max(distance, max_distance); 00083 } 00084 } 00085 00086 // Using static globals for old_elem_id, etc. will fail 00087 // horribly with more than one thread. 00088 libmesh_assert_equal_to (libMesh::n_threads(), 1); 00089 00090 const Real x = point_in(0); 00091 const Real xc = centroid(0); 00092 const Real dx = (x - xc)/max_distance; 00093 00094 // monomials. since they are hierarchic we only need one case block. 00095 switch (i) 00096 { 00097 case 0: 00098 return 1.; 00099 00100 case 1: 00101 return dx; 00102 00103 case 2: 00104 return dx*dx; 00105 00106 case 3: 00107 return dx*dx*dx; 00108 00109 case 4: 00110 return dx*dx*dx*dx; 00111 00112 default: 00113 Real val = 1.; 00114 for (unsigned int index = 0; index != i; ++index) 00115 val *= dx; 00116 return val; 00117 } 00118 00119 libmesh_error(); 00120 return 0.; 00121 00122 } 00123 00124 00125 00126 template <> 00127 Real FE<1,XYZ>::shape_deriv(const ElemType, 00128 const Order, 00129 const unsigned int, 00130 const unsigned int, 00131 const Point&) 00132 { 00133 libMesh::err << "XYZ polynomials require the element\n" 00134 << "because the centroid is needed." 00135 << std::endl; 00136 00137 libmesh_error(); 00138 return 0.; 00139 } 00140 00141 00142 00143 template <> 00144 Real FE<1,XYZ>::shape_deriv(const Elem* elem, 00145 const Order libmesh_dbg_var(order), 00146 const unsigned int i, 00147 const unsigned int libmesh_dbg_var(j), 00148 const Point& point_in) 00149 { 00150 libmesh_assert(elem); 00151 libmesh_assert_less_equal (i, order + elem->p_level()); 00152 00153 // only d()/dxi in 1D! 00154 00155 libmesh_assert_equal_to (j, 0); 00156 00157 // Only recompute the centroid if the element 00158 // has changed from the last one we computed. 00159 // This avoids repeated centroid calculations 00160 // when called in succession with the same element. 00161 if (elem->id() != old_elem_id) 00162 { 00163 centroid = elem->centroid(); 00164 old_elem_id = elem->id(); 00165 max_distance = 0.; 00166 for (unsigned int p = 0; p < elem->n_nodes(); p++) 00167 { 00168 const Real distance = std::abs(centroid(0) - elem->point(p)(0)); 00169 max_distance = std::max(distance, max_distance); 00170 } 00171 } 00172 00173 // Using static globals for old_elem_id, etc. will fail 00174 // horribly with more than one thread. 00175 libmesh_assert_equal_to (libMesh::n_threads(), 1); 00176 00177 const Real x = point_in(0); 00178 const Real xc = centroid(0); 00179 const Real dx = (x - xc)/max_distance; 00180 00181 // monomials. since they are hierarchic we only need one case block. 00182 switch (i) 00183 { 00184 case 0: 00185 return 0.; 00186 00187 case 1: 00188 return 1.; 00189 00190 case 2: 00191 return 2.*dx/max_distance; 00192 00193 case 3: 00194 return 3.*dx*dx/max_distance; 00195 00196 case 4: 00197 return 4.*dx*dx*dx/max_distance; 00198 00199 default: 00200 Real val = i; 00201 for (unsigned int index = 1; index != i; ++index) 00202 val *= dx; 00203 return val/max_distance; 00204 } 00205 00206 libmesh_error(); 00207 return 0.; 00208 } 00209 00210 00211 00212 template <> 00213 Real FE<1,XYZ>::shape_second_deriv(const ElemType, 00214 const Order, 00215 const unsigned int, 00216 const unsigned int, 00217 const Point&) 00218 { 00219 libMesh::err << "XYZ polynomials require the element\n" 00220 << "because the centroid is needed." 00221 << std::endl; 00222 00223 libmesh_error(); 00224 return 0.; 00225 } 00226 00227 00228 00229 template <> 00230 Real FE<1,XYZ>::shape_second_deriv(const Elem* elem, 00231 const Order libmesh_dbg_var(order), 00232 const unsigned int i, 00233 const unsigned int libmesh_dbg_var(j), 00234 const Point& point_in) 00235 { 00236 libmesh_assert(elem); 00237 libmesh_assert_less_equal (i, order + elem->p_level()); 00238 00239 // only d2()/dxi2 in 1D! 00240 00241 libmesh_assert_equal_to (j, 0); 00242 00243 // Only recompute the centroid if the element 00244 // has changed from the last one we computed. 00245 // This avoids repeated centroid calculations 00246 // when called in succession with the same element. 00247 if (elem->id() != old_elem_id) 00248 { 00249 centroid = elem->centroid(); 00250 old_elem_id = elem->id(); 00251 max_distance = 0.; 00252 for (unsigned int p = 0; p < elem->n_nodes(); p++) 00253 { 00254 const Real distance = std::abs(centroid(0) - elem->point(p)(0)); 00255 max_distance = std::max(distance, max_distance); 00256 } 00257 } 00258 00259 // Using static globals for old_elem_id, etc. will fail 00260 // horribly with more than one thread. 00261 libmesh_assert_equal_to (libMesh::n_threads(), 1); 00262 00263 const Real x = point_in(0); 00264 const Real xc = centroid(0); 00265 const Real dx = (x - xc)/max_distance; 00266 const Real dist2 = pow(max_distance,2.); 00267 00268 // monomials. since they are hierarchic we only need one case block. 00269 switch (i) 00270 { 00271 case 0: 00272 case 1: 00273 return 0.; 00274 00275 case 2: 00276 return 2./dist2; 00277 00278 case 3: 00279 return 6.*dx/dist2; 00280 00281 case 4: 00282 return 12.*dx*dx/dist2; 00283 00284 default: 00285 Real val = 2.; 00286 for (unsigned int index = 2; index != i; ++index) 00287 val *= (index+1) * dx; 00288 return val/dist2; 00289 } 00290 00291 libmesh_error(); 00292 return 0.; 00293 } 00294 00295 } // namespace libMesh
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:47 UTC
Hosted By: