fe_szabab_shape_1D.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2018 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 {
84  libmesh_assert(elem);
85 
86  return FE<1,SZABAB>::shape(elem->type(), static_cast<Order>(order + elem->p_level()), i, p);
87 }
88 
89 
90 
91 template <>
93  const Order libmesh_dbg_var(order),
94  const unsigned int i,
95  const unsigned int libmesh_dbg_var(j),
96  const Point & p)
97 {
98  // only d()/dxi in 1D!
99  libmesh_assert_equal_to (j, 0);
100 
101  const Real xi = p(0);
102  const Real xi2 = xi*xi;
103 
104  // Use this libmesh_assert rather than a switch with a single entry...
105  // It will go away in optimized mode, essentially has the same effect.
106  libmesh_assert_less_equal (order, SEVENTH);
107 
108  // switch (order)
109  // {
110  // case FIRST:
111  // case SECOND:
112  // case THIRD:
113  // case FOURTH:
114  // case FIFTH:
115  // case SIXTH:
116  // case SEVENTH:
117 
118  switch(i)
119  {
120  case 0:return -1./2.;
121  case 1:return 1./2.;
122  case 2:return 1./2.*2.4494897427831780982*xi;
123  case 3:return -1./4.*3.1622776601683793320+3./4.*3.1622776601683793320*xi2;
124  case 4:return 1./16.*3.7416573867739413856*(-12.+20*xi2)*xi;
125  case 5:return 9./16.*1.4142135623730950488+(-45./8.*1.4142135623730950488+105./16.*1.4142135623730950488*xi2)*xi2;
126  case 6:return 1./32.*4.6904157598234295546*(30.+(-140.+126.*xi2)*xi2)*xi;
127  case 7:return -5./32.*5.0990195135927848300+(105./32.*5.0990195135927848300+(-315./32.*5.0990195135927848300+231./32.*5.0990195135927848300*xi2)*xi2)*xi2;
128  case 8:return 1./256.*5.4772255750516611346*(-280.+(2520.+(-5544.+3432.*xi2)*xi2)*xi2)*xi;
129 
130  default:
131  libmesh_error_msg("Invalid shape function index!");
132  }
133 }
134 
135 
136 
137 template <>
139  const Order order,
140  const unsigned int i,
141  const unsigned int j,
142  const Point & p)
143 {
144  libmesh_assert(elem);
145 
146  return FE<1,SZABAB>::shape_deriv(elem->type(),
147  static_cast<Order>(order + elem->p_level()), i, j, p);
148 }
149 
150 
151 
152 template <>
154  const Order,
155  const unsigned int,
156  const unsigned int,
157  const Point &)
158 {
159  static bool warning_given = false;
160 
161  if (!warning_given)
162  libMesh::err << "Second derivatives for Szabab elements "
163  << " are not yet implemented!"
164  << std::endl;
165 
166  warning_given = true;
167  return 0.;
168 }
169 
170 
171 
172 template <>
174  const Order,
175  const unsigned int,
176  const unsigned int,
177  const Point &)
178 {
179  static bool warning_given = false;
180 
181  if (!warning_given)
182  libMesh::err << "Second derivatives for Szabab elements "
183  << " are not yet implemented!"
184  << std::endl;
185 
186  warning_given = true;
187  return 0.;
188 }
189 
190 } // namespace libMesh
191 
192 
193 #endif //LIBMESH_ENABLE_HIGHER_ORDER_SHAPES
static OutputShape shape(const ElemType t, const Order o, const unsigned int i, const Point &p)
The base class for all geometric element types.
Definition: elem.h:100
static OutputShape shape_deriv(const ElemType t, const Order o, const unsigned int i, const unsigned int j, const Point &p)
unsigned int p_level() const
Definition: elem.h:2555
OStreamProxy err(std::cerr)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual ElemType type() const =0
A geometric point in (x,y,z) space.
Definition: point.h:38
static OutputShape shape_second_deriv(const ElemType t, const Order o, const unsigned int i, const unsigned int j, const Point &p)