libMesh
fe_bernstein.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 
23 #ifdef LIBMESH_ENABLE_HIGHER_ORDER_SHAPES
24 
25 #include "libmesh/elem.h"
26 #include "libmesh/enum_to_string.h"
27 #include "libmesh/fe.h"
28 #include "libmesh/fe_interface.h"
29 #include "libmesh/fe_macro.h"
30 
31 namespace libMesh
32 {
33 
34 // ------------------------------------------------------------
35 // Bernstein-specific implementations, Steffen Petersen 2005
36 
37 // Anonymous namespace for local helper functions
38 namespace {
39 
40 void bernstein_nodal_soln(const Elem * elem,
41  const Order order,
42  const std::vector<Number> & elem_soln,
43  std::vector<Number> & nodal_soln,
44  const bool add_p_level)
45 {
46  const unsigned int n_nodes = elem->n_nodes();
47 
48  const ElemType elem_type = elem->type();
49 
50  nodal_soln.resize(n_nodes);
51 
52  const Order totalorder = order + add_p_level*elem->p_level();
53 
54  // FEType object to be passed to various FEInterface functions below.
55  FEType fe_type(order, BERNSTEIN);
56 
57  switch (totalorder)
58  {
59  // Constant shape functions
60  case CONSTANT:
61  {
62  libmesh_assert_equal_to (elem_soln.size(), 1);
63 
64  std::fill(nodal_soln.begin(), nodal_soln.end(), elem_soln[0]);
65 
66  return;
67  }
68 
69 
70  // For other bases do interpolation at the nodes
71  // explicitly.
72  case FIRST:
73  case SECOND:
74  case THIRD:
75  case FOURTH:
76  case FIFTH:
77  case SIXTH:
78  {
79  const unsigned int n_sf =
80  FEInterface::n_shape_functions(fe_type, elem);
81 
82  std::vector<Point> refspace_nodes;
83  FEBase::get_refspace_nodes(elem_type,refspace_nodes);
84  libmesh_assert_equal_to (refspace_nodes.size(), n_nodes);
85  libmesh_assert_equal_to (elem_soln.size(), n_sf);
86 
87  // Zero before summation
88  std::fill(nodal_soln.begin(), nodal_soln.end(), 0);
89 
90  for (unsigned int n=0; n<n_nodes; n++)
91  // u_i = Sum (alpha_i phi_i)
92  for (unsigned int i=0; i<n_sf; i++)
93  nodal_soln[n] += elem_soln[i] *
94  FEInterface::shape(fe_type, elem, i, refspace_nodes[n]);
95 
96  return;
97  }
98 
99  default:
100  libmesh_error_msg("ERROR: Invalid total order " << totalorder);
101  }
102 } // bernstein_nodal_soln()
103 
104 
105 
106 unsigned int BERNSTEIN_n_dofs(const ElemType t, const Order o)
107 {
108  switch (t)
109  {
110  case NODEELEM:
111  return 1;
112  case EDGE2:
113  case EDGE3:
114  return (o+1);
115  case QUAD4:
116  case QUADSHELL4:
117  libmesh_assert_less (o, 2);
118  libmesh_fallthrough();
119  case QUAD8:
120  case QUADSHELL8:
121  {
122  if (o == 1)
123  return 4;
124  else if (o == 2)
125  return 8;
126  else
127  libmesh_error_msg("ERROR: Invalid Order " << Utility::enum_to_string(o) << " selected for BERNSTEIN FE family!");
128  }
129  case QUAD9:
130  case QUADSHELL9:
131  return ((o+1)*(o+1));
132  case HEX8:
133  libmesh_assert_less (o, 2);
134  libmesh_fallthrough();
135  case HEX20:
136  {
137  if (o == 1)
138  return 8;
139  else if (o == 2)
140  return 20;
141  else
142  libmesh_error_msg("ERROR: Invalid Order " << Utility::enum_to_string(o) << " selected for BERNSTEIN FE family!");
143  }
144  case HEX27:
145  return ((o+1)*(o+1)*(o+1));
146  case TRI3:
147  case TRISHELL3:
148  libmesh_assert_less (o, 2);
149  libmesh_fallthrough();
150  case TRI6:
151  case TRI7:
152  return ((o+1)*(o+2)/2);
153  case TET4:
154  libmesh_assert_less (o, 2);
155  libmesh_fallthrough();
156  case TET10:
157  case TET14:
158  libmesh_assert_less (o, 3);
159  return ((o+1)*(o+2)*(o+3)/6);
160  case INVALID_ELEM:
161  return 0;
162  default:
163  libmesh_error_msg("ERROR: Invalid ElemType " << Utility::enum_to_string(t) << " selected for BERNSTEIN FE family!");
164  }
165 } // BERNSTEIN_n_dofs()
166 
167 
168 
169 unsigned int BERNSTEIN_n_dofs(const Elem * e, const Order o)
170 {
171  libmesh_assert(e);
172  return BERNSTEIN_n_dofs(e->type(), o);
173 }
174 
175 
176 
177 unsigned int BERNSTEIN_n_dofs_at_node(const ElemType t,
178  const Order o,
179  const unsigned int n)
180 {
181  switch (t)
182  {
183  case NODEELEM:
184  return 1;
185 
186  case EDGE2:
187  case EDGE3:
188  switch (n)
189  {
190  case 0:
191  case 1:
192  return 1;
193  case 2:
194  libmesh_assert (o>1);
195  return (o-1);
196  default:
197  libmesh_error_msg("ERROR: Invalid node ID " << n << " selected for EDGE2/3!");
198  }
199  case TRI3:
200  libmesh_assert_less (n, 3);
201  libmesh_assert_less (o, 2);
202  libmesh_fallthrough();
203  case TRI6:
204  // Internal DoFs are associated with the elem on a Tri6, or node 6 on a Tri7
205  libmesh_assert_less (n, 6);
206  libmesh_fallthrough();
207  case TRI7:
208  switch (n)
209  {
210  case 0:
211  case 1:
212  case 2:
213  return 1;
214 
215  case 3:
216  case 4:
217  case 5:
218  return (o-1);
219 
220  case 6:
221  return ((o-1)*(o-2)/2);
222  default:
223  libmesh_error_msg("ERROR: Invalid node ID " << n << " selected for TRI!");
224  }
225  case QUAD4:
226  libmesh_assert_less (n, 4);
227  libmesh_assert_less (o, 2);
228  libmesh_fallthrough();
229  case QUAD8:
230  libmesh_assert_less (n, 8);
231  libmesh_assert_less (o, 3);
232  libmesh_fallthrough();
233  case QUAD9:
234  {
235  switch (n)
236  {
237  case 0:
238  case 1:
239  case 2:
240  case 3:
241  return 1;
242 
243  case 4:
244  case 5:
245  case 6:
246  case 7:
247  return (o-1);
248 
249  case 8:
250  return ((o-1)*(o-1));
251 
252  default:
253  libmesh_error_msg("ERROR: Invalid node ID " << n << " selected for QUAD9!");
254  }
255  }
256  case HEX8:
257  libmesh_assert_less (n, 8);
258  libmesh_assert_less (o, 2);
259  libmesh_fallthrough();
260  case HEX20:
261  libmesh_assert_less (n, 20);
262  libmesh_assert_less (o, 3);
263  libmesh_fallthrough();
264  case HEX27:
265  switch (n)
266  {
267  case 0:
268  case 1:
269  case 2:
270  case 3:
271  case 4:
272  case 5:
273  case 6:
274  case 7:
275  return 1;
276 
277  case 8:
278  case 9:
279  case 10:
280  case 11:
281  case 12:
282  case 13:
283  case 14:
284  case 15:
285  case 16:
286  case 17:
287  case 18:
288  case 19:
289  return (o-1);
290 
291  case 20:
292  case 21:
293  case 22:
294  case 23:
295  case 24:
296  case 25:
297  return ((o-1)*(o-1));
298 
299  case 26:
300  return ((o-1)*(o-1)*(o-1));
301 
302  default:
303  libmesh_error_msg("ERROR: Invalid node ID " << n << " selected for HEX27!");
304  }
305  case TET4:
306  libmesh_assert_less (n, 4);
307  libmesh_assert_less (o, 2);
308  libmesh_fallthrough();
309  case TET10:
310  libmesh_assert_less (n, 10);
311  libmesh_fallthrough();
312  case TET14:
313  libmesh_assert_less (o, 3);
314  switch (n)
315  {
316  case 0:
317  case 1:
318  case 2:
319  case 3:
320  return 1;
321 
322  case 4:
323  case 5:
324  case 6:
325  case 7:
326  case 8:
327  case 9:
328  return (o-1);
329 
330  case 10:
331  case 11:
332  case 12:
333  case 13:
334  return 0;
335 
336  default:
337  libmesh_error_msg("ERROR: Invalid node ID " << n << " selected for TET4/10/14!");
338  }
339  case INVALID_ELEM:
340  return 0;
341  default:
342  libmesh_error_msg("ERROR: Invalid ElemType " << Utility::enum_to_string(t) << " selected for BERNSTEIN FE family!");
343  }
344 } // BERNSTEIN_n_dofs_at_node()
345 
346 
347 
348 unsigned int BERNSTEIN_n_dofs_at_node(const Elem & e,
349  const Order o,
350  const unsigned int n)
351 {
352  return BERNSTEIN_n_dofs_at_node(e.type(), o, n);
353 }
354 
355 
356 
357 unsigned int BERNSTEIN_n_dofs_per_elem(const ElemType t, const Order o)
358 {
359  switch (t)
360  {
361  case NODEELEM:
362  case EDGE2:
363  case EDGE3:
364  return 0;
365  case TRI3:
366  case TRISHELL3:
367  case QUAD4:
368  case QUADSHELL4:
369  return 0;
370  case TRI6:
371  return ((o-1)*(o-2)/2);
372  case TRI7:
373  return 0;
374  case QUAD8:
375  case QUADSHELL8:
376  if (o <= 2)
377  return 0;
378  libmesh_fallthrough();
379  case QUAD9:
380  case QUADSHELL9:
381  return 0;
382  case HEX8:
383  libmesh_assert_less (o, 2);
384  return 0;
385  case HEX20:
386  libmesh_assert_less (o, 3);
387  return 0;
388  case HEX27:
389  return 0;
390  case TET4:
391  libmesh_assert_less (o, 2);
392  libmesh_fallthrough();
393  case TET10:
394  case TET14:
395  libmesh_assert_less (o, 3);
396  return 0;
397  case INVALID_ELEM:
398  return 0;
399  default:
400  libmesh_error_msg("ERROR: Invalid ElemType " << Utility::enum_to_string(t) << " selected for BERNSTEIN FE family!");
401  }
402 } // BERNSTEIN_n_dofs_per_elem
403 
404 
405 
406 unsigned int BERNSTEIN_n_dofs_per_elem(const Elem & e, const Order o)
407 {
408  return BERNSTEIN_n_dofs_per_elem(e.type(), o);
409 }
410 
411 } // anonymous namespace
412 
413 
414 // Instantiate (side_) nodal_soln() function for every dimension
415 LIBMESH_FE_NODAL_SOLN(BERNSTEIN, bernstein_nodal_soln)
417 
418 // Instantiate n_dofs*() functions for every dimension
419 LIBMESH_DEFAULT_NDOFS(BERNSTEIN)
420 
421 // Bernstein FEMs are C^0 continuous
422 template <> FEContinuity FE<0,BERNSTEIN>::get_continuity() const { return C_ZERO; }
423 template <> FEContinuity FE<1,BERNSTEIN>::get_continuity() const { return C_ZERO; }
424 template <> FEContinuity FE<2,BERNSTEIN>::get_continuity() const { return C_ZERO; }
425 template <> FEContinuity FE<3,BERNSTEIN>::get_continuity() const { return C_ZERO; }
426 
427 // Bernstein FEMs are not hierarchic
428 template <> bool FE<0,BERNSTEIN>::is_hierarchic() const { return false; }
429 template <> bool FE<1,BERNSTEIN>::is_hierarchic() const { return false; }
430 template <> bool FE<2,BERNSTEIN>::is_hierarchic() const { return false; }
431 template <> bool FE<3,BERNSTEIN>::is_hierarchic() const { return false; }
432 
433 #ifdef LIBMESH_ENABLE_AMR
434 // compute_constraints() specializations are only needed for 2 and 3D
435 template <>
437  DofMap & dof_map,
438  const unsigned int variable_number,
439  const Elem * elem)
440 { compute_proj_constraints(constraints, dof_map, variable_number, elem); }
441 
442 template <>
444  DofMap & dof_map,
445  const unsigned int variable_number,
446  const Elem * elem)
447 { compute_proj_constraints(constraints, dof_map, variable_number, elem); }
448 #endif // #ifdef LIBMESH_ENABLE_AMR
449 
450 // Bernstein shapes need reinit only for approximation orders >= 3,
451 // but we might reach that with p refinement
452 template <> bool FE<0,BERNSTEIN>::shapes_need_reinit() const { return true; }
453 template <> bool FE<1,BERNSTEIN>::shapes_need_reinit() const { return true; }
454 template <> bool FE<2,BERNSTEIN>::shapes_need_reinit() const { return true; }
455 template <> bool FE<3,BERNSTEIN>::shapes_need_reinit() const { return true; }
456 
457 } // namespace libMesh
458 
459 #endif //LIBMESH_ENABLE_HIGHER_ORDER_SHAPES
ElemType
Defines an enum for geometric element types.
Order
defines an enum for polynomial orders.
Definition: enum_order.h:40
This is the base class from which all geometric element types are derived.
Definition: elem.h:94
virtual bool shapes_need_reinit() const override
The libMesh namespace provides an interface to certain functionality in the library.
LIBMESH_FE_NODAL_SOLN(LIBMESH_FE_SIDE_NODAL_SOLN() LIBMESH_DEFAULT_NDOFS(BERNSTEIN) template<> FEContinuity FE< 0 BERNSTEIN, bernstein_nodal_soln)
Definition: fe_bernstein.C:415
virtual bool is_hierarchic() const override
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:179
LIBMESH_FE_SIDE_NODAL_SOLN(HIERARCHIC_VEC)
const dof_id_type n_nodes
Definition: tecplot_io.C:67
static unsigned int n_shape_functions(const unsigned int dim, const FEType &fe_t, const ElemType t)
Definition: fe_interface.C:272
static Real shape(const unsigned int dim, const FEType &fe_t, const ElemType t, const unsigned int i, const Point &p)
Definition: fe_interface.C:760
static void get_refspace_nodes(const ElemType t, std::vector< Point > &nodes)
Definition: fe_abstract.C:400
libmesh_assert(ctx)
virtual FEContinuity get_continuity() const override
std::string enum_to_string(const T e)
static void compute_constraints(DofConstraints &constraints, DofMap &dof_map, const unsigned int variable_number, const Elem *elem)
Computes the constraint matrix contributions (for non-conforming adapted meshes) corresponding to var...
FEContinuity
defines an enum for finite element types to libmesh_assert a certain level (or type? Hcurl?) of continuity.
The constraint matrix storage format.
Definition: dof_map.h:108