libMesh
fe_rational_shape_3D.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2025 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 // Local includes
20 #include "libmesh/libmesh_config.h"
21 #ifdef LIBMESH_ENABLE_HIGHER_ORDER_SHAPES
22 
23 #include "libmesh/elem.h"
24 #include "libmesh/fe.h"
25 #include "libmesh/fe_interface.h"
26 
27 
28 namespace {
29 using namespace libMesh;
30 
31 static const FEFamily _underlying_fe_family = BERNSTEIN;
32 
33 } // anonymous namespace
34 
35 
36 
37 namespace libMesh
38 {
39 
40 
41 template <>
43  const Order order,
44  const unsigned int i,
45  const Point & p,
46  const bool add_p_level)
47 {
48  libmesh_assert(elem);
49 
50  // FEType object for the non-rational basis underlying this one
51  FEType fe_type(order, _underlying_fe_family);
52 
53  return rational_fe_shape(*elem, fe_type, i, p, add_p_level);
54 }
55 
56 
57 
58 template <>
60  const Order,
61  const unsigned int,
62  const Point &)
63 {
64  libmesh_error_msg("Rational bases require the real element \nto query nodal weighting.");
65  return 0.;
66 }
67 
68 template <>
70  const Elem * elem,
71  const unsigned int i,
72  const Point & p,
73  const bool add_p_level)
74 {
76  (elem, fet.order, i, p, add_p_level);
77 }
78 
79 
80 template <>
82  const Order order,
83  const unsigned int i,
84  const unsigned int j,
85  const Point & p,
86  const bool add_p_level)
87 {
88  libmesh_assert(elem);
89 
90  FEType underlying_fe_type(order, _underlying_fe_family);
91 
92  return rational_fe_shape_deriv(*elem, underlying_fe_type, i, j, p,
93  add_p_level);
94 }
95 
96 
97 template <>
99  const Order,
100  const unsigned int,
101  const unsigned int,
102  const Point &)
103 {
104  libmesh_error_msg("Rational bases require the real element \nto query nodal weighting.");
105  return 0.;
106 }
107 
108 
109 template <>
111  const Elem * elem,
112  const unsigned int i,
113  const unsigned int j,
114  const Point & p,
115  const bool add_p_level)
116 {
118  (elem, fet.order, i, j, p, add_p_level);
119 }
120 
121 
122 
123 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
124 
125 template <>
127  const Order order,
128  const unsigned int i,
129  const unsigned int j,
130  const Point & p,
131  const bool add_p_level)
132 {
133  libmesh_assert(elem);
134 
135  // FEType object to be passed to various FEInterface functions below.
136  FEType underlying_fe_type(order, _underlying_fe_family);
137 
138  return rational_fe_shape_second_deriv(*elem, underlying_fe_type, i,
139  j, p, add_p_level);
140 }
141 
142 
143 
144 template <>
146  const Order,
147  const unsigned int,
148  const unsigned int,
149  const Point &)
150 {
151  libmesh_error_msg("Rational bases require the real element \nto query nodal weighting.");
152  return 0.;
153 }
154 
155 
156 template <>
158  const Elem * elem,
159  const unsigned int i,
160  const unsigned int j,
161  const Point & p,
162  const bool add_p_level)
163 {
165  (elem, fet.order, i, j, p, add_p_level);
166 }
167 
168 
169 #endif
170 
171 
172 template<>
174  (const Elem * elem,
175  const Order o,
176  const unsigned int i,
177  const std::vector<Point> & p,
178  std::vector<OutputShape> & vi,
179  const bool add_p_level)
180 {
181  libmesh_assert_equal_to(p.size(), vi.size());
182  for (auto j : index_range(vi))
183  vi[j] = FE<3,RATIONAL_BERNSTEIN>::shape (elem, o, i, p[j], add_p_level);
184 }
185 
186 template<>
188  (const Elem * elem,
189  const Order o,
190  const std::vector<Point> & p,
191  std::vector<std::vector<OutputShape>> & v,
192  const bool add_p_level)
193 {
194  FEType underlying_fe_type(o, _underlying_fe_family);
195 
196  rational_all_shapes(*elem, underlying_fe_type, p, v, add_p_level);
197 }
198 
199 
200 template<>
202  (const Elem * elem,
203  const Order o,
204  const unsigned int i,
205  const unsigned int j,
206  const std::vector<Point> & p,
207  std::vector<OutputShape> & v,
208  const bool add_p_level)
209 {
210  libmesh_assert_equal_to(p.size(), v.size());
211  for (auto vi : index_range(v))
212  v[vi] = FE<3,RATIONAL_BERNSTEIN>::shape_deriv (elem, o, i, j, p[vi], add_p_level);
213 }
214 
215 
216 template <>
217 void
219  const Order o,
220  const std::vector<Point> & p,
221  std::vector<std::vector<Real>> * comps[3],
222  const bool add_p_level)
223 {
224  FEType underlying_fe_type(o, _underlying_fe_family);
225 
226  rational_all_shape_derivs (*elem, underlying_fe_type, p,
227  comps, add_p_level);
228 }
229 
230 
231 } // namespace libMesh
232 
233 
234 #endif// LIBMESH_ENABLE_HIGHER_ORDER_SHAPES
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:196
ElemType
Defines an enum for geometric element types.
Order
defines an enum for polynomial orders.
Definition: enum_order.h:40
static OutputShape shape(const ElemType t, const Order o, const unsigned int i, const Point &p)
static void all_shape_derivs(const Elem *elem, const Order o, const std::vector< Point > &p, std::vector< std::vector< OutputShape >> *comps[3], const bool add_p_level=true)
Fills comps with dphidxi (and in higher dimensions, eta/zeta) derivative component values for all sha...
This is the base class from which all geometric element types are derived.
Definition: elem.h:94
static OutputShape shape_deriv(const ElemType t, const Order o, const unsigned int i, const unsigned int j, const Point &p)
OrderWrapper order
The approximation order of the element.
Definition: fe_type.h:215
The libMesh namespace provides an interface to certain functionality in the library.
A specific instantiation of the FEBase class.
Definition: fe.h:127
libmesh_assert(ctx)
void rational_all_shape_derivs(const Elem &elem, const FEType underlying_fe_type, const std::vector< Point > &p, std::vector< std::vector< OutputShape >> *comps[3], const bool add_p_level)
Definition: fe.C:1338
Real rational_fe_shape(const Elem &elem, const FEType underlying_fe_type, const unsigned int i, const Point &p, const bool add_p_level)
Definition: fe.C:1119
Real rational_fe_shape_second_deriv(const Elem &elem, const FEType underlying_fe_type, const unsigned int i, const unsigned int j, const Point &p, const bool add_p_level)
Definition: fe.C:1202
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
static void shape_derivs(const Elem *elem, const Order o, const unsigned int i, const unsigned int j, const std::vector< Point > &p, std::vector< OutputShape > &v, const bool add_p_level=true)
Fills v with the derivative of the shape function, evaluated at all points p.
static void shapes(const Elem *elem, const Order o, const unsigned int i, const std::vector< Point > &p, std::vector< OutputShape > &v, const bool add_p_level=true)
Fills v with the values of the shape function, evaluated at all points p.
FEFamily
defines an enum for finite element families.
void rational_all_shapes(const Elem &elem, const FEType underlying_fe_type, const std::vector< Point > &p, std::vector< std::vector< Real >> &v, const bool add_p_level)
Definition: fe.C:1308
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:39
auto index_range(const T &sizable)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
Definition: int_range.h:117
static OutputShape shape_second_deriv(const ElemType t, const Order o, const unsigned int i, const unsigned int j, const Point &p)
Real rational_fe_shape_deriv(const Elem &elem, const FEType underlying_fe_type, const unsigned int i, const unsigned int j, const Point &p, const bool add_p_level)
Definition: fe.C:1153
static void all_shapes(const Elem *elem, const Order o, const std::vector< Point > &p, std::vector< std::vector< OutputShape > > &v, const bool add_p_level=true)
Fills v[i][qp] with the values of the shape functions, evaluated at all points in p...