fe_monomial.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 00020 // Local includes 00021 #include "libmesh/fe.h" 00022 #include "libmesh/elem.h" 00023 #include "libmesh/fe_interface.h" 00024 00025 namespace libMesh 00026 { 00027 00028 // ------------------------------------------------------------ 00029 // Monomials-specific implementations 00030 00031 00032 // Anonymous namespace for local helper functions 00033 namespace { 00034 00035 void monomial_nodal_soln(const Elem* elem, 00036 const Order order, 00037 const std::vector<Number>& elem_soln, 00038 std::vector<Number>& nodal_soln, 00039 const unsigned Dim) 00040 { 00041 const unsigned int n_nodes = elem->n_nodes(); 00042 00043 const ElemType elem_type = elem->type(); 00044 00045 nodal_soln.resize(n_nodes); 00046 00047 const Order totalorder = static_cast<Order>(order+elem->p_level()); 00048 00049 switch (totalorder) 00050 { 00051 // Constant shape functions 00052 case CONSTANT: 00053 { 00054 libmesh_assert_equal_to (elem_soln.size(), 1); 00055 00056 const Number val = elem_soln[0]; 00057 00058 for (unsigned int n=0; n<n_nodes; n++) 00059 nodal_soln[n] = val; 00060 00061 return; 00062 } 00063 00064 00065 // For other orders, do interpolation at the nodes 00066 // explicitly. 00067 default: 00068 { 00069 // FEType object to be passed to various FEInterface functions below. 00070 FEType fe_type(totalorder, MONOMIAL); 00071 00072 const unsigned int n_sf = 00073 // FE<Dim,T>::n_shape_functions(elem_type, totalorder); 00074 FEInterface::n_shape_functions(Dim, fe_type, elem_type); 00075 00076 std::vector<Point> refspace_nodes; 00077 FEBase::get_refspace_nodes(elem_type,refspace_nodes); 00078 libmesh_assert_equal_to (refspace_nodes.size(), n_nodes); 00079 00080 for (unsigned int n=0; n<n_nodes; n++) 00081 { 00082 libmesh_assert_equal_to (elem_soln.size(), n_sf); 00083 00084 // Zero before summation 00085 nodal_soln[n] = 0; 00086 00087 // u_i = Sum (alpha_i phi_i) 00088 for (unsigned int i=0; i<n_sf; i++) 00089 nodal_soln[n] += elem_soln[i] * 00090 // FE<Dim,T>::shape(elem, order, i, mapped_point); 00091 FEInterface::shape(Dim, fe_type, elem, i, refspace_nodes[n]); 00092 } 00093 00094 return; 00095 } // default 00096 } // switch 00097 } // monomial_nodal_soln() 00098 00099 00100 00101 00102 unsigned int monomial_n_dofs(const ElemType t, const Order o) 00103 { 00104 switch (o) 00105 { 00106 00107 // constant shape functions 00108 // no matter what shape there is only one DOF. 00109 case CONSTANT: 00110 return 1; 00111 00112 00113 // Discontinuous linear shape functions 00114 // expressed in the monomials. 00115 case FIRST: 00116 { 00117 switch (t) 00118 { 00119 case NODEELEM: 00120 return 1; 00121 00122 case EDGE2: 00123 case EDGE3: 00124 case EDGE4: 00125 return 2; 00126 00127 case TRI3: 00128 case TRI6: 00129 case QUAD4: 00130 case QUAD8: 00131 case QUAD9: 00132 return 3; 00133 00134 case TET4: 00135 case TET10: 00136 case HEX8: 00137 case HEX20: 00138 case HEX27: 00139 case PRISM6: 00140 case PRISM15: 00141 case PRISM18: 00142 case PYRAMID5: 00143 return 4; 00144 00145 default: 00146 { 00147 #ifdef DEBUG 00148 libMesh::err << "ERROR: Bad ElemType = " << t 00149 << " for " << o << "th order approximation!" 00150 << std::endl; 00151 #endif 00152 libmesh_error(); 00153 } 00154 } 00155 } 00156 00157 00158 // Discontinuous quadratic shape functions 00159 // expressed in the monomials. 00160 case SECOND: 00161 { 00162 switch (t) 00163 { 00164 case NODEELEM: 00165 return 1; 00166 00167 case EDGE2: 00168 case EDGE3: 00169 case EDGE4: 00170 return 3; 00171 00172 case TRI3: 00173 case TRI6: 00174 case QUAD4: 00175 case QUAD8: 00176 case QUAD9: 00177 return 6; 00178 00179 case TET4: 00180 case TET10: 00181 case HEX8: 00182 case HEX20: 00183 case HEX27: 00184 case PRISM6: 00185 case PRISM15: 00186 case PRISM18: 00187 case PYRAMID5: 00188 return 10; 00189 00190 default: 00191 { 00192 #ifdef DEBUG 00193 libMesh::err << "ERROR: Bad ElemType = " << t 00194 << " for " << o << "th order approximation!" 00195 << std::endl; 00196 #endif 00197 libmesh_error(); 00198 } 00199 } 00200 } 00201 00202 00203 // Discontinuous cubic shape functions 00204 // expressed in the monomials. 00205 case THIRD: 00206 { 00207 switch (t) 00208 { 00209 case NODEELEM: 00210 return 1; 00211 00212 case EDGE2: 00213 case EDGE3: 00214 case EDGE4: 00215 return 4; 00216 00217 case TRI3: 00218 case TRI6: 00219 case QUAD4: 00220 case QUAD8: 00221 case QUAD9: 00222 return 10; 00223 00224 case TET4: 00225 case TET10: 00226 case HEX8: 00227 case HEX20: 00228 case HEX27: 00229 case PRISM6: 00230 case PRISM15: 00231 case PRISM18: 00232 case PYRAMID5: 00233 return 20; 00234 00235 default: 00236 { 00237 #ifdef DEBUG 00238 libMesh::err << "ERROR: Bad ElemType = " << t 00239 << " for " << o << "th order approximation!" 00240 << std::endl; 00241 #endif 00242 libmesh_error(); 00243 } 00244 } 00245 } 00246 00247 00248 // Discontinuous quartic shape functions 00249 // expressed in the monomials. 00250 case FOURTH: 00251 { 00252 switch (t) 00253 { 00254 case NODEELEM: 00255 return 1; 00256 00257 case EDGE2: 00258 case EDGE3: 00259 return 5; 00260 00261 case TRI3: 00262 case TRI6: 00263 case QUAD4: 00264 case QUAD8: 00265 case QUAD9: 00266 return 15; 00267 00268 case TET4: 00269 case TET10: 00270 case HEX8: 00271 case HEX20: 00272 case HEX27: 00273 case PRISM6: 00274 case PRISM15: 00275 case PRISM18: 00276 case PYRAMID5: 00277 return 35; 00278 00279 default: 00280 { 00281 #ifdef DEBUG 00282 libMesh::err << "ERROR: Bad ElemType = " << t 00283 << " for " << o << "th order approximation!" 00284 << std::endl; 00285 #endif 00286 libmesh_error(); 00287 } 00288 } 00289 } 00290 00291 00292 default: 00293 { 00294 const unsigned int order = static_cast<unsigned int>(o); 00295 switch (t) 00296 { 00297 case NODEELEM: 00298 return 1; 00299 00300 case EDGE2: 00301 case EDGE3: 00302 return (order+1); 00303 00304 case TRI3: 00305 case TRI6: 00306 case QUAD4: 00307 case QUAD8: 00308 case QUAD9: 00309 return (order+1)*(order+2)/2; 00310 00311 case TET4: 00312 case TET10: 00313 case HEX8: 00314 case HEX20: 00315 case HEX27: 00316 case PRISM6: 00317 case PRISM15: 00318 case PRISM18: 00319 case PYRAMID5: 00320 return (order+1)*(order+2)*(order+3)/6; 00321 00322 default: 00323 { 00324 #ifdef DEBUG 00325 libMesh::err << "ERROR: Bad ElemType = " << t 00326 << " for " << o << "th order approximation!" 00327 << std::endl; 00328 #endif 00329 libmesh_error(); 00330 } 00331 } 00332 } 00333 } 00334 00335 libmesh_error(); 00336 00337 return 0; 00338 } // monomial_n_dofs() 00339 00340 00341 } // anonymous namespace 00342 00343 00344 00345 00346 00347 // Do full-specialization for every dimension, instead 00348 // of explicit instantiation at the end of this file. 00349 // This could be macro-ified so that it fits on one line... 00350 template <> 00351 void FE<0,MONOMIAL>::nodal_soln(const Elem* elem, 00352 const Order order, 00353 const std::vector<Number>& elem_soln, 00354 std::vector<Number>& nodal_soln) 00355 { monomial_nodal_soln(elem, order, elem_soln, nodal_soln, /*Dim=*/0); } 00356 00357 template <> 00358 void FE<1,MONOMIAL>::nodal_soln(const Elem* elem, 00359 const Order order, 00360 const std::vector<Number>& elem_soln, 00361 std::vector<Number>& nodal_soln) 00362 { monomial_nodal_soln(elem, order, elem_soln, nodal_soln, /*Dim=*/1); } 00363 00364 template <> 00365 void FE<2,MONOMIAL>::nodal_soln(const Elem* elem, 00366 const Order order, 00367 const std::vector<Number>& elem_soln, 00368 std::vector<Number>& nodal_soln) 00369 { monomial_nodal_soln(elem, order, elem_soln, nodal_soln, /*Dim=*/2); } 00370 00371 template <> 00372 void FE<3,MONOMIAL>::nodal_soln(const Elem* elem, 00373 const Order order, 00374 const std::vector<Number>& elem_soln, 00375 std::vector<Number>& nodal_soln) 00376 { monomial_nodal_soln(elem, order, elem_soln, nodal_soln, /*Dim=*/3); } 00377 00378 00379 // Full specialization of n_dofs() function for every dimension 00380 template <> unsigned int FE<0,MONOMIAL>::n_dofs(const ElemType t, const Order o) { return monomial_n_dofs(t, o); } 00381 template <> unsigned int FE<1,MONOMIAL>::n_dofs(const ElemType t, const Order o) { return monomial_n_dofs(t, o); } 00382 template <> unsigned int FE<2,MONOMIAL>::n_dofs(const ElemType t, const Order o) { return monomial_n_dofs(t, o); } 00383 template <> unsigned int FE<3,MONOMIAL>::n_dofs(const ElemType t, const Order o) { return monomial_n_dofs(t, o); } 00384 00385 // Full specialization of n_dofs_at_node() function for every dimension. 00386 // Monomials have no dofs at nodes, only element dofs. 00387 template <> unsigned int FE<0,MONOMIAL>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { return 0; } 00388 template <> unsigned int FE<1,MONOMIAL>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { return 0; } 00389 template <> unsigned int FE<2,MONOMIAL>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { return 0; } 00390 template <> unsigned int FE<3,MONOMIAL>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { return 0; } 00391 00392 // Full specialization of n_dofs_per_elem() function for every dimension. 00393 template <> unsigned int FE<0,MONOMIAL>::n_dofs_per_elem(const ElemType t, const Order o) { return monomial_n_dofs(t, o); } 00394 template <> unsigned int FE<1,MONOMIAL>::n_dofs_per_elem(const ElemType t, const Order o) { return monomial_n_dofs(t, o); } 00395 template <> unsigned int FE<2,MONOMIAL>::n_dofs_per_elem(const ElemType t, const Order o) { return monomial_n_dofs(t, o); } 00396 template <> unsigned int FE<3,MONOMIAL>::n_dofs_per_elem(const ElemType t, const Order o) { return monomial_n_dofs(t, o); } 00397 00398 00399 // Full specialization of get_continuity() function for every dimension. 00400 template <> FEContinuity FE<0,MONOMIAL>::get_continuity() const { return DISCONTINUOUS; } 00401 template <> FEContinuity FE<1,MONOMIAL>::get_continuity() const { return DISCONTINUOUS; } 00402 template <> FEContinuity FE<2,MONOMIAL>::get_continuity() const { return DISCONTINUOUS; } 00403 template <> FEContinuity FE<3,MONOMIAL>::get_continuity() const { return DISCONTINUOUS; } 00404 00405 // Full specialization of is_hierarchic() function for every dimension. 00406 // The monomials are hierarchic! 00407 template <> bool FE<0,MONOMIAL>::is_hierarchic() const { return true; } 00408 template <> bool FE<1,MONOMIAL>::is_hierarchic() const { return true; } 00409 template <> bool FE<2,MONOMIAL>::is_hierarchic() const { return true; } 00410 template <> bool FE<3,MONOMIAL>::is_hierarchic() const { return true; } 00411 00412 #ifdef LIBMESH_ENABLE_AMR 00413 00414 // Full specialization of compute_constraints() function for 2D and 00415 // 3D only. There are no constraints for discontinuous elements, so 00416 // we do nothing. 00417 template <> void FE<2,MONOMIAL>::compute_constraints (DofConstraints &, DofMap &, const unsigned int, const Elem*) {} 00418 template <> void FE<3,MONOMIAL>::compute_constraints (DofConstraints &, DofMap &, const unsigned int, const Elem*) {} 00419 00420 #endif // #ifdef LIBMESH_ENABLE_AMR 00421 00422 // Full specialization of shapes_need_reinit() function for every dimension. 00423 template <> bool FE<0,MONOMIAL>::shapes_need_reinit() const { return false; } 00424 template <> bool FE<1,MONOMIAL>::shapes_need_reinit() const { return false; } 00425 template <> bool FE<2,MONOMIAL>::shapes_need_reinit() const { return false; } 00426 template <> bool FE<3,MONOMIAL>::shapes_need_reinit() const { return false; } 00427 00428 } // namespace libMesh
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:46 UTC
Hosted By: