$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 namespace libMesh 00026 { 00027 00028 00029 00030 00031 template <> 00032 Real FE<1,MONOMIAL>::shape(const ElemType, 00033 const Order libmesh_dbg_var(order), 00034 const unsigned int i, 00035 const Point& p) 00036 { 00037 const Real xi = p(0); 00038 00039 libmesh_assert_less_equal (i, static_cast<unsigned int>(order)); 00040 00041 // monomials. since they are hierarchic we only need one case block. 00042 switch (i) 00043 { 00044 case 0: 00045 return 1.; 00046 00047 case 1: 00048 return xi; 00049 00050 case 2: 00051 return xi*xi; 00052 00053 case 3: 00054 return xi*xi*xi; 00055 00056 case 4: 00057 return xi*xi*xi*xi; 00058 00059 default: 00060 Real val = 1.; 00061 for (unsigned int index = 0; index != i; ++index) 00062 val *= xi; 00063 return val; 00064 } 00065 00066 libmesh_error_msg("We'll never get here!"); 00067 return 0.; 00068 } 00069 00070 00071 00072 template <> 00073 Real FE<1,MONOMIAL>::shape(const Elem* elem, 00074 const Order order, 00075 const unsigned int i, 00076 const Point& p) 00077 { 00078 libmesh_assert(elem); 00079 00080 return FE<1,MONOMIAL>::shape(elem->type(), static_cast<Order>(order + elem->p_level()), i, p); 00081 } 00082 00083 00084 00085 template <> 00086 Real FE<1,MONOMIAL>::shape_deriv(const ElemType, 00087 const Order libmesh_dbg_var(order), 00088 const unsigned int i, 00089 const unsigned int libmesh_dbg_var(j), 00090 const Point& p) 00091 { 00092 // only d()/dxi in 1D! 00093 00094 libmesh_assert_equal_to (j, 0); 00095 00096 const Real xi = p(0); 00097 00098 libmesh_assert_less_equal (i, static_cast<unsigned int>(order)); 00099 00100 // monomials. since they are hierarchic we only need one case block. 00101 switch (i) 00102 { 00103 case 0: 00104 return 0.; 00105 00106 case 1: 00107 return 1.; 00108 00109 case 2: 00110 return 2.*xi; 00111 00112 case 3: 00113 return 3.*xi*xi; 00114 00115 case 4: 00116 return 4.*xi*xi*xi; 00117 00118 default: 00119 Real val = i; 00120 for (unsigned int index = 1; index != i; ++index) 00121 val *= xi; 00122 return val; 00123 } 00124 00125 libmesh_error_msg("We'll never get here!"); 00126 return 0.; 00127 } 00128 00129 00130 00131 template <> 00132 Real FE<1,MONOMIAL>::shape_deriv(const Elem* elem, 00133 const Order order, 00134 const unsigned int i, 00135 const unsigned int j, 00136 const Point& p) 00137 { 00138 libmesh_assert(elem); 00139 00140 return FE<1,MONOMIAL>::shape_deriv(elem->type(), 00141 static_cast<Order>(order + elem->p_level()), i, j, p); 00142 } 00143 00144 00145 00146 template <> 00147 Real FE<1,MONOMIAL>::shape_second_deriv(const ElemType, 00148 const Order libmesh_dbg_var(order), 00149 const unsigned int i, 00150 const unsigned int libmesh_dbg_var(j), 00151 const Point& p) 00152 { 00153 // only d()/dxi in 1D! 00154 00155 libmesh_assert_equal_to (j, 0); 00156 00157 const Real xi = p(0); 00158 00159 libmesh_assert_less_equal (i, static_cast<unsigned int>(order)); 00160 00161 switch (i) 00162 { 00163 case 0: 00164 case 1: 00165 return 0.; 00166 00167 case 2: 00168 return 2.; 00169 00170 case 3: 00171 return 6.*xi; 00172 00173 case 4: 00174 return 12.*xi*xi; 00175 00176 default: 00177 Real val = 2.; 00178 for (unsigned int index = 2; index != i; ++index) 00179 val *= (index+1) * xi; 00180 return val; 00181 } 00182 00183 libmesh_error_msg("We'll never get here!"); 00184 return 0.; 00185 } 00186 00187 00188 00189 template <> 00190 Real FE<1,MONOMIAL>::shape_second_deriv(const Elem* elem, 00191 const Order order, 00192 const unsigned int i, 00193 const unsigned int j, 00194 const Point& p) 00195 { 00196 libmesh_assert(elem); 00197 00198 return FE<1,MONOMIAL>::shape_second_deriv(elem->type(), 00199 static_cast<Order>(order + elem->p_level()), i, j, p); 00200 } 00201 00202 } // namespace libMesh