libMesh
fe_xyz_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 // C++ includes
20 
21 // Local includes
22 #include "libmesh/fe.h"
23 #include "libmesh/elem.h"
24 
25 
26 
27 namespace libMesh
28 {
29 
30 
31 template <>
33  const Order,
34  const unsigned int,
35  const Point &)
36 {
37  libmesh_error_msg("XYZ polynomials require the element \n because the centroid is needed.");
38  return 0.;
39 }
40 
41 
42 
43 template <>
44 Real FE<1,XYZ>::shape(const Elem * elem,
45  const Order libmesh_dbg_var(order),
46  const unsigned int i,
47  const Point & point_in,
48  const bool libmesh_dbg_var(add_p_level))
49 {
50  libmesh_assert(elem);
51  libmesh_assert_less_equal (i, order + add_p_level * elem->p_level());
52 
53  Point centroid = elem->centroid();
54  Real max_distance = 0.;
55  for (const Point & p : elem->node_ref_range())
56  {
57  const Real distance = std::abs(centroid(0) - p(0));
58  max_distance = std::max(distance, max_distance);
59  }
60 
61  const Real x = point_in(0);
62  const Real xc = centroid(0);
63  const Real dx = (x - xc)/max_distance;
64 
65  // monomials. since they are hierarchic we only need one case block.
66  switch (i)
67  {
68  case 0:
69  return 1.;
70 
71  case 1:
72  return dx;
73 
74  case 2:
75  return dx*dx;
76 
77  case 3:
78  return dx*dx*dx;
79 
80  case 4:
81  return dx*dx*dx*dx;
82 
83  default:
84  Real val = 1.;
85  for (unsigned int index = 0; index != i; ++index)
86  val *= dx;
87  return val;
88  }
89 }
90 
91 
92 
93 template <>
95  const Order,
96  const unsigned int,
97  const unsigned int,
98  const Point &)
99 {
100  libmesh_error_msg("XYZ polynomials require the element \nbecause the centroid is needed.");
101  return 0.;
102 }
103 
104 
105 
106 template <>
108  const Order libmesh_dbg_var(order),
109  const unsigned int i,
110  const unsigned int libmesh_dbg_var(j),
111  const Point & point_in,
112  const bool libmesh_dbg_var(add_p_level))
113 {
114  libmesh_assert(elem);
115  libmesh_assert_less_equal (i, order + add_p_level * elem->p_level());
116 
117  // only d()/dxi in 1D!
118 
119  libmesh_assert_equal_to (j, 0);
120 
121  Point centroid = elem->centroid();
122  Real max_distance = 0.;
123  for (const Point & p : elem->node_ref_range())
124  {
125  const Real distance = std::abs(centroid(0) - p(0));
126  max_distance = std::max(distance, max_distance);
127  }
128 
129  const Real x = point_in(0);
130  const Real xc = centroid(0);
131  const Real dx = (x - xc)/max_distance;
132 
133  // monomials. since they are hierarchic we only need one case block.
134  switch (i)
135  {
136  case 0:
137  return 0.;
138 
139  case 1:
140  return 1./max_distance;
141 
142  case 2:
143  return 2.*dx/max_distance;
144 
145  case 3:
146  return 3.*dx*dx/max_distance;
147 
148  case 4:
149  return 4.*dx*dx*dx/max_distance;
150 
151  default:
152  Real val = i;
153  for (unsigned int index = 1; index != i; ++index)
154  val *= dx;
155  return val/max_distance;
156  }
157 }
158 
159 
160 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
161 
162 template <>
164  const Order,
165  const unsigned int,
166  const unsigned int,
167  const Point &)
168 {
169  libmesh_error_msg("XYZ polynomials require the element \nbecause the centroid is needed.");
170  return 0.;
171 }
172 
173 
174 
175 template <>
177  const Order libmesh_dbg_var(order),
178  const unsigned int i,
179  const unsigned int libmesh_dbg_var(j),
180  const Point & point_in,
181  const bool libmesh_dbg_var(add_p_level))
182 {
183  libmesh_assert(elem);
184  libmesh_assert_less_equal (i, order + add_p_level * elem->p_level());
185 
186  // only d2()/dxi2 in 1D!
187 
188  libmesh_assert_equal_to (j, 0);
189 
190  Point centroid = elem->centroid();
191  Real max_distance = 0.;
192  for (const Point & p : elem->node_ref_range())
193  {
194  const Real distance = std::abs(centroid(0) - p(0));
195  max_distance = std::max(distance, max_distance);
196  }
197 
198  const Real x = point_in(0);
199  const Real xc = centroid(0);
200  const Real dx = (x - xc)/max_distance;
201  const Real dist2 = pow(max_distance,2.);
202 
203  // monomials. since they are hierarchic we only need one case block.
204  switch (i)
205  {
206  case 0:
207  case 1:
208  return 0.;
209 
210  case 2:
211  return 2./dist2;
212 
213  case 3:
214  return 6.*dx/dist2;
215 
216  case 4:
217  return 12.*dx*dx/dist2;
218 
219  default:
220  Real val = 2.;
221  for (unsigned int index = 2; index != i; ++index)
222  val *= (index+1) * dx;
223  return val/dist2;
224  }
225 }
226 #endif
227 
228 } // namespace libMesh
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::Elem::centroid
virtual Point centroid() const
Definition: elem.C:345
libMesh::libmesh_assert
libmesh_assert(ctx)
std::abs
MetaPhysicL::DualNumber< T, D > abs(const MetaPhysicL::DualNumber< T, D > &in)
libMesh::FE::shape
static OutputShape shape(const ElemType t, const Order o, const unsigned int i, const Point &p)
libMesh::Elem::node_ref_range
SimpleRange< NodeRefIter > node_ref_range()
Returns a range with all nodes of an element, usable in range-based for loops.
Definition: elem.h:2152
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
std::pow
double pow(double a, int b)
Definition: libmesh_augment_std_namespace.h:58
distance
Real distance(const Point &p)
Definition: subdomains_ex3.C:50
libMesh::Elem
This is the base class from which all geometric element types are derived.
Definition: elem.h:100
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::ElemType
ElemType
Defines an enum for geometric element types.
Definition: enum_elem_type.h:33