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