libMesh
fe_rational_shape_1D.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 
20 // Local includes
21 #include "libmesh/libmesh_config.h"
22 #ifdef LIBMESH_ENABLE_HIGHER_ORDER_SHAPES
23 
24 #include "libmesh/libmesh_common.h"
25 #include "libmesh/elem.h"
26 #include "libmesh/fe.h"
27 #include "libmesh/fe_interface.h"
28 
29 
30 namespace {
31 using namespace libMesh;
32 
33 static const FEFamily _underlying_fe_family = BERNSTEIN;
34 
35 } // anonymous namespace
36 
37 
38 
39 namespace libMesh
40 {
41 
42 
43 template <>
45  const Order order,
46  const unsigned int i,
47  const Point & p,
48  const bool add_p_level)
49 {
50  libmesh_assert(elem);
51 
52  // FEType object for the non-rational basis underlying this one
53  FEType fe_type(order, _underlying_fe_family);
54 
55  return rational_fe_shape(*elem, fe_type, i, p, add_p_level);
56 }
57 
58 
59 
60 template <>
62  const Order,
63  const unsigned int,
64  const Point &)
65 {
66  libmesh_error_msg("Rational bases require the real element \nto query nodal weighting.");
67  return 0.;
68 }
69 
70 
71 template <>
73  const Elem * elem,
74  const unsigned int i,
75  const Point & p,
76  const bool add_p_level)
77 {
78  libmesh_assert(elem);
79  return FE<1,RATIONAL_BERNSTEIN>::shape(elem, fet.order, i, p, add_p_level);
80 }
81 
82 
83 template <>
85  const Order order,
86  const unsigned int i,
87  const unsigned int j,
88  const Point & p,
89  const bool add_p_level)
90 {
91  libmesh_assert(elem);
92 
93  FEType underlying_fe_type(order, _underlying_fe_family);
94 
95  return rational_fe_shape_deriv(*elem, underlying_fe_type, i, j, p,
96  add_p_level);
97 }
98 
99 
100 template <>
102  const Order,
103  const unsigned int,
104  const unsigned int,
105  const Point &)
106 {
107  libmesh_error_msg("Rational bases require the real element \nto query nodal weighting.");
108  return 0.;
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 {
120  libmesh_assert(elem);
121  return FE<1,RATIONAL_BERNSTEIN>::shape_deriv(elem, fet.order, i, j, p, add_p_level);
122 }
123 
124 
125 
126 
127 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
128 
129 
130 template <>
132  const Order order,
133  const unsigned int i,
134  const unsigned int j,
135  const Point & p,
136  const bool add_p_level)
137 {
138  libmesh_assert(elem);
139 
140  // FEType object to be passed to various FEInterface functions below.
141  FEType underlying_fe_type(order, _underlying_fe_family);
142 
143  return rational_fe_shape_second_deriv(*elem, underlying_fe_type, i,
144  j, p, add_p_level);
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 
161 
162 template <>
164  const Elem * elem,
165  const unsigned int i,
166  const unsigned int j,
167  const Point & p,
168  const bool add_p_level)
169 {
170  libmesh_assert(elem);
171  return FE<1,RATIONAL_BERNSTEIN>::shape_second_deriv(elem, fet.order, i, j, p, add_p_level);
172 }
173 
174 
175 #endif
176 
177 
178 template<>
180  (const Elem * elem,
181  const Order o,
182  const unsigned int i,
183  const std::vector<Point> & p,
184  std::vector<OutputShape> & vi,
185  const bool add_p_level)
186 {
187  libmesh_assert_equal_to(p.size(), vi.size());
188  for (auto j : index_range(vi))
189  vi[j] = FE<1,RATIONAL_BERNSTEIN>::shape (elem, o, i, p[j], add_p_level);
190 }
191 
192 template<>
194  (const Elem * elem,
195  const Order o,
196  const std::vector<Point> & p,
197  std::vector<std::vector<OutputShape>> & v,
198  const bool add_p_level)
199 {
200  FEType underlying_fe_type(o, _underlying_fe_family);
201 
202  rational_all_shapes(*elem, underlying_fe_type, p, v, add_p_level);
203 }
204 
205 
206 template<>
208  (const Elem * elem,
209  const Order o,
210  const unsigned int i,
211  const unsigned int j,
212  const std::vector<Point> & p,
213  std::vector<OutputShape> & v,
214  const bool add_p_level)
215 {
216  libmesh_assert_equal_to(p.size(), v.size());
217  for (auto vi : index_range(v))
218  v[vi] = FE<1,RATIONAL_BERNSTEIN>::shape_deriv (elem, o, i, j, p[vi], add_p_level);
219 }
220 
221 
222 template <>
223 void
225  const Order o,
226  const std::vector<Point> & p,
227  std::vector<std::vector<Real>> * comps[3],
228  const bool add_p_level)
229 {
230  FEType underlying_fe_type(o, _underlying_fe_family);
231 
232  rational_all_shape_derivs (*elem, underlying_fe_type, p,
233  comps, add_p_level);
234 }
235 
236 
237 } // namespace libMesh
238 
239 
240 #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...