libMesh
fe_raviart.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/dof_map.h"
22 #include "libmesh/elem.h"
23 #include "libmesh/enum_to_string.h"
24 #include "libmesh/fe.h"
25 #include "libmesh/fe_interface.h"
26 #include "libmesh/fe_macro.h"
27 #include "libmesh/tensor_value.h"
28 
29 
30 namespace libMesh
31 {
32 
33 
42 
43 
44 // Anonymous namespace for local helper functions
45 namespace {
46 void raviart_thomas_nodal_soln(const Elem * elem,
47  const Order order,
48  const std::vector<Number> & elem_soln,
49  const int dim,
50  const int vdim,
51  std::vector<Number> & nodal_soln,
52  const bool add_p_level)
53 {
54  const unsigned int n_nodes = elem->n_nodes();
55  const ElemType elem_type = elem->type();
56 
57  const Order totalorder = order + add_p_level*elem->p_level();
58 
59  nodal_soln.resize(n_nodes*vdim);
60 
61  FEType p_refined_fe_type(totalorder, RAVIART_THOMAS);
62 
63  if (elem_type != TRI6 && elem_type != TRI7 &&
64  elem_type != QUAD8 && elem_type != QUAD9 &&
65  elem_type != TET14 && elem_type != HEX27)
66  libmesh_error_msg("ERROR: Invalid ElemType " << Utility::enum_to_string(elem_type) << " selected for RAVIART_THOMAS FE family!");
67 
68  const unsigned int n_sf = FEInterface::n_shape_functions(p_refined_fe_type, elem, false);
69 
70  std::vector<Point> refspace_nodes;
71  FEVectorBase::get_refspace_nodes(elem_type,refspace_nodes);
72  libmesh_assert_equal_to (refspace_nodes.size(), n_nodes);
73  libmesh_assert_equal_to (elem_soln.size(), n_sf);
74 
75  // Need to create new fe object so the shape function has the FETransformation
76  // applied to it.
77  std::unique_ptr<FEVectorBase> vis_fe = FEVectorBase::build(dim, p_refined_fe_type);
78 
79  const std::vector<std::vector<RealGradient>> & vis_phi = vis_fe->get_phi();
80 
81  vis_fe->reinit(elem,&refspace_nodes);
82 
83  // Zero before summation
84  std::fill(nodal_soln.begin(), nodal_soln.end(), 0);
85 
86  for (unsigned int n = 0; n < n_nodes; n++)
87  // u = Sum (u_i phi_i)
88  for (unsigned int i=0; i<n_sf; i++)
89  for (int d = 0; d < vdim; d++)
90  nodal_soln[vdim*n+d] += elem_soln[i]*(vis_phi[i][n](d));
91 
92  return;
93 } // raviart_thomas_nodal_soln
94 
95 
96 
97 unsigned int raviart_thomas_n_dofs(const ElemType t, const Order o)
98 {
99  libmesh_assert_greater (o, 0);
100  switch (t)
101  {
102  case TRI6:
103  case TRI7:
104  return o*(o+2);
105  case QUAD8:
106  case QUAD9:
107  return 2*o*(o+1);
108  case TET14:
109  return o*(o+1)*(o+3)/2;
110  case HEX27:
111  return 3*o*o*(o+1);
112  case INVALID_ELEM:
113  return 0;
114  default:
115  libmesh_error_msg("ERROR: Invalid ElemType " << Utility::enum_to_string(t) << " selected for RAVIART_THOMAS FE family!");
116  }
117 }
118 
119 
120 
121 unsigned int raviart_thomas_n_dofs(const Elem * e, const Order o)
122 {
123  libmesh_assert(e);
124  return raviart_thomas_n_dofs(e->type(), o);
125 }
126 
127 
128 
129 unsigned int raviart_thomas_n_dofs_at_node(const ElemType t,
130  const Order o,
131  const unsigned int n)
132 {
133  libmesh_assert_greater (o, 0);
134  switch (t)
135  {
136  case TRI6:
137  case TRI7:
138  {
139  switch (n)
140  {
141  case 0:
142  case 1:
143  case 2:
144  return 0;
145  case 3:
146  case 4:
147  case 5:
148  return o;
149  case 6:
150  libmesh_assert_equal_to(t, TRI7);
151  return 0;
152  default:
153  libmesh_error_msg("ERROR: Invalid node ID " << n);
154  }
155  }
156  case QUAD8:
157  case QUAD9:
158  {
159  switch (n)
160  {
161  case 0:
162  case 1:
163  case 2:
164  case 3:
165  return 0;
166  case 4:
167  case 5:
168  case 6:
169  case 7:
170  return o;
171  case 8:
172  libmesh_assert_equal_to(t, QUAD9);
173  return 0;
174  default:
175  libmesh_error_msg("ERROR: Invalid node ID " << n);
176  }
177  }
178  case TET14:
179  {
180  switch (n)
181  {
182  case 0:
183  case 1:
184  case 2:
185  case 3:
186  case 4:
187  case 5:
188  case 6:
189  case 7:
190  case 8:
191  case 9:
192  return 0;
193  case 10:
194  case 11:
195  case 12:
196  case 13:
197  return o*(o+1)/2;
198  default:
199  libmesh_error_msg("ERROR: Invalid node ID " << n);
200  }
201  }
202  case HEX27:
203  {
204  switch (n)
205  {
206  case 0:
207  case 1:
208  case 2:
209  case 3:
210  case 4:
211  case 5:
212  case 6:
213  case 7:
214  case 8:
215  case 9:
216  case 10:
217  case 11:
218  case 12:
219  case 13:
220  case 14:
221  case 15:
222  case 16:
223  case 17:
224  case 18:
225  case 19:
226  return 0;
227  case 20:
228  case 21:
229  case 22:
230  case 23:
231  case 24:
232  case 25:
233  return o*o;
234  case 26:
235  return 0;
236  default:
237  libmesh_error_msg("ERROR: Invalid node ID " << n);
238  }
239  }
240  case INVALID_ELEM:
241  return 0;
242  default:
243  libmesh_error_msg("ERROR: Invalid ElemType " << Utility::enum_to_string(t) << " selected for RAVIART_THOMAS FE family!");
244  }
245 }
246 
247 
248 
249 unsigned int raviart_thomas_n_dofs_at_node(const Elem & e,
250  const Order o,
251  const unsigned int n)
252 {
253  return raviart_thomas_n_dofs_at_node(e.type(), o, n);
254 }
255 
256 
257 
258 unsigned int raviart_thomas_n_dofs_per_elem(const ElemType t,
259  const Order o)
260 {
261  libmesh_assert_greater (o, 0);
262  switch (t)
263  {
264  case TRI6:
265  case TRI7:
266  return o*(o-1);
267  case QUAD8:
268  case QUAD9:
269  return 2*o*(o-1);
270  case TET14:
271  return (o+1)*o*(o-1)/2;
272  case HEX27:
273  return 3*o*o*(o-1);
274  case INVALID_ELEM:
275  return 0;
276  default:
277  libmesh_error_msg("ERROR: Invalid ElemType " << Utility::enum_to_string(t) << " selected for RAVIART_THOMAS FE family!");
278  }
279 }
280 
281 
282 
283 unsigned int raviart_thomas_n_dofs_per_elem(const Elem & e,
284  const Order o)
285 {
286  return raviart_thomas_n_dofs_per_elem(e.type(), o);
287 }
288 
289 
290 #ifdef LIBMESH_ENABLE_AMR
291 void raviart_thomas_compute_constraints (DofConstraints & /*constraints*/,
292  DofMap & /*dof_map*/,
293  const unsigned int /*variable_number*/,
294  const Elem * libmesh_dbg_var(elem),
295  const unsigned Dim)
296 {
297  // Only constrain elements in 2,3D.
298  if (Dim == 1)
299  return;
300 
301  libmesh_assert(elem);
302 
303  libmesh_not_implemented();
304 } // raviart_thomas_compute_constraints()
305 #endif // #ifdef LIBMESH_ENABLE_AMR
306 
307 } // anonymous namespace
308 
309 #define RAVIART_LOW_D_ERROR_MESSAGE \
310  libmesh_error_msg("ERROR: This method makes no sense for low-D elements!");
311 
312 
313 // Do full-specialization for every dimension, instead
314 // of explicit instantiation at the end of this file.
315 template <>
317  const Order,
318  const std::vector<Number> &,
319  std::vector<Number> &,
320  bool,
321  const unsigned)
322 { RAVIART_LOW_D_ERROR_MESSAGE }
323 
324 template <>
326  const Order,
327  const std::vector<Number> &,
328  std::vector<Number> &,
329  bool,
330  const unsigned)
331 { RAVIART_LOW_D_ERROR_MESSAGE }
332 
333 template <>
335  const Order order,
336  const std::vector<Number> & elem_soln,
337  std::vector<Number> & nodal_soln,
338  const bool add_p_level,
339  const unsigned vdim)
340 { raviart_thomas_nodal_soln(elem, order, elem_soln, 2 /*dim*/, vdim, nodal_soln, add_p_level); }
341 
342 template <>
344  const Order order,
345  const std::vector<Number> & elem_soln,
346  std::vector<Number> & nodal_soln,
347  const bool add_p_level,
348  const unsigned)
349 { raviart_thomas_nodal_soln(elem, order, elem_soln, 3 /*dim*/, 3 /*vdim*/, nodal_soln, add_p_level); }
350 
351 template <>
353  const Order,
354  const std::vector<Number> &,
355  std::vector<Number> &,
356  bool,
357  const unsigned)
358 { RAVIART_LOW_D_ERROR_MESSAGE }
359 
360 template <>
362  const Order,
363  const std::vector<Number> &,
364  std::vector<Number> &,
365  bool,
366  const unsigned)
367 { RAVIART_LOW_D_ERROR_MESSAGE }
368 
369 template <>
371  const Order order,
372  const std::vector<Number> & elem_soln,
373  std::vector<Number> & nodal_soln,
374  const bool add_p_level,
375  const unsigned vdim)
376 { raviart_thomas_nodal_soln(elem, order, elem_soln, 2 /*dim*/, vdim, nodal_soln, add_p_level); }
377 
378 template <>
380  const Order order,
381  const std::vector<Number> & elem_soln,
382  std::vector<Number> & nodal_soln,
383  const bool add_p_level,
384  const unsigned)
385 { raviart_thomas_nodal_soln(elem, order, elem_soln, 3 /*dim*/, 3 /*vdim*/, nodal_soln, add_p_level); }
386 
389 
390 
391 // Do full-specialization for every dimension, instead
392 // of explicit instantiation at the end of this function.
393 // This could be macro-ified.
394 template <> unsigned int FE<0,RAVIART_THOMAS>::n_dofs(const ElemType, const Order) { RAVIART_LOW_D_ERROR_MESSAGE }
395 template <> unsigned int FE<1,RAVIART_THOMAS>::n_dofs(const ElemType, const Order) { RAVIART_LOW_D_ERROR_MESSAGE }
396 template <> unsigned int FE<2,RAVIART_THOMAS>::n_dofs(const ElemType t, const Order o) { return raviart_thomas_n_dofs(t, o); }
397 template <> unsigned int FE<3,RAVIART_THOMAS>::n_dofs(const ElemType t, const Order o) { return raviart_thomas_n_dofs(t, o); }
398 
399 template <> unsigned int FE<0,RAVIART_THOMAS>::n_dofs(const Elem *, const Order) { RAVIART_LOW_D_ERROR_MESSAGE }
400 template <> unsigned int FE<1,RAVIART_THOMAS>::n_dofs(const Elem *, const Order) { RAVIART_LOW_D_ERROR_MESSAGE }
401 template <> unsigned int FE<2,RAVIART_THOMAS>::n_dofs(const Elem * e, const Order o) { return raviart_thomas_n_dofs(e, o); }
402 template <> unsigned int FE<3,RAVIART_THOMAS>::n_dofs(const Elem * e, const Order o) { return raviart_thomas_n_dofs(e, o); }
403 
404 template <> unsigned int FE<0,L2_RAVIART_THOMAS>::n_dofs(const ElemType, const Order) { RAVIART_LOW_D_ERROR_MESSAGE }
405 template <> unsigned int FE<1,L2_RAVIART_THOMAS>::n_dofs(const ElemType, const Order) { RAVIART_LOW_D_ERROR_MESSAGE }
406 template <> unsigned int FE<2,L2_RAVIART_THOMAS>::n_dofs(const ElemType t, const Order o) { return raviart_thomas_n_dofs(t, o); }
407 template <> unsigned int FE<3,L2_RAVIART_THOMAS>::n_dofs(const ElemType t, const Order o) { return raviart_thomas_n_dofs(t, o); }
408 
409 template <> unsigned int FE<0,L2_RAVIART_THOMAS>::n_dofs(const Elem *, const Order) { RAVIART_LOW_D_ERROR_MESSAGE }
410 template <> unsigned int FE<1,L2_RAVIART_THOMAS>::n_dofs(const Elem *, const Order) { RAVIART_LOW_D_ERROR_MESSAGE }
411 template <> unsigned int FE<2,L2_RAVIART_THOMAS>::n_dofs(const Elem * e, const Order o) { return raviart_thomas_n_dofs(e, o); }
412 template <> unsigned int FE<3,L2_RAVIART_THOMAS>::n_dofs(const Elem * e, const Order o) { return raviart_thomas_n_dofs(e, o); }
413 
414 template <> unsigned int FE<0,RAVIART_THOMAS>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { RAVIART_LOW_D_ERROR_MESSAGE }
415 template <> unsigned int FE<1,RAVIART_THOMAS>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { RAVIART_LOW_D_ERROR_MESSAGE }
416 template <> unsigned int FE<2,RAVIART_THOMAS>::n_dofs_at_node(const ElemType t, const Order o, const unsigned int n) { return raviart_thomas_n_dofs_at_node(t, o, n); }
417 template <> unsigned int FE<3,RAVIART_THOMAS>::n_dofs_at_node(const ElemType t, const Order o, const unsigned int n) { return raviart_thomas_n_dofs_at_node(t, o, n); }
418 
419 template <> unsigned int FE<0,RAVIART_THOMAS>::n_dofs_at_node(const Elem &, const Order, const unsigned int) { RAVIART_LOW_D_ERROR_MESSAGE }
420 template <> unsigned int FE<1,RAVIART_THOMAS>::n_dofs_at_node(const Elem &, const Order, const unsigned int) { RAVIART_LOW_D_ERROR_MESSAGE }
421 template <> unsigned int FE<2,RAVIART_THOMAS>::n_dofs_at_node(const Elem & e, const Order o, const unsigned int n) { return raviart_thomas_n_dofs_at_node(e, o, n); }
422 template <> unsigned int FE<3,RAVIART_THOMAS>::n_dofs_at_node(const Elem & e, const Order o, const unsigned int n) { return raviart_thomas_n_dofs_at_node(e, o, n); }
423 
424 template <> unsigned int FE<0,L2_RAVIART_THOMAS>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { RAVIART_LOW_D_ERROR_MESSAGE }
425 template <> unsigned int FE<1,L2_RAVIART_THOMAS>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { RAVIART_LOW_D_ERROR_MESSAGE }
426 template <> unsigned int FE<2,L2_RAVIART_THOMAS>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { return 0; }
427 template <> unsigned int FE<3,L2_RAVIART_THOMAS>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { return 0; }
428 
429 template <> unsigned int FE<0,L2_RAVIART_THOMAS>::n_dofs_at_node(const Elem &, const Order, const unsigned int) { RAVIART_LOW_D_ERROR_MESSAGE }
430 template <> unsigned int FE<1,L2_RAVIART_THOMAS>::n_dofs_at_node(const Elem &, const Order, const unsigned int) { RAVIART_LOW_D_ERROR_MESSAGE }
431 template <> unsigned int FE<2,L2_RAVIART_THOMAS>::n_dofs_at_node(const Elem &, const Order, const unsigned int) { return 0; }
432 template <> unsigned int FE<3,L2_RAVIART_THOMAS>::n_dofs_at_node(const Elem &, const Order, const unsigned int) { return 0; }
433 
434 template <> unsigned int FE<0,RAVIART_THOMAS>::n_dofs_per_elem(const ElemType, const Order) { RAVIART_LOW_D_ERROR_MESSAGE }
435 template <> unsigned int FE<1,RAVIART_THOMAS>::n_dofs_per_elem(const ElemType, const Order) { RAVIART_LOW_D_ERROR_MESSAGE }
436 template <> unsigned int FE<2,RAVIART_THOMAS>::n_dofs_per_elem(const ElemType t, const Order o) { return raviart_thomas_n_dofs_per_elem(t, o); }
437 template <> unsigned int FE<3,RAVIART_THOMAS>::n_dofs_per_elem(const ElemType t, const Order o) { return raviart_thomas_n_dofs_per_elem(t, o); }
438 
439 template <> unsigned int FE<0,RAVIART_THOMAS>::n_dofs_per_elem(const Elem &, const Order) { RAVIART_LOW_D_ERROR_MESSAGE }
440 template <> unsigned int FE<1,RAVIART_THOMAS>::n_dofs_per_elem(const Elem &, const Order) { RAVIART_LOW_D_ERROR_MESSAGE }
441 template <> unsigned int FE<2,RAVIART_THOMAS>::n_dofs_per_elem(const Elem & e, const Order o) { return raviart_thomas_n_dofs_per_elem(e, o); }
442 template <> unsigned int FE<3,RAVIART_THOMAS>::n_dofs_per_elem(const Elem & e, const Order o) { return raviart_thomas_n_dofs_per_elem(e, o); }
443 
444 // L2 Raviart-Thomas elements have all their dofs per element
445 template <> unsigned int FE<0,L2_RAVIART_THOMAS>::n_dofs_per_elem(const ElemType, const Order) { RAVIART_LOW_D_ERROR_MESSAGE }
446 template <> unsigned int FE<1,L2_RAVIART_THOMAS>::n_dofs_per_elem(const ElemType, const Order) { RAVIART_LOW_D_ERROR_MESSAGE }
447 template <> unsigned int FE<2,L2_RAVIART_THOMAS>::n_dofs_per_elem(const ElemType t, const Order o) { return n_dofs(t, o); }
448 template <> unsigned int FE<3,L2_RAVIART_THOMAS>::n_dofs_per_elem(const ElemType t, const Order o) { return n_dofs(t, o); }
449 
450 template <> unsigned int FE<0,L2_RAVIART_THOMAS>::n_dofs_per_elem(const Elem &, const Order) { RAVIART_LOW_D_ERROR_MESSAGE }
451 template <> unsigned int FE<1,L2_RAVIART_THOMAS>::n_dofs_per_elem(const Elem &, const Order) { RAVIART_LOW_D_ERROR_MESSAGE }
452 template <> unsigned int FE<2,L2_RAVIART_THOMAS>::n_dofs_per_elem(const Elem & e, const Order o) { return n_dofs(&e, o); }
453 template <> unsigned int FE<3,L2_RAVIART_THOMAS>::n_dofs_per_elem(const Elem & e, const Order o) { return n_dofs(&e, o); }
454 
455 // Raviart-Thomas FEMs are always normally continuous
456 template <> FEContinuity FE<0,RAVIART_THOMAS>::get_continuity() const { RAVIART_LOW_D_ERROR_MESSAGE }
457 template <> FEContinuity FE<1,RAVIART_THOMAS>::get_continuity() const { RAVIART_LOW_D_ERROR_MESSAGE }
460 
461 // L2 Raviart-Thomas FEMs are discontinuous
462 template <> FEContinuity FE<0,L2_RAVIART_THOMAS>::get_continuity() const { RAVIART_LOW_D_ERROR_MESSAGE }
463 template <> FEContinuity FE<1,L2_RAVIART_THOMAS>::get_continuity() const { RAVIART_LOW_D_ERROR_MESSAGE }
466 
467 // Raviart-Thomas FEMs are not hierarchic
468 template <> bool FE<0,RAVIART_THOMAS>::is_hierarchic() const { RAVIART_LOW_D_ERROR_MESSAGE }
469 template <> bool FE<1,RAVIART_THOMAS>::is_hierarchic() const { RAVIART_LOW_D_ERROR_MESSAGE }
470 template <> bool FE<2,RAVIART_THOMAS>::is_hierarchic() const { return false; }
471 template <> bool FE<3,RAVIART_THOMAS>::is_hierarchic() const { return false; }
472 template <> bool FE<0,L2_RAVIART_THOMAS>::is_hierarchic() const { RAVIART_LOW_D_ERROR_MESSAGE }
473 template <> bool FE<1,L2_RAVIART_THOMAS>::is_hierarchic() const { RAVIART_LOW_D_ERROR_MESSAGE }
474 template <> bool FE<2,L2_RAVIART_THOMAS>::is_hierarchic() const { return false; }
475 template <> bool FE<3,L2_RAVIART_THOMAS>::is_hierarchic() const { return false; }
476 
477 // Raviart-Thomas FEM shapes always need to be reinit'ed (because of orientation dependence)
478 template <> bool FE<0,RAVIART_THOMAS>::shapes_need_reinit() const { RAVIART_LOW_D_ERROR_MESSAGE }
479 template <> bool FE<1,RAVIART_THOMAS>::shapes_need_reinit() const { RAVIART_LOW_D_ERROR_MESSAGE }
480 template <> bool FE<2,RAVIART_THOMAS>::shapes_need_reinit() const { return true; }
481 template <> bool FE<3,RAVIART_THOMAS>::shapes_need_reinit() const { return true; }
482 
483 // Do L2 Raviart-Thomas FEM shapes always need to be reinit'ed because of orientation dependence?
484 template <> bool FE<0,L2_RAVIART_THOMAS>::shapes_need_reinit() const { RAVIART_LOW_D_ERROR_MESSAGE }
485 template <> bool FE<1,L2_RAVIART_THOMAS>::shapes_need_reinit() const { RAVIART_LOW_D_ERROR_MESSAGE }
486 template <> bool FE<2,L2_RAVIART_THOMAS>::shapes_need_reinit() const { return true; }
487 template <> bool FE<3,L2_RAVIART_THOMAS>::shapes_need_reinit() const { return true; }
488 
489 #ifdef LIBMESH_ENABLE_AMR
490 template <>
492  DofMap &,
493  const unsigned int,
494  const Elem *)
495 { RAVIART_LOW_D_ERROR_MESSAGE }
496 
497 template <>
499  DofMap &,
500  const unsigned int,
501  const Elem *)
502 { RAVIART_LOW_D_ERROR_MESSAGE }
503 
504 template <>
506  DofMap & dof_map,
507  const unsigned int variable_number,
508  const Elem * elem)
509 { raviart_thomas_compute_constraints(constraints, dof_map, variable_number, elem, /*Dim=*/2); }
510 
511 template <>
513  DofMap & dof_map,
514  const unsigned int variable_number,
515  const Elem * elem)
516 { raviart_thomas_compute_constraints(constraints, dof_map, variable_number, elem, /*Dim=*/3); }
517 
518 template <>
520  DofMap &,
521  const unsigned int,
522  const Elem *)
523 { RAVIART_LOW_D_ERROR_MESSAGE }
524 
525 template <>
527  DofMap &,
528  const unsigned int,
529  const Elem *)
530 { RAVIART_LOW_D_ERROR_MESSAGE }
531 
532 template <>
534  DofMap & dof_map,
535  const unsigned int variable_number,
536  const Elem * elem)
537 { raviart_thomas_compute_constraints(constraints, dof_map, variable_number, elem, /*Dim=*/2); }
538 
539 template <>
541  DofMap & dof_map,
542  const unsigned int variable_number,
543  const Elem * elem)
544 { raviart_thomas_compute_constraints(constraints, dof_map, variable_number, elem, /*Dim=*/3); }
545 #endif // LIBMESH_ENABLE_AMR
546 
547 // Specialize useless shape function methods
548 template <>
549 RealGradient FE<0,RAVIART_THOMAS>::shape(const ElemType, const Order,const unsigned int,const Point &)
550 { RAVIART_LOW_D_ERROR_MESSAGE }
551 template <>
552 RealGradient FE<0,RAVIART_THOMAS>::shape(const Elem *,const Order,const unsigned int,const Point &,const bool)
553 { RAVIART_LOW_D_ERROR_MESSAGE }
554 template <>
555 RealGradient FE<0,L2_RAVIART_THOMAS>::shape(const ElemType, const Order,const unsigned int,const Point &)
556 { RAVIART_LOW_D_ERROR_MESSAGE }
557 template <>
558 RealGradient FE<0,L2_RAVIART_THOMAS>::shape(const Elem *,const Order,const unsigned int,const Point &,const bool)
559 { RAVIART_LOW_D_ERROR_MESSAGE }
560 template <>
562  const unsigned int,const Point &)
563 { RAVIART_LOW_D_ERROR_MESSAGE }
564 template <>
565 RealGradient FE<0,RAVIART_THOMAS>::shape_deriv(const Elem *,const Order,const unsigned int,
566  const unsigned int,const Point &,const bool)
567 { RAVIART_LOW_D_ERROR_MESSAGE }
568 template <>
570  const unsigned int,const Point &)
571 { RAVIART_LOW_D_ERROR_MESSAGE }
572 template <>
574  const unsigned int,const Point &,const bool)
575 { RAVIART_LOW_D_ERROR_MESSAGE }
576 
577 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
578 template <>
580  const unsigned int,const Point &)
581 { RAVIART_LOW_D_ERROR_MESSAGE }
582 template <>
584  const unsigned int,const Point &,const bool)
585 { RAVIART_LOW_D_ERROR_MESSAGE }
586 template <>
588  const unsigned int,const Point &)
589 { RAVIART_LOW_D_ERROR_MESSAGE }
590 template <>
592  const unsigned int,const Point &,const bool)
593 { RAVIART_LOW_D_ERROR_MESSAGE }
594 
595 #endif
596 
597 template <>
598 RealGradient FE<1,RAVIART_THOMAS>::shape(const ElemType, const Order,const unsigned int,const Point &)
599 { RAVIART_LOW_D_ERROR_MESSAGE }
600 template <>
601 RealGradient FE<1,RAVIART_THOMAS>::shape(const Elem *,const Order,const unsigned int,const Point &,const bool)
602 { RAVIART_LOW_D_ERROR_MESSAGE }
603 template <>
604 RealGradient FE<1,L2_RAVIART_THOMAS>::shape(const ElemType, const Order,const unsigned int,const Point &)
605 { RAVIART_LOW_D_ERROR_MESSAGE }
606 template <>
607 RealGradient FE<1,L2_RAVIART_THOMAS>::shape(const Elem *,const Order,const unsigned int,const Point &,const bool)
608 { RAVIART_LOW_D_ERROR_MESSAGE }
609 template <>
611  const unsigned int,const Point &)
612 { RAVIART_LOW_D_ERROR_MESSAGE }
613 template <>
614 RealGradient FE<1,RAVIART_THOMAS>::shape_deriv(const Elem *,const Order,const unsigned int,
615  const unsigned int,const Point &,const bool)
616 { RAVIART_LOW_D_ERROR_MESSAGE }
617 template <>
619  const unsigned int,const Point &)
620 { RAVIART_LOW_D_ERROR_MESSAGE }
621 template <>
623  const unsigned int,const Point &,const bool)
624 { RAVIART_LOW_D_ERROR_MESSAGE }
625 
626 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
627 template <>
629  const unsigned int,const Point &)
630 { RAVIART_LOW_D_ERROR_MESSAGE }
631 template <>
633  const unsigned int,const Point &,const bool)
634 { RAVIART_LOW_D_ERROR_MESSAGE }
635 template <>
637  const unsigned int,const Point &)
638 { RAVIART_LOW_D_ERROR_MESSAGE }
639 template <>
641  const unsigned int,const Point &,const bool)
642 { RAVIART_LOW_D_ERROR_MESSAGE }
643 
644 #endif
645 
646 } // namespace libMesh
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:196
static unsigned int n_dofs(const ElemType t, const Order o)
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)
unsigned int dim
This is the base class from which all geometric element types are derived.
Definition: elem.h:94
static unsigned int n_dofs_at_node(const ElemType t, const Order o, const unsigned int n)
static OutputShape shape_deriv(const ElemType t, const Order o, const unsigned int i, const unsigned int j, const Point &p)
unsigned int p_level() const
Definition: elem.h:3108
virtual bool shapes_need_reinit() const override
The libMesh namespace provides an interface to certain functionality in the library.
virtual bool is_hierarchic() const override
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:179
LIBMESH_DEFAULT_VECTORIZED_FE(template<>Real FE< 0, BERNSTEIN)
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
virtual unsigned int n_nodes() const =0
static void get_refspace_nodes(const ElemType t, std::vector< Point > &nodes)
Definition: fe_abstract.C:400
static std::unique_ptr< FEGenericBase > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type.
libmesh_assert(ctx)
static unsigned int n_dofs_per_elem(const ElemType t, const Order o)
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...
static void nodal_soln(const Elem *elem, const Order o, const std::vector< Number > &elem_soln, std::vector< Number > &nodal_soln, bool add_p_level=true, const unsigned vdim=1)
Build the nodal soln from the element soln.
FEContinuity
defines an enum for finite element types to libmesh_assert a certain level (or type? Hcurl?) of continuity.
virtual ElemType type() const =0
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:39
The constraint matrix storage format.
Definition: dof_map.h:108
static OutputShape shape_second_deriv(const ElemType t, const Order o, const unsigned int i, const unsigned int j, const Point &p)