fe_clough.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/elem.h" 00022 #include "libmesh/fe.h" 00023 #include "libmesh/fe_interface.h" 00024 00025 00026 namespace libMesh 00027 { 00028 00029 // ------------------------------------------------------------ 00030 // Clough-specific implementations 00031 00032 // Anonymous namespace for local helper functions 00033 namespace { 00034 00035 void clough_nodal_soln(const Elem* elem, 00036 const Order order, 00037 const std::vector<Number>& elem_soln, 00038 std::vector<Number>& nodal_soln, 00039 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 // FEType object to be passed to various FEInterface functions below. 00050 FEType fe_type(totalorder, CLOUGH); 00051 00052 switch (totalorder) 00053 { 00054 // Piecewise cubic shape functions with linear flux edges 00055 case SECOND: 00056 // Piecewise cubic shape functions 00057 case THIRD: 00058 { 00059 00060 const unsigned int n_sf = 00061 // FE<Dim,T>::n_shape_functions(elem_type, totalorder); 00062 FEInterface::n_shape_functions(Dim, fe_type, elem_type); 00063 00064 std::vector<Point> refspace_nodes; 00065 FEBase::get_refspace_nodes(elem_type,refspace_nodes); 00066 libmesh_assert_equal_to (refspace_nodes.size(), n_nodes); 00067 00068 for (unsigned int n=0; n<n_nodes; n++) 00069 { 00070 libmesh_assert_equal_to (elem_soln.size(), n_sf); 00071 00072 // Zero before summation 00073 nodal_soln[n] = 0; 00074 00075 // u_i = Sum (alpha_i phi_i) 00076 for (unsigned int i=0; i<n_sf; i++) 00077 nodal_soln[n] += elem_soln[i] * 00078 // FE<Dim,T>::shape(elem, order, i, mapped_point); 00079 FEInterface::shape(Dim, fe_type, elem, i, refspace_nodes[n]); 00080 } 00081 00082 return; 00083 } 00084 00085 default: 00086 { 00087 libmesh_error(); 00088 } 00089 } 00090 } // clough_nodal_soln() 00091 00092 00093 00094 00095 unsigned int clough_n_dofs(const ElemType t, const Order o) 00096 { 00097 switch (o) 00098 { 00099 // Piecewise cubic shape functions with linear flux edges 00100 case SECOND: 00101 { 00102 switch (t) 00103 { 00104 case TRI6: 00105 return 9; 00106 00107 default: 00108 { 00109 #ifdef DEBUG 00110 libMesh::err << "ERROR: Bad ElemType = " << t 00111 << " for " << o << "th order approximation!" 00112 << std::endl; 00113 #endif 00114 libmesh_error(); 00115 } 00116 } 00117 } 00118 // Piecewise cubic Clough-Tocher element 00119 case THIRD: 00120 { 00121 switch (t) 00122 { 00123 case TRI6: 00124 return 12; 00125 00126 default: 00127 { 00128 #ifdef DEBUG 00129 libMesh::err << "ERROR: Bad ElemType = " << t 00130 << " for " << o << "th order approximation!" 00131 << std::endl; 00132 #endif 00133 libmesh_error(); 00134 } 00135 } 00136 } 00137 00138 default: 00139 { 00140 libmesh_error(); 00141 } 00142 } 00143 00144 libmesh_error(); 00145 return 0; 00146 } // clough_n_dofs() 00147 00148 00149 00150 00151 00152 unsigned int clough_n_dofs_at_node(const ElemType t, 00153 const Order o, 00154 const unsigned int n) 00155 { 00156 switch (o) 00157 { 00158 // Piecewise cubic shape functions with linear flux edges 00159 case SECOND: 00160 { 00161 switch (t) 00162 { 00163 // The 2D Clough-Tocher defined on a 6-noded triangle 00164 case TRI6: 00165 { 00166 switch (n) 00167 { 00168 case 0: 00169 case 1: 00170 case 2: 00171 return 3; 00172 00173 case 3: 00174 case 4: 00175 case 5: 00176 return 0; 00177 00178 default: 00179 libmesh_error(); 00180 } 00181 } 00182 00183 default: 00184 { 00185 #ifdef DEBUG 00186 libMesh::err << "ERROR: Bad ElemType = " << t 00187 << " for " << o << "th order approximation!" 00188 << std::endl; 00189 #endif 00190 libmesh_error(); 00191 } 00192 00193 } 00194 } 00195 // The third-order clough shape functions 00196 case THIRD: 00197 { 00198 switch (t) 00199 { 00200 // The 2D Clough-Tocher defined on a 6-noded triangle 00201 case TRI6: 00202 { 00203 switch (n) 00204 { 00205 case 0: 00206 case 1: 00207 case 2: 00208 return 3; 00209 00210 case 3: 00211 case 4: 00212 case 5: 00213 return 1; 00214 00215 default: 00216 libmesh_error(); 00217 } 00218 } 00219 00220 default: 00221 { 00222 #ifdef DEBUG 00223 libMesh::err << "ERROR: Bad ElemType = " << t 00224 << " for " << o << "th order approximation!" 00225 << std::endl; 00226 #endif 00227 libmesh_error(); 00228 } 00229 00230 } 00231 } 00232 default: 00233 { 00234 libmesh_error(); 00235 } 00236 } 00237 00238 libmesh_error(); 00239 00240 return 0; 00241 } // clough_n_dofs_at_node() 00242 00243 00244 00245 00246 unsigned int clough_n_dofs_per_elem(const ElemType t, const Order o) 00247 { 00248 switch (o) 00249 { 00250 // Piecewise cubic shape functions with linear flux edges 00251 case SECOND: 00252 // The third-order Clough-Tocher shape functions 00253 case THIRD: 00254 { 00255 switch (t) 00256 { 00257 // The 2D clough defined on a 6-noded triangle 00258 case TRI6: 00259 return 0; 00260 00261 default: 00262 { 00263 #ifdef DEBUG 00264 libMesh::err << "ERROR: Bad ElemType = " << t 00265 << " for " << o << "th order approximation!" 00266 << std::endl; 00267 #endif 00268 libmesh_error(); 00269 } 00270 00271 } 00272 } 00273 // Otherwise no DOFS per element 00274 default: 00275 libmesh_error(); 00276 return 0; 00277 } 00278 } // clough_n_dofs_per_elem() 00279 00280 } // anonymous 00281 00282 00283 00284 00285 // Do full-specialization of nodal_soln() function for every 00286 // dimension, instead of explicit instantiation at the end of this 00287 // file. 00288 // This could be macro-ified so that it fits on one line... 00289 template <> 00290 void FE<0,CLOUGH>::nodal_soln(const Elem* elem, 00291 const Order order, 00292 const std::vector<Number>& elem_soln, 00293 std::vector<Number>& nodal_soln) 00294 { clough_nodal_soln(elem, order, elem_soln, nodal_soln, /*Dim=*/0); } 00295 00296 template <> 00297 void FE<1,CLOUGH>::nodal_soln(const Elem* elem, 00298 const Order order, 00299 const std::vector<Number>& elem_soln, 00300 std::vector<Number>& nodal_soln) 00301 { clough_nodal_soln(elem, order, elem_soln, nodal_soln, /*Dim=*/1); } 00302 00303 template <> 00304 void FE<2,CLOUGH>::nodal_soln(const Elem* elem, 00305 const Order order, 00306 const std::vector<Number>& elem_soln, 00307 std::vector<Number>& nodal_soln) 00308 { clough_nodal_soln(elem, order, elem_soln, nodal_soln, /*Dim=*/2); } 00309 00310 template <> 00311 void FE<3,CLOUGH>::nodal_soln(const Elem* elem, 00312 const Order order, 00313 const std::vector<Number>& elem_soln, 00314 std::vector<Number>& nodal_soln) 00315 { clough_nodal_soln(elem, order, elem_soln, nodal_soln, /*Dim=*/3); } 00316 00317 00318 // Full specialization of n_dofs() function for every dimension 00319 template <> unsigned int FE<0,CLOUGH>::n_dofs(const ElemType t, const Order o) { return clough_n_dofs(t, o); } 00320 template <> unsigned int FE<1,CLOUGH>::n_dofs(const ElemType t, const Order o) { return clough_n_dofs(t, o); } 00321 template <> unsigned int FE<2,CLOUGH>::n_dofs(const ElemType t, const Order o) { return clough_n_dofs(t, o); } 00322 template <> unsigned int FE<3,CLOUGH>::n_dofs(const ElemType t, const Order o) { return clough_n_dofs(t, o); } 00323 00324 00325 // Full specialization of n_dofs_at_node() function for every dimension. 00326 template <> unsigned int FE<0,CLOUGH>::n_dofs_at_node(const ElemType t, const Order o, const unsigned int n) { return clough_n_dofs_at_node(t, o, n); } 00327 template <> unsigned int FE<1,CLOUGH>::n_dofs_at_node(const ElemType t, const Order o, const unsigned int n) { return clough_n_dofs_at_node(t, o, n); } 00328 template <> unsigned int FE<2,CLOUGH>::n_dofs_at_node(const ElemType t, const Order o, const unsigned int n) { return clough_n_dofs_at_node(t, o, n); } 00329 template <> unsigned int FE<3,CLOUGH>::n_dofs_at_node(const ElemType t, const Order o, const unsigned int n) { return clough_n_dofs_at_node(t, o, n); } 00330 00331 // Full specialization of n_dofs_per_elem() function for every dimension. 00332 template <> unsigned int FE<0,CLOUGH>::n_dofs_per_elem(const ElemType t, const Order o) { return clough_n_dofs_per_elem(t, o); } 00333 template <> unsigned int FE<1,CLOUGH>::n_dofs_per_elem(const ElemType t, const Order o) { return clough_n_dofs_per_elem(t, o); } 00334 template <> unsigned int FE<2,CLOUGH>::n_dofs_per_elem(const ElemType t, const Order o) { return clough_n_dofs_per_elem(t, o); } 00335 template <> unsigned int FE<3,CLOUGH>::n_dofs_per_elem(const ElemType t, const Order o) { return clough_n_dofs_per_elem(t, o); } 00336 00337 // Clough FEMs are C^1 continuous 00338 template <> FEContinuity FE<0,CLOUGH>::get_continuity() const { return C_ONE; } 00339 template <> FEContinuity FE<1,CLOUGH>::get_continuity() const { return C_ONE; } 00340 template <> FEContinuity FE<2,CLOUGH>::get_continuity() const { return C_ONE; } 00341 template <> FEContinuity FE<3,CLOUGH>::get_continuity() const { return C_ONE; } 00342 00343 // Clough FEMs are not (currently) hierarchic 00344 template <> bool FE<0,CLOUGH>::is_hierarchic() const { return false; } // FIXME - this will be changed 00345 template <> bool FE<1,CLOUGH>::is_hierarchic() const { return false; } // FIXME - this will be changed 00346 template <> bool FE<2,CLOUGH>::is_hierarchic() const { return false; } // FIXME - this will be changed 00347 template <> bool FE<3,CLOUGH>::is_hierarchic() const { return false; } // FIXME - this will be changed 00348 00349 #ifdef LIBMESH_ENABLE_AMR 00350 // compute_constraints() specializations are only needed for 2 and 3D 00351 template <> 00352 void FE<2,CLOUGH>::compute_constraints (DofConstraints &constraints, 00353 DofMap &dof_map, 00354 const unsigned int variable_number, 00355 const Elem* elem) 00356 { compute_proj_constraints(constraints, dof_map, variable_number, elem); } 00357 00358 template <> 00359 void FE<3,CLOUGH>::compute_constraints (DofConstraints &constraints, 00360 DofMap &dof_map, 00361 const unsigned int variable_number, 00362 const Elem* elem) 00363 { compute_proj_constraints(constraints, dof_map, variable_number, elem); } 00364 #endif // #ifdef LIBMESH_ENABLE_AMR 00365 00366 // Clough FEM shapes need reinit 00367 template <> bool FE<0,CLOUGH>::shapes_need_reinit() const { return true; } 00368 template <> bool FE<1,CLOUGH>::shapes_need_reinit() const { return true; } 00369 template <> bool FE<2,CLOUGH>::shapes_need_reinit() const { return true; } 00370 template <> bool FE<3,CLOUGH>::shapes_need_reinit() const { return true; } 00371 00372 } // namespace libMesh
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:46 UTC
Hosted By: