libMesh
fe_rational_shape_2D.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 
69 
70 template <>
72  const Elem * elem,
73  const unsigned int i,
74  const Point & p,
75  const bool add_p_level)
76 {
78  (elem, fet.order, i, p, add_p_level);
79 }
80 
81 
82 template <>
84  const Order order,
85  const unsigned int i,
86  const unsigned int j,
87  const Point & p,
88  const bool add_p_level)
89 {
90  libmesh_assert(elem);
91 
92  FEType underlying_fe_type(order, _underlying_fe_family);
93 
94  return rational_fe_shape_deriv(*elem, underlying_fe_type, i, j, p,
95  add_p_level);
96 }
97 
98 
99 template <>
101  const Order,
102  const unsigned int,
103  const unsigned int,
104  const Point &)
105 {
106  libmesh_error_msg("Rational bases require the real element \nto query nodal weighting.");
107  return 0.;
108 }
109 
110 
111 
112 template <>
114  const Elem * elem,
115  const unsigned int i,
116  const unsigned int j,
117  const Point & p,
118  const bool add_p_level)
119 {
121  (elem, fet.order, i, j, p, add_p_level);
122 }
123 
124 
125 
126 
127 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
128 
129 template <>
131  const Order order,
132  const unsigned int i,
133  const unsigned int j,
134  const Point & p,
135  const bool add_p_level)
136 {
137  libmesh_assert(elem);
138 
139  // FEType object to be passed to various FEInterface functions below.
140  FEType underlying_fe_type(order, _underlying_fe_family);
141 
142  return rational_fe_shape_second_deriv(*elem, underlying_fe_type, i,
143  j, p, add_p_level);
144 }
145 
146 
147 
148 template <>
150  const Order,
151  const unsigned int,
152  const unsigned int,
153  const Point &)
154 {
155  libmesh_error_msg("Rational bases require the real element \nto query nodal weighting.");
156  return 0.;
157 }
158 
159 
160 template <>
162  const Elem * elem,
163  const unsigned int i,
164  const unsigned int j,
165  const Point & p,
166  const bool add_p_level)
167 {
169  (elem, fet.order, i, j, p, add_p_level);
170 }
171 
172 
173 #endif
174 
175 
176 template<>
178  (const Elem * elem,
179  const Order o,
180  const unsigned int i,
181  const std::vector<Point> & p,
182  std::vector<OutputShape> & vi,
183  const bool add_p_level)
184 {
185  libmesh_assert_equal_to(p.size(), vi.size());
186  for (auto j : index_range(vi))
187  vi[j] = FE<2,RATIONAL_BERNSTEIN>::shape (elem, o, i, p[j], add_p_level);
188 }
189 
190 template<>
192  (const Elem * elem,
193  const Order o,
194  const std::vector<Point> & p,
195  std::vector<std::vector<OutputShape>> & v,
196  const bool add_p_level)
197 {
198  FEType underlying_fe_type(o, _underlying_fe_family);
199 
200  rational_all_shapes(*elem, underlying_fe_type, p, v, add_p_level);
201 }
202 
203 
204 template<>
206  (const Elem * elem,
207  const Order o,
208  const unsigned int i,
209  const unsigned int j,
210  const std::vector<Point> & p,
211  std::vector<OutputShape> & v,
212  const bool add_p_level)
213 {
214  libmesh_assert_equal_to(p.size(), v.size());
215  for (auto vi : index_range(v))
216  v[vi] = FE<2,RATIONAL_BERNSTEIN>::shape_deriv (elem, o, i, j, p[vi], add_p_level);
217 }
218 
219 
220 template <>
221 void
223  const Order o,
224  const std::vector<Point> & p,
225  std::vector<std::vector<Real>> * comps[3],
226  const bool add_p_level)
227 {
228  FEType underlying_fe_type(o, _underlying_fe_family);
229 
230  rational_all_shape_derivs (*elem, underlying_fe_type, p,
231  comps, add_p_level);
232 }
233 
234 
235 } // namespace libMesh
236 
237 
238 #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...