libMesh
inf_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_INF_FE_H
21 #define LIBMESH_INF_FE_H
22 
23 #include "libmesh/libmesh_config.h"
24 
25 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
26 
27 // Local includes
28 #include "libmesh/fe_base.h"
29 #include "libmesh/inf_fe_map.h"
30 
31 // C++ includes
32 #include <cstddef>
33 
34 namespace libMesh
35 {
36 
37 
38 // forward declarations
39 class Elem;
40 class FEComputeData;
41 
42 
54 {
55 private:
56 
61 
62 public:
63 
68  static Real decay (const unsigned int dim, const Real v);
69 
74  static Real decay_deriv (const unsigned int dim, const Real);
75 
82  static Real D (const Real v) { return (1.-v)*(1.-v)/4.; }
83 
84  static Real Dxr_sq (const Real /*v*/) { return 1.; }
85 
90  static Real D_deriv (const Real v) { return (v-1.)/2.; }
91 
92 
97  static Order mapping_order() { return FIRST; }
98 
113  static unsigned int n_dofs (const Order o_radial)
114  { return static_cast<unsigned int>(o_radial)+1; }
115 
125  static unsigned int n_dofs_at_node (const Order o_radial,
126  const unsigned int n_onion);
127 
136  static unsigned int n_dofs_per_elem (const Order o_radial)
137  { return static_cast<unsigned int>(o_radial)+1; }
138 
139 };
140 
141 
142 
152 {
153 private:
154 
159 
160 public:
161 
168  static std::unique_ptr<const Elem> build_elem (const Elem * inf_elem);
169 
175  static ElemType get_elem_type (const ElemType type);
176 
182  static unsigned int n_base_mapping_sf (const Elem & base_elem,
183  const Order base_mapping_order);
184 
185 };
186 
187 
188 
223 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map>
224 class InfFE : public FEBase
225 {
226 
227  /*
228  * Protect the nested class
229  */
230 protected:
231 
232 
233 public:
234 
235  // InfFE continued
236 
250  explicit
251  InfFE(const FEType & fet);
252  ~InfFE() = default;
253 
254  // The static public members for access from FEInterface etc
255 
268  static Real shape(const FEType & fet,
269  const ElemType t,
270  const unsigned int i,
271  const Point & p);
272 
285  static Real shape(const FEType & fet,
286  const Elem * elem,
287  const unsigned int i,
288  const Point & p);
289 
302  static Real shape(const FEType fet,
303  const Elem * elem,
304  const unsigned int i,
305  const Point & p,
306  const bool add_p_level);
307 
320  static Real shape_deriv (const FEType & fet,
321  const Elem * inf_elem,
322  const unsigned int i,
323  const unsigned int j,
324  const Point & p);
325 
338  static Real shape_deriv (const FEType fet,
339  const Elem * inf_elem,
340  const unsigned int i,
341  const unsigned int j,
342  const Point & p,
343  const bool add_p_level);
344 
357  static Real shape_deriv (const FEType & fet,
358  const ElemType inf_elem_type,
359  const unsigned int i,
360  const unsigned int j,
361  const Point & p);
362 
373  static void compute_data(const FEType & fe_t,
374  const Elem * inf_elem,
375  FEComputeData & data);
376 
381  static unsigned int n_shape_functions (const FEType & fet,
382  const Elem * inf_elem)
383  { return n_dofs(fet, inf_elem); }
384 
385 #ifdef LIBMESH_ENABLE_DEPRECATED
386  /*
387  * \deprecated Call the version of this function that takes a
388  * pointer-to-Elem instead.
389  */
390  static unsigned int n_shape_functions (const FEType & fet,
391  const ElemType t)
392  { return n_dofs(fet, t); }
393 
398  static unsigned int n_dofs(const FEType & fet,
399  const ElemType inf_elem_type);
400 #endif // LIBMESH_ENABLE_DEPRECATED
401 
408  static unsigned int n_dofs(const FEType & fet,
409  const Elem * inf_elem);
410 
411 #ifdef LIBMESH_ENABLE_DEPRECATED
412 
416  static unsigned int n_dofs_at_node(const FEType & fet,
417  const ElemType inf_elem_type,
418  const unsigned int n);
419 #endif // LIBMESH_ENABLE_DEPRECATED
420 
425  static unsigned int n_dofs_at_node(const FEType & fet,
426  const Elem * inf_elem,
427  const unsigned int n);
428 
429 #ifdef LIBMESH_ENABLE_DEPRECATED
430 
437  static unsigned int n_dofs_per_elem(const FEType & fet,
438  const ElemType inf_elem_type);
439 #endif // LIBMESH_ENABLE_DEPRECATED
440 
445  static unsigned int n_dofs_per_elem(const FEType & fet,
446  const Elem * inf_elem);
447 
451  virtual FEContinuity get_continuity() const override
452  { return C_ZERO; } // FIXME - is this true??
453 
458  virtual bool is_hierarchic() const override
459  { return false; } // FIXME - Inf FEs don't handle p elevation yet
460 
468  static void nodal_soln(const FEType & fet,
469  const Elem * elem,
470  const std::vector<Number> & elem_soln,
471  std::vector<Number> & nodal_soln);
472 
473 
474  static Point map (const Elem * inf_elem,
475  const Point & reference_point)
476  {
477  // libmesh_deprecated(); // soon
478  return InfFEMap::map(Dim, inf_elem, reference_point);
479  }
480 
481 
482  static Point inverse_map (const Elem * elem,
483  const Point & p,
484  const Real tolerance = TOLERANCE,
485  const bool secure = true)
486  {
487  // libmesh_deprecated(); // soon
488  return InfFEMap::inverse_map(Dim, elem, p, tolerance, secure);
489  }
490 
491 
492  static void inverse_map (const Elem * elem,
493  const std::vector<Point> & physical_points,
494  std::vector<Point> & reference_points,
495  const Real tolerance = TOLERANCE,
496  const bool secure = true)
497  {
498  // libmesh_deprecated(); // soon
499  return InfFEMap::inverse_map(Dim, elem, physical_points,
500  reference_points, tolerance, secure);
501  }
502 
503 
504  // The workhorses of InfFE. These are often used during matrix assembly.
505 
513  virtual void reinit (const Elem * elem,
514  const std::vector<Point> * const pts = nullptr,
515  const std::vector<Real> * const weights = nullptr) override;
516 
522  virtual void reinit (const Elem * inf_elem,
523  const unsigned int s,
524  const Real tolerance = TOLERANCE,
525  const std::vector<Point> * const pts = nullptr,
526  const std::vector<Real> * const weights = nullptr) override;
527 
533  virtual void edge_reinit (const Elem * elem,
534  const unsigned int edge,
535  const Real tolerance = TOLERANCE,
536  const std::vector<Point> * const pts = nullptr,
537  const std::vector<Real> * const weights = nullptr) override;
538 
543  virtual void side_map (const Elem * /* elem */,
544  const Elem * /* side */,
545  const unsigned int /* s */,
546  const std::vector<Point> & /* reference_side_points */,
547  std::vector<Point> & /* reference_points */) override
548  {
549  libmesh_not_implemented();
550  }
551 
563  virtual void attach_quadrature_rule (QBase * q) override;
564 
569  virtual unsigned int n_shape_functions () const override
570  { return _n_total_approx_sf; }
571 
577  virtual unsigned int n_quadrature_points () const override
578  { libmesh_assert(radial_qrule); return this->_n_total_qp; }
579 
584  virtual const std::vector<Point> & get_xyz () const override
586  calculate_xyz = true; return xyz; }
587 
594  virtual const std::vector<Real> & get_JxW () const override
595  {
597  calculate_jxw = true;
598  return this->JxW;
599  }
600 
609  virtual const std::vector<Real> & get_JxWxdecay_sq () const override
611  calculate_map_scaled = true; return this->JxWxdecay;}
612 
613 
624  virtual const std::vector<std::vector<OutputShape>> & get_phi_over_decayxR () const override
626  calculate_phi_scaled = true; return phixr; }
627 
628 
634  virtual const std::vector<std::vector<OutputGradient>> & get_dphi_over_decayxR () const override
636  calculate_dphi_scaled = true; return dphixr; }
637 
638 
646  virtual const std::vector<std::vector<OutputGradient>> & get_dphi_over_decay () const override
648  calculate_dphi_scaled = true; return dphixr_sq; }
649 
650 
655  virtual const std::vector<RealGradient> & get_dxyzdxi() const override
656  { calculate_map = true; libmesh_not_implemented();}
657 
658 
663  virtual const std::vector<RealGradient> & get_dxyzdeta() const override
664  { calculate_map = true; libmesh_not_implemented();}
665 
666 
671  virtual const std::vector<RealGradient> & get_dxyzdzeta() const override
672  { calculate_map = true; libmesh_not_implemented();}
673 
674 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
675 
678  virtual const std::vector<RealGradient> & get_d2xyzdxi2() const override
679  { calculate_map = true; libmesh_not_implemented();}
683  virtual const std::vector<RealGradient> & get_d2xyzdeta2() const override
684  { calculate_map = true; libmesh_not_implemented();}
688  virtual const std::vector<RealGradient> & get_d2xyzdzeta2() const override
689  { calculate_map = true; libmesh_not_implemented();}
693  virtual const std::vector<RealGradient> & get_d2xyzdxideta() const override
694  { calculate_map = true; libmesh_not_implemented();}
698  virtual const std::vector<RealGradient> & get_d2xyzdxidzeta() const override
699  { calculate_map = true; libmesh_not_implemented();}
703  virtual const std::vector<RealGradient> & get_d2xyzdetadzeta() const override
704  { calculate_map = true; libmesh_not_implemented();}
705 #endif
706 
711  virtual const std::vector<Real> & get_dxidx() const override
713  calculate_map = true; return dxidx_map;}
714 
715 
720  virtual const std::vector<Real> & get_dxidy() const override
722  calculate_map = true; return dxidy_map;}
723 
724 
729  virtual const std::vector<Real> & get_dxidz() const override
731  calculate_map = true; return dxidz_map;}
732 
733 
738  virtual const std::vector<Real> & get_detadx() const override
740  calculate_map = true; return detadx_map;}
741 
742 
747  virtual const std::vector<Real> & get_detady() const override
749  calculate_map = true; return detady_map;}
750 
751 
756  virtual const std::vector<Real> & get_detadz() const override
758  calculate_map = true; return detadz_map;}
759 
760 
765  virtual const std::vector<Real> & get_dzetadx() const override
767  calculate_map = true; return dzetadx_map;}
768 
769 
774  virtual const std::vector<Real> & get_dzetady() const override
776  calculate_map = true; return dzetady_map;}
777 
778 
783  virtual const std::vector<Real> & get_dzetadz() const override
785  calculate_map = true; return dzetadz_map;}
786 
787 
795  virtual const std::vector<Real> & get_Sobolev_weight() const override
797  calculate_phi = true; return weight; }
798 
799 
805  virtual const std::vector<RealGradient> & get_Sobolev_dweight() const override
807  calculate_dphi = true; return dweight; }
808 
809 
810 
814  virtual const std::vector<std::vector<Point>> & get_tangents() const override
816  calculate_map = true; return tangents; }
817 
821  virtual const std::vector<Point> & get_normals() const override
823  calculate_map = true; return normals; }
824 
825 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
826 
829  virtual const std::vector<Real> & get_curvatures() const override
830  { calculate_map = true; libmesh_not_implemented();}
831 #endif
832 
838  virtual const std::vector<Real> & get_Sobolev_weightxR_sq() const override
840  calculate_phi_scaled = true; return weightxr_sq; }
841 
842 
849  virtual const std::vector<RealGradient> & get_Sobolev_dweightxR_sq() const override
851  calculate_dphi_scaled = true; return dweightxr_sq; }
852 
853 
854 #ifdef LIBMESH_ENABLE_AMR
855 
861  static void inf_compute_constraints (DofConstraints & constraints,
862  DofMap & dof_map,
863  const unsigned int variable_number,
864  const Elem * child_elem);
865 #endif
866 
867 #ifdef LIBMESH_ENABLE_NODE_CONSTRAINTS
868  static void inf_compute_node_constraints (NodeConstraints & constraints,
869  const Elem * elem);
870 #endif
871 
872 
873 protected:
874 
875  // static members used by the workhorses
876 
893  static Real eval(Real v,
894  Order o_radial,
895  unsigned int i);
896 
902  static Real eval_deriv(Real v,
903  Order o_radial,
904  unsigned int i);
905 
906 
907 
908  // Non-static members used by the workhorses
909 
914  void update_base_elem (const Elem * inf_elem);
915 
919  virtual void init_base_shape_functions(const std::vector<Point> &,
920  const Elem *) override
921  { libmesh_not_implemented(); }
922 
927  virtual void determine_calculations() override;
928 
934  void init_radial_shape_functions(const Elem * inf_elem,
935  const std::vector<Point> * radial_pts = nullptr);
936 
943  void init_shape_functions(const std::vector<Point> & radial_qp,
944  const std::vector<Point> & base_qp,
945  const Elem * inf_elem);
946 
951  void init_face_shape_functions (const std::vector<Point> &,
952  const Elem * inf_side);
953 
964  void compute_shape_functions(const Elem * inf_elem,
965  const std::vector<Point> & base_qp,
966  const std::vector<Point> & radial_qp);
967 
968  void compute_face_functions();
969 
974  virtual void compute_shape_functions(const Elem *, const std::vector<Point> & ) override
975  {
976  //FIXME: it seems this function cannot be left out because
977  // it is pure virtual in \p FEBase
978  libmesh_not_implemented();
979  }
980 
984  mutable bool calculate_map_scaled;
985 
989  mutable bool calculate_phi_scaled;
990 
994  mutable bool calculate_dphi_scaled;
995 
996 
1000  mutable bool calculate_xyz;
1001 
1002 
1007  mutable bool calculate_jxw;
1008 
1009  // Miscellaneous static members
1010 
1017  static void compute_node_indices (const ElemType inf_elem_type,
1018  const unsigned int outer_node_index,
1019  unsigned int & base_node,
1020  unsigned int & radial_node);
1021 
1030  static void compute_node_indices_fast (const ElemType inf_elem_type,
1031  const unsigned int outer_node_index,
1032  unsigned int & base_node,
1033  unsigned int & radial_node);
1034 
1035 #ifdef LIBMESH_ENABLE_DEPRECATED
1036  /*
1037  * \deprecated Call the version of this function that takes an Elem * instead.
1038  */
1039  static void compute_shape_indices (const FEType & fet,
1040  const ElemType inf_elem_type,
1041  const unsigned int i,
1042  unsigned int & base_shape,
1043  unsigned int & radial_shape);
1044 #endif // LIBMESH_ENABLE_DEPRECATED
1045 
1052  static void compute_shape_indices (const FEType & fet,
1053  const Elem * inf_elem,
1054  const unsigned int i,
1055  unsigned int & base_shape,
1056  unsigned int & radial_shape);
1057 
1063  std::vector<Point> xyz;
1064 
1065  std::vector<Real> weightxr_sq;
1072  std::vector<Real> dweightdv;
1073 
1074  std::vector<RealGradient> dweightxr_sq;
1075 
1083  std::vector<Real> som;
1088  std::vector<Real> dsomdv;
1089 
1094  std::vector<std::vector<Real>> mode;
1095 
1100  std::vector<std::vector<Real>> dmodedv;
1101 
1102  // mapping of reference element to physical element
1103  // These vectors usually belong to \p this->fe_map
1104  // but for infinite elements, \p FEMap cannot
1105  // compute them.
1106  std::vector<Real> dxidx_map;
1107  std::vector<Real> dxidy_map;
1108  std::vector<Real> dxidz_map;
1109  std::vector<Real> detadx_map;
1110  std::vector<Real> detady_map;
1111  std::vector<Real> detadz_map;
1112  std::vector<Real> dzetadx_map;
1113  std::vector<Real> dzetady_map;
1114  std::vector<Real> dzetadz_map;
1115 
1116  // scaled mapping: Similar to the above
1117  // vectors, but scaled by radial (infinite) coordinate.
1118  std::vector<Real> dxidx_map_scaled;
1119  std::vector<Real> dxidy_map_scaled;
1120  std::vector<Real> dxidz_map_scaled;
1121  std::vector<Real> detadx_map_scaled;
1122  std::vector<Real> detady_map_scaled;
1123  std::vector<Real> detadz_map_scaled;
1124  std::vector<Real> dzetadx_map_scaled;
1125  std::vector<Real> dzetady_map_scaled;
1126  std::vector<Real> dzetadz_map_scaled;
1127 
1128 
1129  // respectively weighted shape functions.
1130  //FIXME: these names are correct only in 3D
1131  // but avoid long and clumbsy names...
1132  std::vector<std::vector<Real>> phixr;
1133  std::vector<std::vector<RealGradient>> dphixr;
1134  std::vector<std::vector<RealGradient>> dphixr_sq;
1135 
1136  std::vector<Real> JxWxdecay;
1137  std::vector<Real> JxW;
1138 
1139  std::vector<Point> normals;
1140  std::vector<std::vector<Point>> tangents;
1141 
1142  // numbering scheme maps
1143 
1152  std::vector<unsigned int> _radial_node_index;
1153 
1162  std::vector<unsigned int> _base_node_index;
1163 
1172  std::vector<unsigned int> _radial_shape_index;
1173 
1182  std::vector<unsigned int> _base_shape_index;
1183 
1184  // some more protected members
1185 
1190  unsigned int _n_total_approx_sf;
1191 
1196  std::vector<Real> _total_qrule_weights;
1197 
1202  std::unique_ptr<QBase> base_qrule;
1203 
1208  std::unique_ptr<QBase> radial_qrule;
1209 
1215  std::unique_ptr<const Elem> base_elem;
1216 
1223  std::unique_ptr<FEBase> base_fe;
1224 
1234 
1235 
1236 private:
1237 
1241  virtual bool shapes_need_reinit() const override;
1242 
1251 
1252 
1253 #ifdef DEBUG
1254 
1259  static bool _warned_for_shape;
1260  static bool _warned_for_dshape;
1261 
1262 #endif
1263 
1269  template <unsigned int friend_Dim, FEFamily friend_T_radial, InfMapType friend_T_map>
1270  friend class InfFE;
1271 
1272  friend class InfFEMap;
1273 };
1274 
1275 
1276 
1277 // InfFERadial class inline members
1278 // FIXME: decay in 3D is a/r
1279 // thus, this function fixes the mapping v<->r explicitly.
1280 // This is consistent with the current InfFEMap, which, however, is
1281 // written such that more general mappings can be implemented.
1282 inline
1283 Real InfFERadial::decay(const unsigned int dim, const Real v)
1284 {
1285  switch (dim)
1286  {
1287  case 3:
1288  return (1.-v)/2.;
1289 
1290  case 2:
1291  // according to P. Bettess "Infinite Elements",
1292  // the 2D case is rather tricky:
1293  // - the sqrt() requires special integration rules
1294  // if not both trial- and test function are involved
1295  // - the analytic solutions contain not only the sqrt, but
1296  // also a polynomial with rather many terms, so convergence
1297  // might be bad.
1298  return sqrt((1.-v)/2.);
1299 
1300  case 1:
1301  return 1.;
1302 
1303  default:
1304  libmesh_error_msg("Invalid dim = " << dim);
1305  }
1306 }
1307 
1308 inline
1309 Real InfFERadial::decay_deriv (const unsigned int dim, const Real v)
1310 {
1311  switch (dim)
1312  {
1313  case 3:
1314  return -0.5;
1315 
1316  case 2:
1317  // according to P. Bettess "Infinite Elements",
1318  // the 2D case is rather tricky:
1319  // - the sqrt() requires special integration rules
1320  // if not both trial- and test function are involved
1321  // - the analytic solutions contain not only the sqrt, but
1322  // also a polynomial with rather many terms, so convergence
1323  // might be bad.
1324 
1325  libmesh_assert (-1.-1.e-5 <= v && v < 1.);
1326  return -0.25/sqrt(1.-v);
1327 
1328  case 1:
1329  return 0.;
1330 
1331  default:
1332  libmesh_error_msg("Invalid dim = " << dim);
1333  }
1334 }
1335 
1336 } // namespace libMesh
1337 
1338 
1339 #endif // ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
1340 
1341 
1342 #endif // LIBMESH_INF_FE_H
std::vector< Real > weightxr_sq
Definition: inf_fe.h:1065
static void compute_node_indices_fast(const ElemType inf_elem_type, const unsigned int outer_node_index, unsigned int &base_node, unsigned int &radial_node)
Does the same as compute_node_indices(), but stores the maps for the current element type...
std::vector< Real > detadz_map
Definition: inf_fe.h:1111
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:196
std::vector< Real > dxidx_map
Definition: inf_fe.h:1106
bool calculate_map_scaled
Are we calculating scaled mapping functions?
Definition: inf_fe.h:984
virtual void determine_calculations() override
Determine which values are to be calculated, for both the FE itself and for the FEMap.
Definition: inf_fe.C:329
static Real shape_deriv(const FEType &fet, const Elem *inf_elem, const unsigned int i, const unsigned int j, const Point &p)
virtual const std::vector< RealGradient > & get_dxyzdxi() const override
Definition: inf_fe.h:655
ElemType
Defines an enum for geometric element types.
static ElemType get_elem_type(const ElemType type)
unsigned int _n_total_qp
The total number of quadrature points for the current configuration.
Definition: fe_abstract.h:774
virtual void init_base_shape_functions(const std::vector< Point > &, const Elem *) override
Do not use this derived member in InfFE<Dim,T_radial,T_map>.
Definition: inf_fe.h:919
virtual const std::vector< std::vector< Point > > & get_tangents() const override
Definition: inf_fe.h:814
static Point map(const Elem *inf_elem, const Point &reference_point)
Definition: inf_fe.h:474
static unsigned int n_dofs_per_elem(const FEType &fet, const ElemType inf_elem_type)
static unsigned int n_base_mapping_sf(const Elem &base_elem, const Order base_mapping_order)
std::vector< Real > dsomdv
the first local derivative of the radial decay in local coordinates.
Definition: inf_fe.h:1088
static bool _warned_for_shape
Definition: inf_fe.h:1259
bool calculations_started
Have calculations with this object already been started? Then all get_* functions should already have...
Definition: fe_abstract.h:666
Order
defines an enum for polynomial orders.
Definition: enum_order.h:40
bool calculate_phi
Should we calculate shape functions?
Definition: fe_abstract.h:691
A specific instantiation of the FEBase class.
Definition: fe.h:42
static unsigned int n_dofs(const Order o_radial)
Definition: inf_fe.h:113
std::vector< Real > detadx_map_scaled
Definition: inf_fe.h:1121
static unsigned int n_dofs_per_elem(const Order o_radial)
Definition: inf_fe.h:136
bool calculate_phi_scaled
Are we calculating scaled shape functions?
Definition: inf_fe.h:989
static unsigned int n_dofs_at_node(const FEType &fet, const ElemType inf_elem_type, const unsigned int n)
static Point map(const unsigned int dim, const Elem *inf_elem, const Point &reference_point)
Definition: inf_fe_map.C:40
std::vector< Real > dzetadz_map
Definition: inf_fe.h:1114
static Real eval_deriv(Real v, Order o_radial, unsigned int i)
std::vector< Real > dweightdv
the additional radial weight in local coordinates, over all quadrature points.
Definition: inf_fe.h:1072
static Real decay_deriv(const unsigned int dim, const Real)
Definition: inf_fe.h:1309
std::vector< Real > dxidx_map_scaled
Definition: inf_fe.h:1118
static constexpr Real TOLERANCE
std::vector< Real > dzetady_map
Definition: inf_fe.h:1113
unsigned int dim
std::vector< std::vector< Real > > phixr
Definition: inf_fe.h:1132
virtual const std::vector< RealGradient > & get_dxyzdeta() const override
Definition: inf_fe.h:663
std::vector< Real > dxidz_map
Definition: inf_fe.h:1108
static Real Dxr_sq(const Real)
Definition: inf_fe.h:84
std::unique_ptr< QBase > radial_qrule
The quadrature rule for the base element associated with the current infinite element.
Definition: inf_fe.h:1208
virtual const std::vector< Real > & get_detadz() const override
Definition: inf_fe.h:756
std::vector< unsigned int > _base_node_index
The internal structure of the InfFE – tensor product of base element times radial nodes – has to be...
Definition: inf_fe.h:1162
void init_shape_functions(const std::vector< Point > &radial_qp, const std::vector< Point > &base_qp, const Elem *inf_elem)
Initialize all the data fields like weight, mode, phi, dphidxi, dphideta, dphidzeta, etc.
Definition: inf_fe.C:482
virtual const std::vector< RealGradient > & get_d2xyzdxi2() const override
Definition: inf_fe.h:678
This is the base class from which all geometric element types are derived.
Definition: elem.h:94
std::vector< Real > dxidy_map_scaled
Definition: inf_fe.h:1119
std::vector< Real > dxidy_map
Definition: inf_fe.h:1107
InfFEBase()
Never use an object of this type.
Definition: inf_fe.h:158
virtual const std::vector< RealGradient > & get_d2xyzdzeta2() const override
Definition: inf_fe.h:688
virtual void compute_shape_functions(const Elem *, const std::vector< Point > &) override
Use compute_shape_functions(const Elem*, const std::vector<Point> &, const std::vector<Point> &) inst...
Definition: inf_fe.h:974
std::vector< Real > detadx_map
Definition: inf_fe.h:1109
std::vector< Real > weight
Used for certain infinite element families: the additional radial weight in local coordinates...
Definition: fe_base.h:767
virtual const std::vector< Real > & get_dxidx() const override
Definition: inf_fe.h:711
Class that encapsulates mapping (i.e.
Definition: inf_fe_map.h:47
The Node constraint storage format.
Definition: dof_map.h:156
static Real shape(const FEType &fet, const ElemType t, const unsigned int i, const Point &p)
std::vector< Real > detady_map_scaled
Definition: inf_fe.h:1122
static unsigned int n_shape_functions(const FEType &fet, const ElemType t)
Definition: inf_fe.h:390
The libMesh namespace provides an interface to certain functionality in the library.
static Point inverse_map(const unsigned int dim, const Elem *elem, const Point &p, const Real tolerance=TOLERANCE, const bool secure=true)
Definition: inf_fe_map.C:96
void init_face_shape_functions(const std::vector< Point > &, const Elem *inf_side)
Initialize all the data fields like weight, phi, etc for the side s.
static void compute_node_indices(const ElemType inf_elem_type, const unsigned int outer_node_index, unsigned int &base_node, unsigned int &radial_node)
Computes the indices in the base base_node and in radial direction radial_node (either 0 or 1) associ...
static std::unique_ptr< const Elem > build_elem(const Elem *inf_elem)
Build the base element of an infinite element.
virtual FEContinuity get_continuity() const override
Definition: inf_fe.h:451
static unsigned int n_dofs_at_node(const Order o_radial, const unsigned int n_onion)
std::vector< std::vector< Real > > mode
the radial approximation shapes in local coordinates Needed when setting up the overall shape functio...
Definition: inf_fe.h:1094
This nested class contains most of the static methods related to the base part of an infinite element...
Definition: inf_fe.h:151
void compute_face_functions()
virtual const std::vector< Real > & get_dzetady() const override
Definition: inf_fe.h:774
std::vector< Real > dzetady_map_scaled
Definition: inf_fe.h:1125
std::vector< RealGradient > dweightxr_sq
Definition: inf_fe.h:1074
std::vector< std::vector< RealGradient > > dphixr
Definition: inf_fe.h:1133
static Real D(const Real v)
Definition: inf_fe.h:82
static bool _warned_for_dshape
Definition: inf_fe.h:1260
~InfFE()=default
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:179
virtual const std::vector< Real > & get_dzetadz() const override
Definition: inf_fe.h:783
bool calculate_xyz
Are we calculating the positions of quadrature points?
Definition: inf_fe.h:1000
std::vector< Real > detadz_map_scaled
Definition: inf_fe.h:1123
static void inf_compute_constraints(DofConstraints &constraints, DofMap &dof_map, const unsigned int variable_number, const Elem *child_elem)
Computes the constraint matrix contributions (for non-conforming adapted meshes) corresponding to var...
std::vector< Real > JxWxdecay
Definition: inf_fe.h:1136
bool calculate_jxw
Are we calculating the unscaled jacobian? We avoid it if not requested explicitly; this has the worst...
Definition: inf_fe.h:1007
friend class InfFE
Make all InfFE<Dim,T_radial,T_map> classes friends of each other, so that the protected eval() may be...
Definition: inf_fe.h:1270
virtual const std::vector< Real > & get_JxWxdecay_sq() const override
Definition: inf_fe.h:609
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: inf_fe.h:492
virtual const std::vector< Real > & get_Sobolev_weight() const override
Definition: inf_fe.h:795
static unsigned int n_dofs(const FEType &fet, const ElemType inf_elem_type)
Definition: inf_fe_static.C:67
virtual unsigned int n_quadrature_points() const override
Definition: inf_fe.h:577
static Order mapping_order()
Definition: inf_fe.h:97
std::unique_ptr< FEBase > base_fe
Have a FE<Dim-1,T_base> handy for base approximation.
Definition: inf_fe.h:1223
std::vector< Real > JxW
Definition: inf_fe.h:1137
virtual const std::vector< Real > & get_dxidz() const override
Definition: inf_fe.h:729
libmesh_assert(ctx)
static Real decay(const unsigned int dim, const Real v)
Definition: inf_fe.h:1283
static Real eval(Real v, Order o_radial, unsigned int i)
FEType current_fe_type
This FEType stores the characteristics for which the data structures phi, phi_map etc are currently i...
Definition: inf_fe.h:1233
virtual const std::vector< Real > & get_detadx() const override
Definition: inf_fe.h:738
virtual const std::vector< RealGradient > & get_d2xyzdxideta() const override
Definition: inf_fe.h:693
bool calculate_dphi_scaled
Are we calculating scaled shape function gradients?
Definition: inf_fe.h:994
FEGenericBase< Real > FEBase
virtual const std::vector< Real > & get_dzetadx() const override
Definition: inf_fe.h:765
virtual bool shapes_need_reinit() const override
Definition: inf_fe.C:1152
unsigned int _n_total_approx_sf
The number of total approximation shape functions for the current configuration.
Definition: inf_fe.h:1190
void init_radial_shape_functions(const Elem *inf_elem, const std::vector< Point > *radial_pts=nullptr)
Some of the member data only depend on the radial part of the infinite element.
Definition: inf_fe.C:418
std::vector< unsigned int > _radial_node_index
The internal structure of the InfFE – tensor product of base element times radial nodes – has to be...
Definition: inf_fe.h:1152
void compute_shape_functions(const Elem *inf_elem, const std::vector< Point > &base_qp, const std::vector< Point > &radial_qp)
After having updated the jacobian and the transformation from local to global coordinates in FEAbstra...
Definition: inf_fe.C:784
virtual const std::vector< Real > & get_dxidy() const override
Definition: inf_fe.h:720
std::vector< Real > som
the radial decay in local coordinates.
Definition: inf_fe.h:1083
static Point inverse_map(const Elem *elem, const Point &p, const Real tolerance=TOLERANCE, const bool secure=true)
Definition: inf_fe.h:482
virtual const std::vector< std::vector< OutputGradient > > & get_dphi_over_decay() const override
Definition: inf_fe.h:646
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
bool calculate_dphi
Should we calculate shape function gradients?
Definition: fe_abstract.h:696
virtual const std::vector< Real > & get_Sobolev_weightxR_sq() const override
Definition: inf_fe.h:838
bool calculate_map
Are we calculating mapping functions?
Definition: fe_abstract.h:686
virtual const std::vector< RealGradient > & get_d2xyzdeta2() const override
Definition: inf_fe.h:683
std::vector< Real > dzetadx_map
Definition: inf_fe.h:1112
FEContinuity
defines an enum for finite element types to libmesh_assert a certain level (or type? Hcurl?) of continuity.
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
Not implemented yet.
std::unique_ptr< QBase > base_qrule
The quadrature rule for the base element associated with the current infinite element.
Definition: inf_fe.h:1202
static void compute_shape_indices(const FEType &fet, const ElemType inf_elem_type, const unsigned int i, unsigned int &base_shape, unsigned int &radial_shape)
std::vector< Real > dxidz_map_scaled
Definition: inf_fe.h:1120
std::vector< unsigned int > _radial_shape_index
The internal structure of the InfFE – tensor product of base element shapes times radial shapes – h...
Definition: inf_fe.h:1172
static Real D_deriv(const Real v)
Definition: inf_fe.h:90
std::vector< Point > xyz
Physical quadrature points.
Definition: inf_fe.h:1063
virtual const std::vector< RealGradient > & get_Sobolev_dweight() const override
Definition: inf_fe.h:805
std::vector< Point > normals
Definition: inf_fe.h:1139
InfFERadial()
Never use an object of this type.
Definition: inf_fe.h:60
std::vector< std::vector< RealGradient > > dphixr_sq
Definition: inf_fe.h:1134
static void inf_compute_node_constraints(NodeConstraints &constraints, const Elem *elem)
static bool _warned_for_nodal_soln
static members that are used to issue warning messages only once.
Definition: inf_fe.h:1258
static void compute_data(const FEType &fe_t, const Elem *inf_elem, FEComputeData &data)
Generalized version of shape(), takes an Elem *.
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: inf_fe.C:120
virtual const std::vector< std::vector< OutputShape > > & get_phi_over_decayxR() const override
Definition: inf_fe.h:624
virtual const std::vector< RealGradient > & get_d2xyzdetadzeta() const override
Definition: inf_fe.h:703
virtual const std::vector< Real > & get_detady() const override
Definition: inf_fe.h:747
static void nodal_soln(const FEType &fet, const Elem *elem, const std::vector< Number > &elem_soln, std::vector< Number > &nodal_soln)
Usually, this method would build the nodal soln from the element soln.
virtual const std::vector< Real > & get_JxW() const override
Definition: inf_fe.h:594
std::vector< unsigned int > _base_shape_index
The internal structure of the InfFE – tensor product of base element shapes times radial shapes – h...
Definition: inf_fe.h:1182
std::vector< Real > dzetadz_map_scaled
Definition: inf_fe.h:1126
virtual const std::vector< std::vector< OutputGradient > > & get_dphi_over_decayxR() const override
Definition: inf_fe.h:634
virtual void attach_quadrature_rule(QBase *q) override
The use of quadrature rules with the InfFE class is somewhat different from the approach of the FE cl...
Definition: inf_fe.C:79
virtual const std::vector< RealGradient > & get_d2xyzdxidzeta() const override
Definition: inf_fe.h:698
virtual void side_map(const Elem *, const Elem *, const unsigned int, const std::vector< Point > &, std::vector< Point > &) override
Computes the reference space quadrature points on the side of an element based on the side quadrature...
Definition: inf_fe.h:543
static ElemType _compute_node_indices_fast_current_elem_type
When compute_node_indices_fast() is used, this static variable remembers the element type for which t...
Definition: inf_fe.h:1250
virtual const std::vector< Point > & get_normals() const override
Definition: inf_fe.h:821
std::vector< std::vector< Real > > dmodedv
the first local derivative of the radial approximation shapes.
Definition: inf_fe.h:1100
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:39
virtual unsigned int n_shape_functions() const override
Definition: inf_fe.h:569
virtual bool is_hierarchic() const override
Definition: inf_fe.h:458
The constraint matrix storage format.
Definition: dof_map.h:108
std::vector< std::vector< Point > > tangents
Definition: inf_fe.h:1140
The QBase class provides the basic functionality from which various quadrature rules can be derived...
Definition: quadrature.h:61
std::vector< Real > detady_map
Definition: inf_fe.h:1110
virtual const std::vector< Point > & get_xyz() const override
Definition: inf_fe.h:584
Infinite elements are in some sense directional, compared to conventional finite elements.
Definition: inf_fe.h:53
std::vector< Real > dzetadx_map_scaled
Definition: inf_fe.h:1124
std::vector< Real > _total_qrule_weights
this vector contains the combined integration weights, so that FEAbstract::compute_map() can still be...
Definition: inf_fe.h:1196
virtual const std::vector< RealGradient > & get_Sobolev_dweightxR_sq() const override
Definition: inf_fe.h:849
virtual const std::vector< Real > & get_curvatures() const override
Definition: inf_fe.h:829
static unsigned int n_shape_functions(const FEType &fet, const Elem *inf_elem)
Definition: inf_fe.h:381
void update_base_elem(const Elem *inf_elem)
Updates the protected member base_elem to the appropriate base element for the given inf_elem...
Definition: inf_fe.C:109
std::vector< RealGradient > dweight
Used for certain infinite element families: the global derivative of the additional radial weight ...
Definition: fe_base.h:760
std::unique_ptr< const Elem > base_elem
The "base" (aka non-infinite) element associated with the current infinite element.
Definition: inf_fe.h:1215
virtual const std::vector< RealGradient > & get_dxyzdzeta() const override
Definition: inf_fe.h:671