libMesh
fe.h
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 #ifndef LIBMESH_FE_H
21 #define LIBMESH_FE_H
22 
23 // Local includes
24 #include "libmesh/fe_base.h"
25 #include "libmesh/int_range.h"
26 #include "libmesh/libmesh.h"
27 
28 // C++ includes
29 #include <cstddef>
30 
31 namespace libMesh
32 {
33 
34 // forward declarations
35 class DofConstraints;
36 class DofMap;
37 class QGauss;
38 
39 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
40 
41 template <unsigned int friend_Dim, FEFamily friend_T_radial, InfMapType friend_T_map>
42 class InfFE;
43 
44 #endif
45 
46 
50 template <FEFamily T>
52 {
53  typedef Real type;
54 };
55 
56 
60 template<>
62 {
64 };
65 
66 template<>
68 {
70 };
71 
72 template<>
74 {
76 };
77 
78 template<>
80 {
82 };
83 
84 template<>
86 {
88 };
89 
90 template<>
92 {
94 };
95 
96 template<>
98 {
100 };
101 
102 template<>
104 {
106 };
107 
108 
126 template <unsigned int Dim, FEFamily T>
127 class FE : public FEGenericBase<typename FEOutputType<T>::type>
128 {
129 public:
130 
134  explicit
135  FE(const FEType & fet);
136 
137  typedef typename
140 
149  static OutputShape shape(const ElemType t,
150  const Order o,
151  const unsigned int i,
152  const Point & p);
153 
164  static OutputShape shape(const Elem * elem,
165  const Order o,
166  const unsigned int i,
167  const Point & p,
168  const bool add_p_level = true);
169 
180  static OutputShape shape(const FEType fet,
181  const Elem * elem,
182  const unsigned int i,
183  const Point & p,
184  const bool add_p_level = true);
185 
196  static void shapes(const Elem * elem,
197  const Order o,
198  const unsigned int i,
199  const std::vector<Point> & p,
200  std::vector<OutputShape> & v,
201  const bool add_p_level = true);
202 
203 
214  static void all_shapes(const Elem * elem,
215  const Order o,
216  const std::vector<Point> & p,
217  std::vector<std::vector<OutputShape> > & v,
218  const bool add_p_level = true);
219 
227  static OutputShape shape_deriv(const ElemType t,
228  const Order o,
229  const unsigned int i,
230  const unsigned int j,
231  const Point & p);
232 
241  static OutputShape shape_deriv(const Elem * elem,
242  const Order o,
243  const unsigned int i,
244  const unsigned int j,
245  const Point & p,
246  const bool add_p_level = true);
247 
257  static OutputShape shape_deriv(const FEType fet,
258  const Elem * elem,
259  const unsigned int i,
260  const unsigned int j,
261  const Point & p,
262  const bool add_p_level = true);
263 
274  static void shape_derivs(const Elem * elem,
275  const Order o,
276  const unsigned int i,
277  const unsigned int j,
278  const std::vector<Point> & p,
279  std::vector<OutputShape> & v,
280  const bool add_p_level = true);
281 
293  static void all_shape_derivs(const Elem * elem,
294  const Order o,
295  const std::vector<Point> & p,
296  std::vector<std::vector<OutputShape>> * comps[3],
297  const bool add_p_level = true);
298 
299 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
300 
320  static OutputShape shape_second_deriv(const ElemType t,
321  const Order o,
322  const unsigned int i,
323  const unsigned int j,
324  const Point & p);
325 
348  static OutputShape shape_second_deriv(const Elem * elem,
349  const Order o,
350  const unsigned int i,
351  const unsigned int j,
352  const Point & p,
353  const bool add_p_level = true);
354 
375  static OutputShape shape_second_deriv(const FEType fet,
376  const Elem * elem,
377  const unsigned int i,
378  const unsigned int j,
379  const Point & p,
380  const bool add_p_level = true);
381 
382 #endif //LIBMESH_ENABLE_SECOND_DERIVATIVES
383 
389  static void nodal_soln(const Elem * elem, const Order o,
390  const std::vector<Number> & elem_soln,
391  std::vector<Number> & nodal_soln,
392  bool add_p_level = true,
393  const unsigned vdim = 1);
394 
401  static void side_nodal_soln(const Elem * elem, const Order o,
402  const unsigned int side,
403  const std::vector<Number> & elem_soln,
404  std::vector<Number> & nodal_soln_on_side,
405  bool add_p_level = true,
406  const unsigned vdim = 1);
407 
412  virtual unsigned int n_shape_functions () const override;
413 
425  static unsigned int n_shape_functions (const ElemType t,
426  const Order o)
427  { return FE<Dim,T>::n_dofs (t,o); }
428 
440  static unsigned int n_dofs(const ElemType t,
441  const Order o);
442 
452  static unsigned int n_dofs(const Elem * e,
453  const Order o);
454 
465  static unsigned int n_dofs_at_node(const ElemType t,
466  const Order o,
467  const unsigned int n);
468 
475  static unsigned int n_dofs_at_node(const Elem & e,
476  const Order o,
477  const unsigned int n);
478 
488  static unsigned int n_dofs_per_elem(const ElemType t,
489  const Order o);
490 
497  static unsigned int n_dofs_per_elem(const Elem & e,
498  const Order o);
499 
503  virtual FEContinuity get_continuity() const override;
504 
509  virtual bool is_hierarchic() const override;
510 
517  static void dofs_on_side(const Elem * const elem,
518  const Order o,
519  unsigned int s,
520  std::vector<unsigned int> & di,
521  bool add_p_level=true);
528  static void dofs_on_edge(const Elem * const elem,
529  const Order o,
530  unsigned int e,
531  std::vector<unsigned int> & di,
532  bool add_p_level=true);
533 
534  static Point inverse_map (const Elem * elem,
535  const Point & p,
536  const Real tolerance = TOLERANCE,
537  const bool secure = true)
538  {
539  // libmesh_deprecated(); // soon
540  return FEMap::inverse_map(Dim, elem, p, tolerance, secure, secure);
541  }
542 
543  static void inverse_map (const Elem * elem,
544  const std::vector<Point> & physical_points,
545  std::vector<Point> & reference_points,
546  const Real tolerance = TOLERANCE,
547  const bool secure = true)
548  {
549  // libmesh_deprecated(); // soon
550  FEMap::inverse_map(Dim, elem, physical_points, reference_points,
551  tolerance, secure, secure);
552  }
553 
564  virtual void reinit (const Elem * elem,
565  const std::vector<Point> * const pts = nullptr,
566  const std::vector<Real> * const weights = nullptr) override;
567 
572  virtual void reinit_dual_shape_coeffs (const Elem * elem,
573  const std::vector<Point> & pts,
574  const std::vector<Real> & JxW) override;
575 
580  virtual void reinit_default_dual_shape_coeffs (const Elem * elem) override;
581 
591  virtual void reinit (const Elem * elem,
592  const unsigned int side,
593  const Real tolerance = TOLERANCE,
594  const std::vector<Point> * const pts = nullptr,
595  const std::vector<Real> * const weights = nullptr) override;
596 
606  virtual void edge_reinit (const Elem * elem,
607  const unsigned int edge,
608  const Real tolerance = TOLERANCE,
609  const std::vector<Point> * const pts = nullptr,
610  const std::vector<Real> * const weights = nullptr) override;
611 
616  virtual void side_map (const Elem * elem,
617  const Elem * side,
618  const unsigned int s,
619  const std::vector<Point> & reference_side_points,
620  std::vector<Point> & reference_points) override;
621 
626  virtual void edge_map (const Elem * elem,
627  const Elem * edge,
628  const unsigned int e,
629  const std::vector<Point> & reference_edge_points,
630  std::vector<Point> & reference_points);
631 
637  virtual void attach_quadrature_rule (QBase * q) override;
638 
639 #ifdef LIBMESH_ENABLE_AMR
640 
646  static void compute_constraints (DofConstraints & constraints,
647  DofMap & dof_map,
648  const unsigned int variable_number,
649  const Elem * elem);
650 #endif // #ifdef LIBMESH_ENABLE_AMR
651 
658  virtual bool shapes_need_reinit() const override;
659 
660  static Point map (const Elem * elem,
661  const Point & reference_point)
662  {
663  // libmesh_deprecated(); // soon
664  return FEMap::map(Dim, elem, reference_point);
665  }
666 
667  static Point map_xi (const Elem * elem,
668  const Point & reference_point)
669  {
670  // libmesh_deprecated(); // soon
671  return FEMap::map_deriv(Dim, elem, 0, reference_point);
672  }
673 
674  static Point map_eta (const Elem * elem,
675  const Point & reference_point)
676  {
677  // libmesh_deprecated(); // soon
678  return FEMap::map_deriv(Dim, elem, 1, reference_point);
679  }
680 
681  static Point map_zeta (const Elem * elem,
682  const Point & reference_point)
683  {
684  // libmesh_deprecated(); // soon
685  return FEMap::map_deriv(Dim, elem, 2, reference_point);
686  }
687 
688 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
689 
693  template <unsigned int friend_Dim, FEFamily friend_T_radial, InfMapType friend_T_map>
694  friend class InfFE;
695 #endif
696 
697 protected:
698 
706  virtual void init_shape_functions(const std::vector<Point> & qp,
707  const Elem * e);
708 
712  static void default_all_shape_derivs (const Elem * elem,
713  const Order o,
714  const std::vector<Point> & p,
715  std::vector<std::vector<OutputShape>> * comps[3],
716  const bool add_p_level = true);
717 
718 
722  void init_dual_shape_functions(unsigned int n_shapes, unsigned int n_qp);
723 
724 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
725 
730  virtual void init_base_shape_functions(const std::vector<Point> & qp,
731  const Elem * e) override;
732 
733 #endif
734 
738  static void default_shapes (const Elem * elem,
739  const Order o,
740  const unsigned int i,
741  const std::vector<Point> & p,
742  std::vector<OutputShape> & v,
743  const bool add_p_level = true)
744  {
745  libmesh_assert_equal_to(p.size(), v.size());
746  for (auto vi : index_range(v))
747  v[vi] = FE<Dim,T>::shape (elem, o, i, p[vi], add_p_level);
748  }
749 
753  static void default_all_shapes (const Elem * elem,
754  const Order o,
755  const std::vector<Point> & p,
756  std::vector<std::vector<OutputShape>> & v,
757  const bool add_p_level = true)
758  {
759  for (auto i : index_range(v))
760  {
761  libmesh_assert_equal_to ( p.size(), v[i].size() );
762  FE<Dim,T>::shapes (elem, o, i, p, v[i], add_p_level);
763  }
764  }
765 
769  static void default_shape_derivs (const Elem * elem,
770  const Order o,
771  const unsigned int i,
772  const unsigned int j,
773  const std::vector<Point> & p,
774  std::vector<OutputShape> & v,
775  const bool add_p_level = true)
776  {
777  libmesh_assert_equal_to(p.size(), v.size());
778  for (auto vi : index_range(v))
779  v[vi] = FE<Dim,T>::shape_deriv (elem, o, i, j, p[vi], add_p_level);
780  }
781 
785  static void default_side_nodal_soln(const Elem * elem, const Order o,
786  const unsigned int side,
787  const std::vector<Number> & elem_soln,
788  std::vector<Number> & nodal_soln_on_side,
789  bool add_p_level = true,
790  const unsigned vdim = 1);
791 
796  std::vector<Point> cached_nodes;
797  std::vector<bool> cached_edges, cached_faces;
798 
803  void cache(const Elem * elem);
804 
809  bool matches_cache(const Elem * elem);
810 
815 
817 };
818 
819 
820 
828 template <unsigned int Dim>
829 class FEClough : public FE<Dim,CLOUGH>
830 {
831 public:
832 
837  explicit
838  FEClough(const FEType & fet) :
839  FE<Dim,CLOUGH> (fet)
840  {}
841 };
842 
843 
844 
852 template <unsigned int Dim>
853 class FEHermite : public FE<Dim,HERMITE>
854 {
855 public:
856 
861  explicit
862  FEHermite(const FEType & fet) :
863  FE<Dim,HERMITE> (fet)
864  {}
865 
866 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
867 
870  static Real hermite_raw_shape_second_deriv(const unsigned int basis_num,
871  const Real xi);
872 #endif
873  static Real hermite_raw_shape_deriv(const unsigned int basis_num,
874  const Real xi);
875  static Real hermite_raw_shape(const unsigned int basis_num,
876  const Real xi);
877 };
878 
879 
880 
887 template <>
888 Real FE<2,SUBDIVISION>::shape(const Elem * elem,
889  const Order order,
890  const unsigned int i,
891  const Point & p,
892  const bool add_p_level);
893 
894 template <>
895 Real FE<2,SUBDIVISION>::shape_deriv(const Elem * elem,
896  const Order order,
897  const unsigned int i,
898  const unsigned int j,
899  const Point & p,
900  const bool add_p_level);
901 
902 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
903 template <>
905  const Order order,
906  const unsigned int i,
907  const unsigned int j,
908  const Point & p,
909  const bool add_p_level);
910 
911 #endif
912 
913 class FESubdivision : public FE<2,SUBDIVISION>
914 {
915 public:
916 
922  FESubdivision(const FEType & fet);
923 
934  virtual void reinit (const Elem * elem,
935  const std::vector<Point> * const pts = nullptr,
936  const std::vector<Real> * const weights = nullptr) override;
937 
942  virtual void reinit (const Elem *,
943  const unsigned int,
944  const Real = TOLERANCE,
945  const std::vector<Point> * const = nullptr,
946  const std::vector<Real> * const = nullptr) override
947  { libmesh_not_implemented(); }
948 
954  virtual void attach_quadrature_rule (QBase * q) override;
955 
963  virtual void init_shape_functions(const std::vector<Point> & qp,
964  const Elem * elem) override;
965 
972  static Real regular_shape(const unsigned int i,
973  const Real v,
974  const Real w);
975 
982  static Real regular_shape_deriv(const unsigned int i,
983  const unsigned int j,
984  const Real v,
985  const Real w);
986 
987 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
988 
994  static Real regular_shape_second_deriv(const unsigned int i,
995  const unsigned int j,
996  const Real v,
997  const Real w);
998 
999 
1000 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVE
1001 
1010  static void loop_subdivision_mask(std::vector<Real> & weights,
1011  const unsigned int valence);
1012 
1013 
1019  unsigned int valence);
1020 };
1021 
1022 
1023 
1031 template <unsigned int Dim>
1032 class FEHierarchic : public FE<Dim,HIERARCHIC>
1033 {
1034 public:
1035 
1040  explicit
1041  FEHierarchic(const FEType & fet) :
1042  FE<Dim,HIERARCHIC> (fet)
1043  {}
1044 };
1045 
1046 
1047 
1055 template <unsigned int Dim>
1056 class FEL2Hierarchic : public FE<Dim,L2_HIERARCHIC>
1057 {
1058 public:
1059 
1064  explicit
1065  FEL2Hierarchic(const FEType & fet) :
1066  FE<Dim,L2_HIERARCHIC> (fet)
1067  {}
1068 };
1069 
1070 
1071 
1079 template <unsigned int Dim>
1080 class FELagrange : public FE<Dim,LAGRANGE>
1081 {
1082 public:
1083 
1088  explicit
1089  FELagrange(const FEType & fet) :
1090  FE<Dim,LAGRANGE> (fet)
1091  {}
1092 };
1093 
1094 
1098 template <unsigned int Dim>
1099 class FEL2Lagrange : public FE<Dim,L2_LAGRANGE>
1100 {
1101 public:
1102 
1107  explicit
1108  FEL2Lagrange(const FEType & fet) :
1109  FE<Dim,L2_LAGRANGE> (fet)
1110  {}
1111 };
1112 
1113 
1121 template <unsigned int Dim>
1122 class FEMonomial : public FE<Dim,MONOMIAL>
1123 {
1124 public:
1125 
1130  explicit
1131  FEMonomial(const FEType & fet) :
1132  FE<Dim,MONOMIAL> (fet)
1133  {}
1134 };
1135 
1136 
1140 template <unsigned int Dim>
1141 class FEScalar : public FE<Dim,SCALAR>
1142 {
1143 public:
1144 
1151  explicit
1152  FEScalar(const FEType & fet) :
1153  FE<Dim,SCALAR> (fet)
1154  {}
1155 };
1156 
1157 
1166 template <unsigned int Dim>
1167 class FEXYZ : public FE<Dim,XYZ>
1168 {
1169 public:
1170 
1175  explicit
1176  FEXYZ(const FEType & fet) :
1177  FE<Dim,XYZ> (fet)
1178  {}
1179 
1185  virtual void reinit (const Elem * elem,
1186  const std::vector<Point> * const pts = nullptr,
1187  const std::vector<Real> * const weights = nullptr) override
1188  { FE<Dim,XYZ>::reinit (elem, pts, weights); }
1189 
1194  virtual void reinit (const Elem * elem,
1195  const unsigned int side,
1196  const Real tolerance = TOLERANCE,
1197  const std::vector<Point> * const pts = nullptr,
1198  const std::vector<Real> * const weights = nullptr) override;
1199 
1200 
1201 protected:
1202 
1210  virtual void init_shape_functions(const std::vector<Point> & qp,
1211  const Elem * e) override;
1212 
1223  virtual void compute_shape_functions(const Elem * elem, const std::vector<Point> & qp) override;
1224 
1228  void compute_face_values (const Elem * elem,
1229  const Elem * side,
1230  const std::vector<Real> & weights);
1231 };
1232 
1233 
1234 
1242 template <unsigned int Dim>
1243 class FELagrangeVec : public FE<Dim,LAGRANGE_VEC>
1244 {
1245 public:
1246 
1251  explicit
1252  FELagrangeVec(const FEType & fet) :
1253  FE<Dim,LAGRANGE_VEC> (fet)
1254  {}
1255 };
1256 
1257 
1265 template <unsigned int Dim>
1266 class FEL2LagrangeVec : public FE<Dim,L2_LAGRANGE_VEC>
1267 {
1268 public:
1269 
1274  explicit
1275  FEL2LagrangeVec(const FEType & fet) :
1276  FE<Dim,L2_LAGRANGE_VEC> (fet)
1277  {}
1278 };
1279 
1280 
1281 
1289 template <unsigned int Dim>
1290 class FEHierarchicVec : public FE<Dim,HIERARCHIC_VEC>
1291 {
1292 public:
1293 
1298  explicit
1299  FEHierarchicVec(const FEType & fet) :
1300  FE<Dim,HIERARCHIC_VEC> (fet)
1301  {}
1302 };
1303 
1304 
1312 template <unsigned int Dim>
1313 class FEL2HierarchicVec : public FE<Dim,L2_HIERARCHIC_VEC>
1314 {
1315 public:
1316 
1321  explicit
1322  FEL2HierarchicVec(const FEType & fet) :
1323  FE<Dim,L2_HIERARCHIC_VEC> (fet)
1324  {}
1325 };
1326 
1327 
1335 template <unsigned int Dim>
1336 class FENedelecOne : public FE<Dim,NEDELEC_ONE>
1337 {
1338 public:
1343  explicit
1344  FENedelecOne(const FEType & fet) :
1345  FE<Dim,NEDELEC_ONE> (fet)
1346  {}
1347 };
1348 
1356 template <unsigned int Dim>
1357 class FEMonomialVec : public FE<Dim,MONOMIAL_VEC>
1358 {
1359 public:
1360 
1365  explicit
1366  FEMonomialVec(const FEType & fet) :
1367  FE<Dim,MONOMIAL_VEC> (fet)
1368  {}
1369 };
1370 
1378 template <unsigned int Dim>
1379 class FERaviartThomas : public FE<Dim,RAVIART_THOMAS>
1380 {
1381 public:
1382 
1387  explicit
1388  FERaviartThomas(const FEType & fet) :
1389  FE<Dim,RAVIART_THOMAS> (fet)
1390  {}
1391 };
1392 
1401 template <unsigned int Dim>
1402 class FEL2RaviartThomas : public FE<Dim,L2_RAVIART_THOMAS>
1403 {
1404 public:
1405 
1410  explicit
1411  FEL2RaviartThomas(const FEType & fet) :
1412  FE<Dim,L2_RAVIART_THOMAS> (fet)
1413  {}
1414 };
1415 
1419 namespace FiniteElements
1420 {
1426 
1432 
1438 
1444 
1445 
1451 
1457 
1463 
1464 
1470 
1476 
1482 
1483 
1489 
1495 
1501 
1502 
1508 
1514 
1520 
1521 }
1522 
1527 template <typename OutputShape>
1528 OutputShape fe_fdm_deriv(const Elem * elem,
1529  const Order order,
1530  const unsigned int i,
1531  const unsigned int j,
1532  const Point & p,
1533  const bool add_p_level,
1534  OutputShape(*shape_func)
1535  (const Elem *, const Order,
1536  const unsigned int, const Point &,
1537  const bool));
1538 
1539 template <typename OutputShape>
1540 OutputShape fe_fdm_deriv(const ElemType type,
1541  const Order order,
1542  const unsigned int i,
1543  const unsigned int j,
1544  const Point & p,
1545  OutputShape(*shape_func)
1546  (const ElemType, const Order,
1547  const unsigned int, const Point &));
1548 
1549 template <typename OutputShape>
1550 OutputShape fe_fdm_deriv(const ElemType type,
1551  const Order order,
1552  const Elem * elem,
1553  const unsigned int i,
1554  const unsigned int j,
1555  const Point & p,
1556  OutputShape(*shape_func)
1557  (const ElemType type, const Order,
1558  const Elem *, const unsigned int,
1559  const Point &));
1560 
1561 
1562 template <typename OutputShape>
1563 OutputShape
1564 fe_fdm_second_deriv(const Elem * elem,
1565  const Order order,
1566  const unsigned int i,
1567  const unsigned int j,
1568  const Point & p,
1569  const bool add_p_level,
1570  OutputShape(*deriv_func)
1571  (const Elem *, const Order,
1572  const unsigned int, const unsigned int,
1573  const Point &, const bool));
1574 
1575 template <typename OutputShape>
1576 OutputShape fe_fdm_second_deriv(const ElemType type,
1577  const Order order,
1578  const unsigned int i,
1579  const unsigned int j,
1580  const Point & p,
1581  OutputShape(*deriv_func)
1582  (const ElemType, const Order,
1583  const unsigned int,
1584  const unsigned int,
1585  const Point &));
1586 
1587 template <typename OutputShape>
1588 OutputShape fe_fdm_second_deriv(const ElemType type,
1589  const Order order,
1590  const Elem * elem,
1591  const unsigned int i,
1592  const unsigned int j,
1593  const Point & p,
1594  OutputShape(*deriv_func)
1595  (const ElemType, const Order,
1596  const Elem *,
1597  const unsigned int,
1598  const unsigned int,
1599  const Point &));
1600 
1604 void lagrange_nodal_soln(const Elem * elem,
1605  const Order order,
1606  const std::vector<Number> & elem_soln,
1607  std::vector<Number> & nodal_soln,
1608  bool add_p_level = true);
1609 
1613 unsigned int monomial_n_dofs(const ElemType t, const Order o);
1614 
1615 unsigned int monomial_n_dofs(const Elem * e, const Order o);
1616 
1620 // shapes[i][j] is shape function phi_i at point p[j]
1621 void rational_fe_weighted_shapes(const Elem * elem,
1622  const FEType underlying_fe_type,
1623  std::vector<std::vector<Real>> & shapes,
1624  const std::vector<Point> & p,
1625  const bool add_p_level);
1626 
1627 // shapes[i][q] is shape function phi_i at point p[q]
1628 // derivs[j][i][q] is dphi_i/dxi_j at p[q]
1629 void rational_fe_weighted_shapes_derivs(const Elem * elem,
1630  const FEType fe_type,
1631  std::vector<std::vector<Real>> & shapes,
1632  std::vector<std::vector<std::vector<Real>>> & derivs,
1633  const std::vector<Point> & p,
1634  const bool add_p_level);
1635 
1636 Real rational_fe_shape(const Elem & elem,
1637  const FEType underlying_fe_type,
1638  const unsigned int i,
1639  const Point & p,
1640  const bool add_p_level);
1641 
1642 Real rational_fe_shape_deriv(const Elem & elem,
1643  const FEType underlying_fe_type,
1644  const unsigned int i,
1645  const unsigned int j,
1646  const Point & p,
1647  const bool add_p_level);
1648 
1650  const FEType underlying_fe_type,
1651  const unsigned int i,
1652  const unsigned int j,
1653  const Point & p,
1654  const bool add_p_level);
1655 
1656 void rational_all_shapes (const Elem & elem,
1657  const FEType underlying_fe_type,
1658  const std::vector<Point> & p,
1659  std::vector<std::vector<Real>> & v,
1660  const bool add_p_level);
1661 
1662 template <typename OutputShape>
1663 void rational_all_shape_derivs (const Elem & elem,
1664  const FEType underlying_fe_type,
1665  const std::vector<Point> & p,
1666  std::vector<std::vector<OutputShape>> * comps[3],
1667  const bool add_p_level);
1668 
1669 } // namespace libMesh
1670 
1671 
1672 // Full specialization of all n_dofs type functions, for every
1673 // dimension, with both original ElemType and new Elem signatures
1674 #define LIBMESH_DEFAULT_NDOFS(MyType) \
1675 template <> unsigned int FE<0,MyType>::n_dofs(const ElemType t, const Order o) { return MyType##_n_dofs(t, o); } \
1676 template <> unsigned int FE<1,MyType>::n_dofs(const ElemType t, const Order o) { return MyType##_n_dofs(t, o); } \
1677 template <> unsigned int FE<2,MyType>::n_dofs(const ElemType t, const Order o) { return MyType##_n_dofs(t, o); } \
1678 template <> unsigned int FE<3,MyType>::n_dofs(const ElemType t, const Order o) { return MyType##_n_dofs(t, o); } \
1679  \
1680 template <> unsigned int FE<0,MyType>::n_dofs(const Elem * e, const Order o) { return MyType##_n_dofs(e, o); } \
1681 template <> unsigned int FE<1,MyType>::n_dofs(const Elem * e, const Order o) { return MyType##_n_dofs(e, o); } \
1682 template <> unsigned int FE<2,MyType>::n_dofs(const Elem * e, const Order o) { return MyType##_n_dofs(e, o); } \
1683 template <> unsigned int FE<3,MyType>::n_dofs(const Elem * e, const Order o) { return MyType##_n_dofs(e, o); } \
1684  \
1685 template <> unsigned int FE<0,MyType>::n_dofs_at_node(const ElemType t, const Order o, const unsigned int n) { return MyType##_n_dofs_at_node(t, o, n); } \
1686 template <> unsigned int FE<1,MyType>::n_dofs_at_node(const ElemType t, const Order o, const unsigned int n) { return MyType##_n_dofs_at_node(t, o, n); } \
1687 template <> unsigned int FE<2,MyType>::n_dofs_at_node(const ElemType t, const Order o, const unsigned int n) { return MyType##_n_dofs_at_node(t, o, n); } \
1688 template <> unsigned int FE<3,MyType>::n_dofs_at_node(const ElemType t, const Order o, const unsigned int n) { return MyType##_n_dofs_at_node(t, o, n); } \
1689  \
1690 template <> unsigned int FE<0,MyType>::n_dofs_at_node(const Elem & e, const Order o, const unsigned int n) { return MyType##_n_dofs_at_node(e, o, n); } \
1691 template <> unsigned int FE<1,MyType>::n_dofs_at_node(const Elem & e, const Order o, const unsigned int n) { return MyType##_n_dofs_at_node(e, o, n); } \
1692 template <> unsigned int FE<2,MyType>::n_dofs_at_node(const Elem & e, const Order o, const unsigned int n) { return MyType##_n_dofs_at_node(e, o, n); } \
1693 template <> unsigned int FE<3,MyType>::n_dofs_at_node(const Elem & e, const Order o, const unsigned int n) { return MyType##_n_dofs_at_node(e, o, n); } \
1694  \
1695 template <> unsigned int FE<0,MyType>::n_dofs_per_elem(const ElemType t, const Order o) { return MyType##_n_dofs_per_elem(t, o); } \
1696 template <> unsigned int FE<1,MyType>::n_dofs_per_elem(const ElemType t, const Order o) { return MyType##_n_dofs_per_elem(t, o); } \
1697 template <> unsigned int FE<2,MyType>::n_dofs_per_elem(const ElemType t, const Order o) { return MyType##_n_dofs_per_elem(t, o); } \
1698 template <> unsigned int FE<3,MyType>::n_dofs_per_elem(const ElemType t, const Order o) { return MyType##_n_dofs_per_elem(t, o); } \
1699  \
1700 template <> unsigned int FE<0,MyType>::n_dofs_per_elem(const Elem & e, const Order o) { return MyType##_n_dofs_per_elem(e, o); } \
1701 template <> unsigned int FE<1,MyType>::n_dofs_per_elem(const Elem & e, const Order o) { return MyType##_n_dofs_per_elem(e, o); } \
1702 template <> unsigned int FE<2,MyType>::n_dofs_per_elem(const Elem & e, const Order o) { return MyType##_n_dofs_per_elem(e, o); } \
1703 template <> unsigned int FE<3,MyType>::n_dofs_per_elem(const Elem & e, const Order o) { return MyType##_n_dofs_per_elem(e, o); }
1704 
1705 
1706 #define LIBMESH_DEFAULT_VEC_NDOFS(MyType) \
1707 template <> unsigned int FE<0,MyType##_VEC>::n_dofs(const ElemType t, const Order o) { return FE<0,MyType>::n_dofs(t, o); } \
1708 template <> unsigned int FE<1,MyType##_VEC>::n_dofs(const ElemType t, const Order o) { return FE<1,MyType>::n_dofs(t, o); } \
1709 template <> unsigned int FE<2,MyType##_VEC>::n_dofs(const ElemType t, const Order o) { return 2*FE<2,MyType>::n_dofs(t, o); } \
1710 template <> unsigned int FE<3,MyType##_VEC>::n_dofs(const ElemType t, const Order o) { return 3*FE<3,MyType>::n_dofs(t, o); } \
1711  \
1712 template <> unsigned int FE<0,MyType##_VEC>::n_dofs(const Elem * e, const Order o) { return FE<0,MyType>::n_dofs(e, o); } \
1713 template <> unsigned int FE<1,MyType##_VEC>::n_dofs(const Elem * e, const Order o) { return FE<1,MyType>::n_dofs(e, o); } \
1714 template <> unsigned int FE<2,MyType##_VEC>::n_dofs(const Elem * e, const Order o) { return 2*FE<2,MyType>::n_dofs(e, o); } \
1715 template <> unsigned int FE<3,MyType##_VEC>::n_dofs(const Elem * e, const Order o) { return 3*FE<3,MyType>::n_dofs(e, o); } \
1716  \
1717 template <> unsigned int FE<0,MyType##_VEC>::n_dofs_at_node(const ElemType t, const Order o, const unsigned int n) { return FE<0,MyType>::n_dofs_at_node(t, o, n); } \
1718 template <> unsigned int FE<1,MyType##_VEC>::n_dofs_at_node(const ElemType t, const Order o, const unsigned int n) { return FE<1,MyType>::n_dofs_at_node(t, o, n); } \
1719 template <> unsigned int FE<2,MyType##_VEC>::n_dofs_at_node(const ElemType t, const Order o, const unsigned int n) { return 2*FE<2,MyType>::n_dofs_at_node(t, o, n); } \
1720 template <> unsigned int FE<3,MyType##_VEC>::n_dofs_at_node(const ElemType t, const Order o, const unsigned int n) { return 3*FE<3,MyType>::n_dofs_at_node(t, o, n); } \
1721  \
1722 template <> unsigned int FE<0,MyType##_VEC>::n_dofs_at_node(const Elem & e, const Order o, const unsigned int n) { return FE<0,MyType>::n_dofs_at_node(e.type(), o, n); } \
1723 template <> unsigned int FE<1,MyType##_VEC>::n_dofs_at_node(const Elem & e, const Order o, const unsigned int n) { return FE<1,MyType>::n_dofs_at_node(e.type(), o, n); } \
1724 template <> unsigned int FE<2,MyType##_VEC>::n_dofs_at_node(const Elem & e, const Order o, const unsigned int n) { return 2*FE<2,MyType>::n_dofs_at_node(e.type(), o, n); } \
1725 template <> unsigned int FE<3,MyType##_VEC>::n_dofs_at_node(const Elem & e, const Order o, const unsigned int n) { return 3*FE<3,MyType>::n_dofs_at_node(e.type(), o, n); } \
1726  \
1727 template <> unsigned int FE<0,MyType##_VEC>::n_dofs_per_elem(const ElemType t, const Order o) { return FE<0,MyType>::n_dofs_per_elem(t, o); } \
1728 template <> unsigned int FE<1,MyType##_VEC>::n_dofs_per_elem(const ElemType t, const Order o) { return FE<1,MyType>::n_dofs_per_elem(t, o); } \
1729 template <> unsigned int FE<2,MyType##_VEC>::n_dofs_per_elem(const ElemType t, const Order o) { return 2*FE<2,MyType>::n_dofs_per_elem(t, o); } \
1730 template <> unsigned int FE<3,MyType##_VEC>::n_dofs_per_elem(const ElemType t, const Order o) { return 3*FE<3,MyType>::n_dofs_per_elem(t, o); } \
1731  \
1732 template <> unsigned int FE<0,MyType##_VEC>::n_dofs_per_elem(const Elem & e, const Order o) { return FE<0,MyType>::n_dofs_per_elem(e.type(), o); } \
1733 template <> unsigned int FE<1,MyType##_VEC>::n_dofs_per_elem(const Elem & e, const Order o) { return FE<1,MyType>::n_dofs_per_elem(e.type(), o); } \
1734 template <> unsigned int FE<2,MyType##_VEC>::n_dofs_per_elem(const Elem & e, const Order o) { return 2*FE<2,MyType>::n_dofs_per_elem(e.type(), o); } \
1735 template <> unsigned int FE<3,MyType##_VEC>::n_dofs_per_elem(const Elem & e, const Order o) { return 3*FE<3,MyType>::n_dofs_per_elem(e.type(), o); }
1736 
1737 
1738 #define LIBMESH_DEFAULT_VECTORIZED_FE(MyDim, MyType) \
1739 template<> \
1740 void FE<MyDim,MyType>::all_shapes \
1741  (const Elem * elem, \
1742  const Order o, \
1743  const std::vector<Point> & p, \
1744  std::vector<std::vector<OutputShape>> & v, \
1745  const bool add_p_level) \
1746 { \
1747  FE<MyDim,MyType>::default_all_shapes \
1748  (elem,o,p,v,add_p_level); \
1749 } \
1750  \
1751 template<> \
1752 void FE<MyDim,MyType>::shapes \
1753  (const Elem * elem, \
1754  const Order o, \
1755  const unsigned int i, \
1756  const std::vector<Point> & p, \
1757  std::vector<OutputShape> & v, \
1758  const bool add_p_level) \
1759 { \
1760  FE<MyDim,MyType>::default_shapes \
1761  (elem,o,i,p,v,add_p_level); \
1762 } \
1763  \
1764 template<> \
1765 void FE<MyDim,MyType>::shape_derivs \
1766  (const Elem * elem, \
1767  const Order o, \
1768  const unsigned int i, \
1769  const unsigned int j, \
1770  const std::vector<Point> & p, \
1771  std::vector<OutputShape> & v, \
1772  const bool add_p_level) \
1773 { \
1774  FE<MyDim,MyType>::default_shape_derivs \
1775  (elem,o,i,j,p,v,add_p_level); \
1776 } \
1777  \
1778 template<> \
1779 void FE<MyDim,MyType>::all_shape_derivs \
1780  (const Elem * elem, \
1781  const Order o, \
1782  const std::vector<Point> & p, \
1783  std::vector<std::vector<OutputShape>> * comps[3], \
1784  const bool add_p_level) \
1785 { \
1786  FE<MyDim,MyType>::default_all_shape_derivs \
1787  (elem,o,p,comps,add_p_level); \
1788 }
1789 
1790 
1791 #endif // LIBMESH_FE_H
FELagrangeVec objects are used for working with vector-valued finite elements.
Definition: fe.h:1243
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:196
XYZ finite elements.
Definition: fe.h:1167
FELagrange(const FEType &fet)
Constructor.
Definition: fe.h:1089
static Point map_xi(const Elem *elem, const Point &reference_point)
Definition: fe.h:667
FESubdivision(const FEType &fet)
Constructor.
FEMonomial(const FEType &fet)
Constructor.
Definition: fe.h:1131
static unsigned int n_dofs(const ElemType t, const Order o)
ElemType
Defines an enum for geometric element types.
static Real regular_shape(const unsigned int i, const Real v, const Real w)
Order
defines an enum for polynomial orders.
Definition: enum_order.h:40
static void dofs_on_edge(const Elem *const elem, const Order o, unsigned int e, std::vector< unsigned int > &di, bool add_p_level=true)
Fills the vector di with the local degree of freedom indices associated with edge e of element elem...
Definition: fe.C:126
A specific instantiation of the FEBase class.
Definition: fe.h:42
FENedelecOne objects are used for working with vector-valued Nedelec finite elements of the first kin...
Definition: fe.h:1336
virtual void side_map(const Elem *elem, const Elem *side, const unsigned int s, const std::vector< Point > &reference_side_points, std::vector< Point > &reference_points) override
Computes the reference space quadrature points on the side of an element based on the side quadrature...
Definition: fe_boundary.C:347
static OutputShape shape(const ElemType t, const Order o, const unsigned int i, const Point &p)
FE< 1, MONOMIAL > FEMonomial1D
Convenient definition for a 1D Monomial finite element.
Definition: fe.h:1507
FE< 2, L2_HIERARCHIC > FEL2Hierarchic2D
Convenient definition for a 2D Discontinuous Hierarchic finite element.
Definition: fe.h:1456
virtual void attach_quadrature_rule(QBase *q) override
Provides the class with the quadrature rule, which provides the locations (on a reference element) wh...
Definition: fe.C:87
FE< 3, L2_HIERARCHIC > FEL2Hierarchic3D
Convenient definition for a 3D Discontinuous Hierarchic finite element.
Definition: fe.h:1462
static void init_subdivision_matrix(DenseMatrix< Real > &A, unsigned int valence)
Builds the subdivision matrix A for the Loop scheme.
void compute_face_values(const Elem *elem, const Elem *side, const std::vector< Real > &weights)
Compute the map & shape functions for this face.
void rational_fe_weighted_shapes_derivs(const Elem *elem, const FEType fe_type, std::vector< std::vector< Real >> &shapes, std::vector< std::vector< std::vector< Real >>> &derivs, const std::vector< Point > &p, const bool add_p_level)
Definition: fe.C:1065
Lagrange finite elements.
Definition: fe.h:1080
Monomial finite elements.
Definition: fe.h:1122
void init_dual_shape_functions(unsigned int n_shapes, unsigned int n_qp)
Init dual_phi and potentially dual_dphi, dual_d2phi.
Definition: fe.C:411
FEClough(const FEType &fet)
Constructor.
Definition: fe.h:838
static void loop_subdivision_mask(std::vector< Real > &weights, const unsigned int valence)
Fills the vector weights with the weight coefficients of the Loop subdivision mask for evaluating the...
static constexpr Real TOLERANCE
std::vector< bool > cached_faces
Definition: fe.h:797
FE< 1, L2_LAGRANGE > FEL2Lagrange1D
Convenient definition for a 1D Discontinuous Lagrange finite element.
Definition: fe.h:1488
static Point inverse_map(const unsigned int dim, const Elem *elem, const Point &p, const Real tolerance=TOLERANCE, const bool secure=true, const bool extra_checks=true)
Definition: fe_map.C:1628
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...
static void default_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)
A default implementation for all_shape_derivs.
Definition: fe.C:736
FEHierarchicVec objects are used for working with vector-valued high-order piecewise-continuous finit...
Definition: fe.h:1313
FERaviartThomas objects are used for working with vector-valued Raviart-Thomas finite elements...
Definition: fe.h:1379
static Point map(const Elem *elem, const Point &reference_point)
Definition: fe.h:660
FERaviartThomas(const FEType &fet)
Constructor.
Definition: fe.h:1388
This is the base class from which all geometric element types are derived.
Definition: elem.h:94
FEGenericBase< typename FEOutputType< T >::type >::OutputShape OutputShape
Definition: fe.h:139
FE(const FEType &fet)
Constructor.
Definition: fe.C:62
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)
FE< 3, L2_LAGRANGE > FEL2Lagrange3D
Convenient definition for a 3D Discontinuous Lagrange finite element.
Definition: fe.h:1500
FEL2RaviartThomas(const FEType &fet)
Constructor.
Definition: fe.h:1411
virtual bool shapes_need_reinit() const override
FE< 1, HIERARCHIC > FEHierarchic1D
Convenient definition for a 1D Hierarchic finite element.
Definition: fe.h:1431
The libMesh namespace provides an interface to certain functionality in the library.
static Real hermite_raw_shape_second_deriv(const unsigned int basis_num, const Real xi)
1D hermite functions on unit interval
FEL2RaviartThomas objects are used for working with vector-valued discontinuous Raviart-Thomas finite...
Definition: fe.h:1402
void cache(const Elem *elem)
Repopulate the element cache with the node locations, edge and face orientations of the element elem...
Definition: fe.C:153
FEXYZ(const FEType &fet)
Constructor.
Definition: fe.h:1176
FE< 2, LAGRANGE > FELagrange2D
Convenient definition for a 2D Lagrange finite element.
Definition: fe.h:1475
FELagrangeVec(const FEType &fet)
Constructor.
Definition: fe.h:1252
static void default_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)
A default implementation for all_shapes.
Definition: fe.h:753
virtual void compute_shape_functions(const Elem *elem, const std::vector< Point > &qp) override
After having updated the jacobian and the transformation from local to global coordinates in FEAbstra...
Definition: fe_xyz.C:205
FEL2LagrangeVec(const FEType &fet)
Constructor.
Definition: fe.h:1275
static Real hermite_raw_shape(const unsigned int basis_num, const Real xi)
virtual void attach_quadrature_rule(QBase *q) override
Provides the class with the quadrature rule, which provides the locations (on a reference element) wh...
FE< 3, LAGRANGE > FELagrange3D
Convenient definition for a 3D Lagrange finite element.
Definition: fe.h:1481
FEClough< 2 > FEClough2D
Convenient definition for a 2D Clough-Tocher finite element.
Definition: fe.h:1425
FEL2Lagrange(const FEType &fet)
Constructor.
Definition: fe.h:1108
FEHierarchicVec objects are used for working with vector-valued high-order finite elements...
Definition: fe.h:1290
virtual bool is_hierarchic() const override
OutputShape fe_fdm_deriv(const Elem *elem, const Order order, const unsigned int i, const unsigned int j, const Point &p, const bool add_p_level, OutputShape(*shape_func)(const Elem *, const Order, const unsigned int, const Point &, const bool))
Helper functions for finite differenced derivatives in cases where analytical calculations haven&#39;t be...
Definition: fe.C:911
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:179
FEMonomialVec(const FEType &fet)
Constructor.
Definition: fe.h:1366
FEHierarchic(const FEType &fet)
Constructor.
Definition: fe.h:1041
virtual void reinit(const Elem *elem, const std::vector< Point > *const pts=nullptr, const std::vector< Real > *const weights=nullptr) override
This is at the core of this class.
Definition: fe.C:197
static Point map_eta(const Elem *elem, const Point &reference_point)
Definition: fe.h:674
virtual unsigned int n_shape_functions() const override
Definition: fe.C:75
A specific instantiation of the FEBase class.
Definition: fe.h:127
static void dofs_on_side(const Elem *const elem, const Order o, unsigned int s, std::vector< unsigned int > &di, bool add_p_level=true)
Fills the vector di with the local degree of freedom indices associated with side s of element elem...
Definition: fe.C:99
static Real regular_shape_second_deriv(const unsigned int i, const unsigned int j, const Real v, const Real w)
FEHermite(const FEType &fet)
Constructor.
Definition: fe.h:862
ElemType last_edge
Definition: fe.h:816
FE< 3, HIERARCHIC > FEHierarchic3D
Convenient definition for a 3D Hierarchic finite element.
Definition: fe.h:1443
Discontinuous Hierarchic finite elements.
Definition: fe.h:1056
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
unsigned int monomial_n_dofs(const ElemType t, const Order o)
Helper functions for Discontinuous-Pn type basis functions.
Definition: fe_monomial.C:37
FEL2Hierarchic(const FEType &fet)
Constructor.
Definition: fe.h:1065
FE< 2, MONOMIAL > FEMonomial2D
Convenient definition for a 2D Monomial finite element.
Definition: fe.h:1513
static Real regular_shape_deriv(const unsigned int i, const unsigned int j, const Real v, const Real w)
FEOutputType< T >::type OutputShape
Convenient typedefs for gradients of output, hessians of output, and potentially-complex-valued versi...
Definition: fe_base.h:119
OutputShape fe_fdm_second_deriv(const Elem *elem, const Order order, const unsigned int i, const unsigned int j, const Point &p, const bool add_p_level, OutputShape(*deriv_func)(const Elem *, const Order, const unsigned int, const unsigned int, const Point &, const bool))
Definition: fe.C:969
static void default_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)
A default implementation for shape_derivs.
Definition: fe.h:769
static unsigned int n_dofs_per_elem(const ElemType t, const Order o)
virtual FEContinuity get_continuity() const override
FENedelecOne(const FEType &fet)
Constructor.
Definition: fe.h:1344
void lagrange_nodal_soln(const Elem *elem, const Order order, const std::vector< Number > &elem_soln, std::vector< Number > &nodal_soln, bool add_p_level=true)
Helper functions for Lagrange-based basis functions.
Definition: fe_lagrange.C:43
FEHierarchicVec(const FEType &fet)
Constructor.
Definition: fe.h:1299
std::vector< bool > cached_edges
Definition: fe.h:797
static void side_nodal_soln(const Elem *elem, const Order o, const unsigned int side, const std::vector< Number > &elem_soln, std::vector< Number > &nodal_soln_on_side, bool add_p_level=true, const unsigned vdim=1)
Build the nodal soln on one side from the (full) element soln.
Discontinuous Lagrange finite elements.
Definition: fe.h:1099
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.
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
static unsigned int n_shape_functions(const ElemType t, const Order o)
Definition: fe.h:425
bool matches_cache(const Elem *elem)
Check if the node locations, edge and face orientations held in the element cache match those of elem...
Definition: fe.C:174
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual void edge_map(const Elem *elem, const Elem *edge, const unsigned int e, const std::vector< Point > &reference_edge_points, std::vector< Point > &reference_points)
Computes the reference space quadrature points on the side of an element based on the edge quadrature...
Definition: fe_boundary.C:403
FE< 2, HIERARCHIC > FEHierarchic2D
Convenient definition for a 2D Hierarchic finite element.
Definition: fe.h:1437
static Point map_zeta(const Elem *elem, const Point &reference_point)
Definition: fe.h:681
Hierarchic finite elements.
Definition: fe.h:1032
FEContinuity
defines an enum for finite element types to libmesh_assert a certain level (or type? Hcurl?) of continuity.
static Point map(const unsigned int dim, const Elem *elem, const Point &reference_point)
Definition: fe_map.C:2069
static void inverse_map(const Elem *elem, const std::vector< Point > &physical_points, std::vector< Point > &reference_points, const Real tolerance=TOLERANCE, const bool secure=true)
Definition: fe.h:543
FE< 3, MONOMIAL > FEMonomial3D
Convenient definition for a 3D Monomial finite element.
Definition: fe.h:1519
virtual void init_shape_functions(const std::vector< Point > &qp, const Elem *e)
Update the various member data fields phi, dphidxi, dphideta, dphidzeta, etc.
Definition: fe.C:441
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.
FEMonomialVec objects are used for working with vector-valued discontinuous finite elements...
Definition: fe.h:1357
static Point map_deriv(const unsigned int dim, const Elem *elem, const unsigned int j, const Point &reference_point)
Definition: fe_map.C:2103
static Real hermite_raw_shape_deriv(const unsigned int basis_num, const Real xi)
FEL2LagrangeVec objects are used for working with vector-valued finite elements.
Definition: fe.h:1266
FE< 2, L2_LAGRANGE > FEL2Lagrange2D
Convenient definition for a 2D Discontinuous Lagrange finite element.
Definition: fe.h:1494
FE< 1, LAGRANGE > FELagrange1D
Convenient definition for a 1D Lagrange finite element.
Definition: fe.h:1469
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.
The FEScalar class is used for working with SCALAR variables.
Definition: fe.h:1141
virtual void reinit_default_dual_shape_coeffs(const Elem *elem) override
This computes the default dual shape function coefficients.
Definition: fe.C:394
FEScalar(const FEType &fet)
Constructor.
Definition: fe.h:1152
Hermite finite elements.
Definition: fe.h:853
virtual void reinit_dual_shape_coeffs(const Elem *elem, const std::vector< Point > &pts, const std::vector< Real > &JxW) override
This re-computes the dual shape function coefficients.
Definition: fe.C:371
FE< 1, L2_HIERARCHIC > FEL2Hierarchic1D
Convenient definition for a 1D Discontinuous Hierarchic finite element.
Definition: fe.h:1450
FEL2HierarchicVec(const FEType &fet)
Constructor.
Definition: fe.h:1322
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
std::vector< Point > cached_nodes
Vectors holding the node locations, edge and face orientations of the last element we cached...
Definition: fe.h:796
virtual void init_shape_functions(const std::vector< Point > &qp, const Elem *elem) override
Update the various member data fields phi, dphidxi, dphideta, dphidzeta, etc.
Clough-Tocher finite elements.
Definition: fe.h:829
ElemType last_side
The last side and last edge we did a reinit on.
Definition: fe.h:814
static Point inverse_map(const Elem *elem, const Point &p, const Real tolerance=TOLERANCE, const bool secure=true)
Definition: fe.h:534
void rational_fe_weighted_shapes(const Elem *elem, const FEType underlying_fe_type, std::vector< std::vector< Real >> &shapes, const std::vector< Point > &p, const bool add_p_level)
Helper functions for rational basis functions.
Definition: fe.C:1027
virtual void edge_reinit(const Elem *elem, const unsigned int edge, const Real tolerance=TOLERANCE, const std::vector< Point > *const pts=nullptr, const std::vector< Real > *const weights=nullptr) override
Reinitializes all the physical element-dependent data based on the edge.
Definition: fe_boundary.C:243
virtual void init_shape_functions(const std::vector< Point > &qp, const Elem *e) override
Update the various member data fields phi, dphidxi, dphideta, dphidzeta, etc.
Definition: fe_xyz.C:98
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
virtual void reinit(const Elem *elem, const std::vector< Point > *const pts=nullptr, const std::vector< Real > *const weights=nullptr) override
This is at the core of this class.
virtual void reinit(const Elem *elem, const std::vector< Point > *const pts=nullptr, const std::vector< Real > *const weights=nullptr) override
Explicitly call base class method.
Definition: fe.h:1185
The QBase class provides the basic functionality from which various quadrature rules can be derived...
Definition: quadrature.h:61
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
This class forms the foundation from which generic finite elements may be derived.
static OutputShape shape_second_deriv(const ElemType t, const Order o, const unsigned int i, const unsigned int j, const Point &p)
Most finite element types in libMesh are scalar-valued.
Definition: fe.h:51
static void default_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)
A default implementation for shapes.
Definition: fe.h:738
static void default_side_nodal_soln(const Elem *elem, const Order o, const unsigned int side, const std::vector< Number > &elem_soln, std::vector< Number > &nodal_soln_on_side, bool add_p_level=true, const unsigned vdim=1)
A default implementation for side_nodal_soln.
Definition: fe.C:754
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
virtual void init_base_shape_functions(const std::vector< Point > &qp, const Elem *e) override
Initialize the data fields for the base of an an infinite element.
Definition: fe.C:778
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...
virtual void reinit(const Elem *, const unsigned int, const Real=TOLERANCE, const std::vector< Point > *const =nullptr, const std::vector< Real > *const =nullptr) override
This prevents some compilers being confused by partially overriding this virtual function.
Definition: fe.h:942