$extrastylesheet
fe_szabab_shape_1D.C
Go to the documentation of this file.
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 
00020 // C++ includes
00021 
00022 // Local includes
00023 #include "libmesh/libmesh_config.h"
00024 
00025 #ifdef LIBMESH_ENABLE_HIGHER_ORDER_SHAPES
00026 
00027 #include "libmesh/fe.h"
00028 #include "libmesh/elem.h"
00029 
00030 namespace libMesh
00031 {
00032 
00033 
00034 template <>
00035 Real FE<1,SZABAB>::shape(const ElemType,
00036                          const Order libmesh_dbg_var(order),
00037                          const unsigned int i,
00038                          const Point& p)
00039 {
00040   const Real xi  = p(0);
00041   const Real xi2 = xi*xi;
00042 
00043 
00044   // Use this libmesh_assert rather than a switch with a single entry...
00045   // It will go away in optimized mode, essentially has the same effect.
00046   libmesh_assert_less_equal (order, SEVENTH);
00047 
00048   //   switch (order)
00049   //     {
00050   //     case FIRST:
00051   //     case SECOND:
00052   //     case THIRD:
00053   //     case FOURTH:
00054   //     case FIFTH:
00055   //     case SIXTH:
00056   //     case SEVENTH:
00057 
00058   switch(i)
00059     {
00060       //nodal shape functions
00061     case 0: return 1./2.-1./2.*xi;
00062     case 1: return 1./2.+1./2.*xi;
00063     case 2: return 1./4.  *2.4494897427831780982*(xi2-1.);
00064     case 3: return 1./4.  *3.1622776601683793320*(xi2-1.)*xi;
00065     case 4: return 1./16. *3.7416573867739413856*((5.*xi2-6.)*xi2+1.);
00066     case 5: return 3./16. *1.4142135623730950488*(3.+(-10.+7.*xi2)*xi2)*xi;
00067     case 6: return 1./32. *4.6904157598234295546*(-1.+(15.+(-35.+21.*xi2)*xi2)*xi2);
00068     case 7: return 1./32. *5.0990195135927848300*(-5.+(35.+(-63.+33.*xi2)*xi2)*xi2)*xi;
00069     case 8: return 1./256.*5.4772255750516611346*(5.+(-140.+(630.+(-924.+429.*xi2)*xi2)*xi2)*xi2);
00070 
00071     default:
00072       libmesh_error_msg("Invalid shape function index!");
00073     }
00074 
00075   libmesh_error_msg("We'll never get here!");
00076   return 0.;
00077 }
00078 
00079 
00080 
00081 template <>
00082 Real FE<1,SZABAB>::shape(const Elem* elem,
00083                          const Order order,
00084                          const unsigned int i,
00085                          const Point& p)
00086 {
00087   libmesh_assert(elem);
00088 
00089   return FE<1,SZABAB>::shape(elem->type(), static_cast<Order>(order + elem->p_level()), i, p);
00090 }
00091 
00092 
00093 
00094 template <>
00095 Real FE<1,SZABAB>::shape_deriv(const ElemType,
00096                                const Order libmesh_dbg_var(order),
00097                                const unsigned int i,
00098                                const unsigned int libmesh_dbg_var(j),
00099                                const Point& p)
00100 {
00101   // only d()/dxi in 1D!
00102   libmesh_assert_equal_to (j, 0);
00103 
00104   const Real xi  = p(0);
00105   const Real xi2 = xi*xi;
00106 
00107   // Use this libmesh_assert rather than a switch with a single entry...
00108   // It will go away in optimized mode, essentially has the same effect.
00109   libmesh_assert_less_equal (order, SEVENTH);
00110 
00111   //   switch (order)
00112   //     {
00113   //     case FIRST:
00114   //     case SECOND:
00115   //     case THIRD:
00116   //     case FOURTH:
00117   //     case FIFTH:
00118   //     case SIXTH:
00119   //     case SEVENTH:
00120 
00121   switch(i)
00122     {
00123     case 0:return -1./2.;
00124     case 1:return 1./2.;
00125     case 2:return 1./2.*2.4494897427831780982*xi;
00126     case 3:return -1./4.*3.1622776601683793320+3./4.*3.1622776601683793320*xi2;
00127     case 4:return 1./16.*3.7416573867739413856*(-12.+20*xi2)*xi;
00128     case 5:return 9./16.*1.4142135623730950488+(-45./8.*1.4142135623730950488+105./16.*1.4142135623730950488*xi2)*xi2;
00129     case 6:return 1./32.*4.6904157598234295546*(30.+(-140.+126.*xi2)*xi2)*xi;
00130     case 7:return -5./32.*5.0990195135927848300+(105./32.*5.0990195135927848300+(-315./32.*5.0990195135927848300+231./32.*5.0990195135927848300*xi2)*xi2)*xi2;
00131     case 8:return 1./256.*5.4772255750516611346*(-280.+(2520.+(-5544.+3432.*xi2)*xi2)*xi2)*xi;
00132 
00133     default:
00134       libmesh_error_msg("Invalid shape function index!");
00135     }
00136 
00137   libmesh_error_msg("We'll never get here!");
00138   return 0.;
00139 }
00140 
00141 
00142 
00143 template <>
00144 Real FE<1,SZABAB>::shape_deriv(const Elem* elem,
00145                                const Order order,
00146                                const unsigned int i,
00147                                const unsigned int j,
00148                                const Point& p)
00149 {
00150   libmesh_assert(elem);
00151 
00152   return FE<1,SZABAB>::shape_deriv(elem->type(),
00153                                    static_cast<Order>(order + elem->p_level()), i, j, p);
00154 }
00155 
00156 
00157 
00158 template <>
00159 Real FE<1,SZABAB>::shape_second_deriv(const ElemType,
00160                                       const Order,
00161                                       const unsigned int,
00162                                       const unsigned int,
00163                                       const Point&)
00164 {
00165   static bool warning_given = false;
00166 
00167   if (!warning_given)
00168     libMesh::err << "Second derivatives for Szabab elements "
00169                  << " are not yet implemented!"
00170                  << std::endl;
00171 
00172   warning_given = true;
00173   return 0.;
00174 }
00175 
00176 
00177 
00178 template <>
00179 Real FE<1,SZABAB>::shape_second_deriv(const Elem*,
00180                                       const Order,
00181                                       const unsigned int,
00182                                       const unsigned int,
00183                                       const Point&)
00184 {
00185   static bool warning_given = false;
00186 
00187   if (!warning_given)
00188     libMesh::err << "Second derivatives for Szabab elements "
00189                  << " are not yet implemented!"
00190                  << std::endl;
00191 
00192   warning_given = true;
00193   return 0.;
00194 }
00195 
00196 } // namespace libMesh
00197 
00198 
00199 #endif //LIBMESH_ENABLE_HIGHER_ORDER_SHAPES