fe_nedelec_one_shape_3D.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 namespace libMesh 00026 { 00027 00028 template <> 00029 RealGradient FE<3,NEDELEC_ONE>::shape(const ElemType, 00030 const Order, 00031 const unsigned int, 00032 const Point&) 00033 { 00034 #if LIBMESH_DIM == 3 00035 libMesh::err << "Nedelec elements require the element type\n" 00036 << "because edge orientation is needed." 00037 << std::endl; 00038 libmesh_error(); 00039 #endif 00040 00041 libmesh_error(); 00042 return RealGradient(); 00043 } 00044 00045 00046 00047 template <> 00048 RealGradient FE<3,NEDELEC_ONE>::shape(const Elem* elem, 00049 const Order order, 00050 const unsigned int i, 00051 const Point& /*p*/) 00052 { 00053 #if LIBMESH_DIM == 3 00054 libmesh_assert(elem); 00055 00056 const Order totalorder = static_cast<Order>(order + elem->p_level()); 00057 00058 switch (totalorder) 00059 { 00060 // linear Lagrange shape functions 00061 case FIRST: 00062 { 00063 switch (elem->type()) 00064 { 00065 case HEX20: 00066 case HEX27: 00067 { 00068 libmesh_assert_less (i, 12); 00069 00070 libmesh_not_implemented(); 00071 return RealGradient(); 00072 } 00073 00074 case TET10: 00075 { 00076 libmesh_assert_less (i, 6); 00077 00078 libmesh_not_implemented(); 00079 return RealGradient(); 00080 } 00081 00082 default: 00083 { 00084 libMesh::err << "ERROR: Unsupported 3D element type!: " << elem->type() 00085 << std::endl; 00086 libmesh_error(); 00087 } 00088 } 00089 } 00090 00091 // unsupported order 00092 default: 00093 { 00094 libMesh::err << "ERROR: Unsupported 3D FE order!: " << totalorder 00095 << std::endl; 00096 00097 libmesh_error(); 00098 } 00099 } 00100 #endif 00101 00102 libmesh_error(); 00103 return RealGradient(); 00104 } 00105 00106 00107 00108 00109 template <> 00110 RealGradient FE<3,NEDELEC_ONE>::shape_deriv(const ElemType, 00111 const Order, 00112 const unsigned int, 00113 const unsigned int, 00114 const Point&) 00115 { 00116 #if LIBMESH_DIM == 3 00117 libMesh::err << "Nedelec elements require the element type\n" 00118 << "because edge orientation is needed." 00119 << std::endl; 00120 libmesh_error(); 00121 #endif 00122 00123 libmesh_error(); 00124 return RealGradient(); 00125 } 00126 00127 template <> 00128 RealGradient FE<3,NEDELEC_ONE>::shape_deriv(const Elem* elem, 00129 const Order order, 00130 const unsigned int i, 00131 const unsigned int j, 00132 const Point& /*p*/) 00133 { 00134 #if LIBMESH_DIM == 3 00135 libmesh_assert(elem); 00136 libmesh_assert_less (j, 3); 00137 00138 const Order totalorder = static_cast<Order>(order + elem->p_level()); 00139 00140 switch (totalorder) 00141 { 00142 case FIRST: 00143 { 00144 switch (elem->type()) 00145 { 00146 case HEX20: 00147 case HEX27: 00148 { 00149 libmesh_assert_less (i, 12); 00150 00151 libmesh_not_implemented(); 00152 return RealGradient(); 00153 } 00154 00155 case TET10: 00156 { 00157 libmesh_assert_less (i, 6); 00158 00159 libmesh_not_implemented(); 00160 return RealGradient(); 00161 } 00162 00163 default: 00164 { 00165 libMesh::err << "ERROR: Unsupported 3D element type!: " << elem->type() 00166 << std::endl; 00167 libmesh_error(); 00168 } 00169 } 00170 } 00171 00172 // unsupported order 00173 default: 00174 { 00175 libMesh::err << "ERROR: Unsupported 3D FE order!: " << totalorder 00176 << std::endl; 00177 libmesh_error(); 00178 } 00179 } 00180 00181 #endif 00182 00183 libmesh_error(); 00184 return RealGradient(); 00185 } 00186 00187 00188 00189 template <> 00190 RealGradient FE<3,NEDELEC_ONE>::shape_second_deriv(const ElemType, 00191 const Order, 00192 const unsigned int, 00193 const unsigned int, 00194 const Point&) 00195 { 00196 #if LIBMESH_DIM == 3 00197 libMesh::err << "Nedelec elements require the element type\n" 00198 << "because edge orientation is needed." 00199 << std::endl; 00200 libmesh_error(); 00201 #endif 00202 00203 libmesh_error(); 00204 return RealGradient(); 00205 } 00206 00207 00208 00209 template <> 00210 RealGradient FE<3,NEDELEC_ONE>::shape_second_deriv(const Elem* elem, 00211 const Order order, 00212 const unsigned int i, 00213 const unsigned int j, 00214 const Point& /*p*/) 00215 { 00216 #if LIBMESH_DIM == 3 00217 00218 libmesh_assert(elem); 00219 00220 libmesh_assert_less (j, 6); 00221 00222 const Order totalorder = static_cast<Order>(order + elem->p_level()); 00223 00224 switch (totalorder) 00225 { 00226 // linear Lagrange shape functions 00227 case FIRST: 00228 { 00229 switch (elem->type()) 00230 { 00231 case HEX20: 00232 case HEX27: 00233 { 00234 libmesh_assert_less (i, 12); 00235 00236 libmesh_not_implemented(); 00237 return RealGradient(); 00238 } 00239 00240 case TET10: 00241 { 00242 libmesh_assert_less (i, 6); 00243 00244 libmesh_not_implemented(); 00245 return RealGradient(); 00246 } 00247 00248 default: 00249 { 00250 libMesh::err << "ERROR: Unsupported 3D element type!: " << elem->type() 00251 << std::endl; 00252 libmesh_error(); 00253 } 00254 00255 } //switch(type) 00256 00257 } // case FIRST: 00258 // unsupported order 00259 default: 00260 { 00261 libMesh::err << "ERROR: Unsupported 3D FE order!: " << totalorder 00262 << std::endl; 00263 libmesh_error(); 00264 } 00265 } 00266 00267 #endif 00268 00269 libmesh_error(); 00270 return RealGradient(); 00271 } 00272 00273 } // namespace libMesh
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:46 UTC
Hosted By: