libMesh
fe_szabab_shape_1D.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2019 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 
20 // C++ includes
21 
22 // Local includes
23 #include "libmesh/libmesh_config.h"
24 
25 #ifdef LIBMESH_ENABLE_HIGHER_ORDER_SHAPES
26 
27 #include "libmesh/fe.h"
28 #include "libmesh/elem.h"
29 
30 namespace libMesh
31 {
32 
33 
34 template <>
36  const Order libmesh_dbg_var(order),
37  const unsigned int i,
38  const Point & p)
39 {
40  const Real xi = p(0);
41  const Real xi2 = xi*xi;
42 
43 
44  // Use this libmesh_assert rather than a switch with a single entry...
45  // It will go away in optimized mode, essentially has the same effect.
46  libmesh_assert_less_equal (order, SEVENTH);
47 
48  // switch (order)
49  // {
50  // case FIRST:
51  // case SECOND:
52  // case THIRD:
53  // case FOURTH:
54  // case FIFTH:
55  // case SIXTH:
56  // case SEVENTH:
57 
58  switch(i)
59  {
60  //nodal shape functions
61  case 0: return 1./2.-1./2.*xi;
62  case 1: return 1./2.+1./2.*xi;
63  case 2: return 1./4. *2.4494897427831780982*(xi2-1.);
64  case 3: return 1./4. *3.1622776601683793320*(xi2-1.)*xi;
65  case 4: return 1./16. *3.7416573867739413856*((5.*xi2-6.)*xi2+1.);
66  case 5: return 3./16. *1.4142135623730950488*(3.+(-10.+7.*xi2)*xi2)*xi;
67  case 6: return 1./32. *4.6904157598234295546*(-1.+(15.+(-35.+21.*xi2)*xi2)*xi2);
68  case 7: return 1./32. *5.0990195135927848300*(-5.+(35.+(-63.+33.*xi2)*xi2)*xi2)*xi;
69  case 8: return 1./256.*5.4772255750516611346*(5.+(-140.+(630.+(-924.+429.*xi2)*xi2)*xi2)*xi2);
70 
71  default:
72  libmesh_error_msg("Invalid shape function index!");
73  }
74 }
75 
76 
77 
78 template <>
80  const Order order,
81  const unsigned int i,
82  const Point & p,
83  const bool add_p_level)
84 {
85  libmesh_assert(elem);
86 
87  return FE<1,SZABAB>::shape(elem->type(), static_cast<Order>(order + add_p_level * add_p_level * elem->p_level()), i, p);
88 }
89 
90 
91 
92 template <>
94  const Order libmesh_dbg_var(order),
95  const unsigned int i,
96  const unsigned int libmesh_dbg_var(j),
97  const Point & p)
98 {
99  // only d()/dxi in 1D!
100  libmesh_assert_equal_to (j, 0);
101 
102  const Real xi = p(0);
103  const Real xi2 = xi*xi;
104 
105  // Use this libmesh_assert rather than a switch with a single entry...
106  // It will go away in optimized mode, essentially has the same effect.
107  libmesh_assert_less_equal (order, SEVENTH);
108 
109  // switch (order)
110  // {
111  // case FIRST:
112  // case SECOND:
113  // case THIRD:
114  // case FOURTH:
115  // case FIFTH:
116  // case SIXTH:
117  // case SEVENTH:
118 
119  switch(i)
120  {
121  case 0:return -1./2.;
122  case 1:return 1./2.;
123  case 2:return 1./2.*2.4494897427831780982*xi;
124  case 3:return -1./4.*3.1622776601683793320+3./4.*3.1622776601683793320*xi2;
125  case 4:return 1./16.*3.7416573867739413856*(-12.+20*xi2)*xi;
126  case 5:return 9./16.*1.4142135623730950488+(-45./8.*1.4142135623730950488+105./16.*1.4142135623730950488*xi2)*xi2;
127  case 6:return 1./32.*4.6904157598234295546*(30.+(-140.+126.*xi2)*xi2)*xi;
128  case 7:return -5./32.*5.0990195135927848300+(105./32.*5.0990195135927848300+(-315./32.*5.0990195135927848300+231./32.*5.0990195135927848300*xi2)*xi2)*xi2;
129  case 8:return 1./256.*5.4772255750516611346*(-280.+(2520.+(-5544.+3432.*xi2)*xi2)*xi2)*xi;
130 
131  default:
132  libmesh_error_msg("Invalid shape function index!");
133  }
134 }
135 
136 
137 
138 template <>
140  const Order order,
141  const unsigned int i,
142  const unsigned int j,
143  const Point & p,
144  const bool add_p_level)
145 {
146  libmesh_assert(elem);
147 
148  return FE<1,SZABAB>::shape_deriv(elem->type(),
149  static_cast<Order>(order + add_p_level * elem->p_level()), i, j, p);
150 }
151 
152 
153 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
154 
155 template <>
157  const Order,
158  const unsigned int,
159  const unsigned int,
160  const Point &)
161 {
162  static bool warning_given = false;
163 
164  if (!warning_given)
165  libMesh::err << "Second derivatives for Szabab elements "
166  << " are not yet implemented!"
167  << std::endl;
168 
169  warning_given = true;
170  return 0.;
171 }
172 
173 
174 
175 template <>
177  const Order,
178  const unsigned int,
179  const unsigned int,
180  const Point &,
181  const bool)
182 {
183  static bool warning_given = false;
184 
185  if (!warning_given)
186  libMesh::err << "Second derivatives for Szabab elements "
187  << " are not yet implemented!"
188  << std::endl;
189 
190  warning_given = true;
191  return 0.;
192 }
193 #endif
194 
195 } // namespace libMesh
196 
197 
198 #endif //LIBMESH_ENABLE_HIGHER_ORDER_SHAPES
libMesh::FE::shape_second_deriv
static OutputShape shape_second_deriv(const ElemType t, const Order o, const unsigned int i, const unsigned int j, const Point &p)
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::Order
Order
Definition: enum_order.h:40
libMesh::Elem::p_level
unsigned int p_level() const
Definition: elem.h:2512
libMesh::libmesh_assert
libmesh_assert(ctx)
libMesh::FE::shape
static OutputShape shape(const ElemType t, const Order o, const unsigned int i, const Point &p)
libMesh::FE::shape_deriv
static OutputShape shape_deriv(const ElemType t, const Order o, const unsigned int i, const unsigned int j, const Point &p)
libMesh::Point
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:38
libMesh::SEVENTH
Definition: enum_order.h:48
libMesh::Elem
This is the base class from which all geometric element types are derived.
Definition: elem.h:100
libMesh::err
OStreamProxy err
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::Elem::type
virtual ElemType type() const =0
libMesh::ElemType
ElemType
Defines an enum for geometric element types.
Definition: enum_elem_type.h:33