libMesh
inf_fe.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2024 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 
220 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map>
221 class InfFE : public FEBase
222 {
223 
224  /*
225  * Protect the nested class
226  */
227 protected:
228 
229 
230 public:
231 
232  // InfFE continued
233 
247  explicit
248  InfFE(const FEType & fet);
249  ~InfFE() = default;
250 
251  // The static public members for access from FEInterface etc
252 
265  static Real shape(const FEType & fet,
266  const ElemType t,
267  const unsigned int i,
268  const Point & p);
269 
282  static Real shape(const FEType & fet,
283  const Elem * elem,
284  const unsigned int i,
285  const Point & p);
286 
299  static Real shape(const FEType fet,
300  const Elem * elem,
301  const unsigned int i,
302  const Point & p,
303  const bool add_p_level);
304 
317  static Real shape_deriv (const FEType & fet,
318  const Elem * inf_elem,
319  const unsigned int i,
320  const unsigned int j,
321  const Point & p);
322 
335  static Real shape_deriv (const FEType fet,
336  const Elem * inf_elem,
337  const unsigned int i,
338  const unsigned int j,
339  const Point & p,
340  const bool add_p_level);
341 
354  static Real shape_deriv (const FEType & fet,
355  const ElemType inf_elem_type,
356  const unsigned int i,
357  const unsigned int j,
358  const Point & p);
359 
370  static void compute_data(const FEType & fe_t,
371  const Elem * inf_elem,
372  FEComputeData & data);
373 
381  static unsigned int n_shape_functions (const FEType & fet,
382  const ElemType t)
383  { return n_dofs(fet, t); }
384 
385  static unsigned int n_shape_functions (const FEType & fet,
386  const Elem * inf_elem)
387  { return n_dofs(fet, inf_elem); }
388 
398  static unsigned int n_dofs(const FEType & fet,
399  const ElemType inf_elem_type);
400 
401  static unsigned int n_dofs(const FEType & fet,
402  const Elem * inf_elem);
403 
411  static unsigned int n_dofs_at_node(const FEType & fet,
412  const ElemType inf_elem_type,
413  const unsigned int n);
414 
415  static unsigned int n_dofs_at_node(const FEType & fet,
416  const Elem * inf_elem,
417  const unsigned int n);
418 
426  static unsigned int n_dofs_per_elem(const FEType & fet,
427  const ElemType inf_elem_type);
428 
429  static unsigned int n_dofs_per_elem(const FEType & fet,
430  const Elem * inf_elem);
431 
435  virtual FEContinuity get_continuity() const override
436  { return C_ZERO; } // FIXME - is this true??
437 
442  virtual bool is_hierarchic() const override
443  { return false; } // FIXME - Inf FEs don't handle p elevation yet
444 
452  static void nodal_soln(const FEType & fet,
453  const Elem * elem,
454  const std::vector<Number> & elem_soln,
455  std::vector<Number> & nodal_soln);
456 
457 
458  static Point map (const Elem * inf_elem,
459  const Point & reference_point)
460  {
461  // libmesh_deprecated(); // soon
462  return InfFEMap::map(Dim, inf_elem, reference_point);
463  }
464 
465 
466  static Point inverse_map (const Elem * elem,
467  const Point & p,
468  const Real tolerance = TOLERANCE,
469  const bool secure = true)
470  {
471  // libmesh_deprecated(); // soon
472  return InfFEMap::inverse_map(Dim, elem, p, tolerance, secure);
473  }
474 
475 
476  static void inverse_map (const Elem * elem,
477  const std::vector<Point> & physical_points,
478  std::vector<Point> & reference_points,
479  const Real tolerance = TOLERANCE,
480  const bool secure = true)
481  {
482  // libmesh_deprecated(); // soon
483  return InfFEMap::inverse_map(Dim, elem, physical_points,
484  reference_points, tolerance, secure);
485  }
486 
487 
488  // The workhorses of InfFE. These are often used during matrix assembly.
489 
497  virtual void reinit (const Elem * elem,
498  const std::vector<Point> * const pts = nullptr,
499  const std::vector<Real> * const weights = nullptr) override;
500 
506  virtual void reinit (const Elem * inf_elem,
507  const unsigned int s,
508  const Real tolerance = TOLERANCE,
509  const std::vector<Point> * const pts = nullptr,
510  const std::vector<Real> * const weights = nullptr) override;
511 
517  virtual void edge_reinit (const Elem * elem,
518  const unsigned int edge,
519  const Real tolerance = TOLERANCE,
520  const std::vector<Point> * const pts = nullptr,
521  const std::vector<Real> * const weights = nullptr) override;
522 
527  virtual void side_map (const Elem * /* elem */,
528  const Elem * /* side */,
529  const unsigned int /* s */,
530  const std::vector<Point> & /* reference_side_points */,
531  std::vector<Point> & /* reference_points */) override
532  {
533  libmesh_not_implemented();
534  }
535 
547  virtual void attach_quadrature_rule (QBase * q) override;
548 
553  virtual unsigned int n_shape_functions () const override
554  { return _n_total_approx_sf; }
555 
561  virtual unsigned int n_quadrature_points () const override
563 
564 
569  virtual const std::vector<Point> & get_xyz () const override
571  calculate_xyz = true; return xyz; }
572 
579  virtual const std::vector<Real> & get_JxW () const override
580  {
582  calculate_jxw = true;
583  return this->JxW;
584  }
585 
594  virtual const std::vector<Real> & get_JxWxdecay_sq () const override
596  calculate_map_scaled = true; return this->JxWxdecay;}
597 
598 
609  virtual const std::vector<std::vector<OutputShape>> & get_phi_over_decayxR () const override
611  calculate_phi_scaled = true; return phixr; }
612 
613 
619  virtual const std::vector<std::vector<OutputGradient>> & get_dphi_over_decayxR () const override
621  calculate_dphi_scaled = true; return dphixr; }
622 
623 
631  virtual const std::vector<std::vector<OutputGradient>> & get_dphi_over_decay () const override
633  calculate_dphi_scaled = true; return dphixr_sq; }
634 
635 
640  virtual const std::vector<RealGradient> & get_dxyzdxi() const override
641  { calculate_map = true; libmesh_not_implemented();}
642 
643 
648  virtual const std::vector<RealGradient> & get_dxyzdeta() const override
649  { calculate_map = true; libmesh_not_implemented();}
650 
651 
656  virtual const std::vector<RealGradient> & get_dxyzdzeta() const override
657  { calculate_map = true; libmesh_not_implemented();}
658 
659 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
660 
663  virtual const std::vector<RealGradient> & get_d2xyzdxi2() const override
664  { calculate_map = true; libmesh_not_implemented();}
668  virtual const std::vector<RealGradient> & get_d2xyzdeta2() const override
669  { calculate_map = true; libmesh_not_implemented();}
673  virtual const std::vector<RealGradient> & get_d2xyzdzeta2() const override
674  { calculate_map = true; libmesh_not_implemented();}
678  virtual const std::vector<RealGradient> & get_d2xyzdxideta() const override
679  { calculate_map = true; libmesh_not_implemented();}
683  virtual const std::vector<RealGradient> & get_d2xyzdxidzeta() const override
684  { calculate_map = true; libmesh_not_implemented();}
688  virtual const std::vector<RealGradient> & get_d2xyzdetadzeta() const override
689  { calculate_map = true; libmesh_not_implemented();}
690 #endif
691 
696  virtual const std::vector<Real> & get_dxidx() const override
698  calculate_map = true; return dxidx_map;}
699 
700 
705  virtual const std::vector<Real> & get_dxidy() const override
707  calculate_map = true; return dxidy_map;}
708 
709 
714  virtual const std::vector<Real> & get_dxidz() const override
716  calculate_map = true; return dxidz_map;}
717 
718 
723  virtual const std::vector<Real> & get_detadx() const override
725  calculate_map = true; return detadx_map;}
726 
727 
732  virtual const std::vector<Real> & get_detady() const override
734  calculate_map = true; return detady_map;}
735 
736 
741  virtual const std::vector<Real> & get_detadz() const override
743  calculate_map = true; return detadz_map;}
744 
745 
750  virtual const std::vector<Real> & get_dzetadx() const override
752  calculate_map = true; return dzetadx_map;}
753 
754 
759  virtual const std::vector<Real> & get_dzetady() const override
761  calculate_map = true; return dzetady_map;}
762 
763 
768  virtual const std::vector<Real> & get_dzetadz() const override
770  calculate_map = true; return dzetadz_map;}
771 
772 
780  virtual const std::vector<Real> & get_Sobolev_weight() const override
782  calculate_phi = true; return weight; }
783 
784 
790  virtual const std::vector<RealGradient> & get_Sobolev_dweight() const override
792  calculate_dphi = true; return dweight; }
793 
794 
795 
799  virtual const std::vector<std::vector<Point>> & get_tangents() const override
801  calculate_map = true; return tangents; }
802 
806  virtual const std::vector<Point> & get_normals() const override
808  calculate_map = true; return normals; }
809 
810 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
811 
814  virtual const std::vector<Real> & get_curvatures() const override
815  { calculate_map = true; libmesh_not_implemented();}
816 #endif
817 
823  virtual const std::vector<Real> & get_Sobolev_weightxR_sq() const override
825  calculate_phi_scaled = true; return weightxr_sq; }
826 
827 
834  virtual const std::vector<RealGradient> & get_Sobolev_dweightxR_sq() const override
836  calculate_dphi_scaled = true; return dweightxr_sq; }
837 
838 
839 #ifdef LIBMESH_ENABLE_AMR
840 
846  static void inf_compute_constraints (DofConstraints & constraints,
847  DofMap & dof_map,
848  const unsigned int variable_number,
849  const Elem * child_elem);
850 #endif
851 
852 #ifdef LIBMESH_ENABLE_NODE_CONSTRAINTS
853  static void inf_compute_node_constraints (NodeConstraints & constraints,
854  const Elem * elem);
855 #endif
856 
857 
858 protected:
859 
860  // static members used by the workhorses
861 
878  static Real eval(Real v,
879  Order o_radial,
880  unsigned int i);
881 
887  static Real eval_deriv(Real v,
888  Order o_radial,
889  unsigned int i);
890 
891 
892 
893  // Non-static members used by the workhorses
894 
899  void update_base_elem (const Elem * inf_elem);
900 
904  virtual void init_base_shape_functions(const std::vector<Point> &,
905  const Elem *) override
906  { libmesh_not_implemented(); }
907 
912  virtual void determine_calculations() override;
913 
919  void init_radial_shape_functions(const Elem * inf_elem,
920  const std::vector<Point> * radial_pts = nullptr);
921 
928  void init_shape_functions(const std::vector<Point> & radial_qp,
929  const std::vector<Point> & base_qp,
930  const Elem * inf_elem);
931 
936  void init_face_shape_functions (const std::vector<Point> &,
937  const Elem * inf_side);
938 
949  void compute_shape_functions(const Elem * inf_elem,
950  const std::vector<Point> & base_qp,
951  const std::vector<Point> & radial_qp);
952 
953  void compute_face_functions();
954 
959  virtual void compute_shape_functions(const Elem *, const std::vector<Point> & ) override
960  {
961  //FIXME: it seems this function cannot be left out because
962  // it is pure virtual in \p FEBase
963  libmesh_not_implemented();
964  }
965 
969  mutable bool calculate_map_scaled;
970 
974  mutable bool calculate_phi_scaled;
975 
979  mutable bool calculate_dphi_scaled;
980 
981 
985  mutable bool calculate_xyz;
986 
987 
992  mutable bool calculate_jxw;
993 
994  // Miscellaneous static members
995 
1002  static void compute_node_indices (const ElemType inf_elem_type,
1003  const unsigned int outer_node_index,
1004  unsigned int & base_node,
1005  unsigned int & radial_node);
1006 
1015  static void compute_node_indices_fast (const ElemType inf_elem_type,
1016  const unsigned int outer_node_index,
1017  unsigned int & base_node,
1018  unsigned int & radial_node);
1019 
1028  static void compute_shape_indices (const FEType & fet,
1029  const ElemType inf_elem_type,
1030  const unsigned int i,
1031  unsigned int & base_shape,
1032  unsigned int & radial_shape);
1033 
1034  static void compute_shape_indices (const FEType & fet,
1035  const Elem * inf_elem,
1036  const unsigned int i,
1037  unsigned int & base_shape,
1038  unsigned int & radial_shape);
1039 
1045  std::vector<Point> xyz;
1046 
1047  std::vector<Real> weightxr_sq;
1054  std::vector<Real> dweightdv;
1055 
1056  std::vector<RealGradient> dweightxr_sq;
1057 
1065  std::vector<Real> som;
1070  std::vector<Real> dsomdv;
1071 
1076  std::vector<std::vector<Real>> mode;
1077 
1082  std::vector<std::vector<Real>> dmodedv;
1083 
1084  // mapping of reference element to physical element
1085  // These vectors usually belong to \p this->fe_map
1086  // but for infinite elements, \p FEMap cannot
1087  // compute them.
1088  std::vector<Real> dxidx_map;
1089  std::vector<Real> dxidy_map;
1090  std::vector<Real> dxidz_map;
1091  std::vector<Real> detadx_map;
1092  std::vector<Real> detady_map;
1093  std::vector<Real> detadz_map;
1094  std::vector<Real> dzetadx_map;
1095  std::vector<Real> dzetady_map;
1096  std::vector<Real> dzetadz_map;
1097 
1098  // scaled mapping: Similar to the above
1099  // vectors, but scaled by radial (infinite) coordinate.
1100  std::vector<Real> dxidx_map_scaled;
1101  std::vector<Real> dxidy_map_scaled;
1102  std::vector<Real> dxidz_map_scaled;
1103  std::vector<Real> detadx_map_scaled;
1104  std::vector<Real> detady_map_scaled;
1105  std::vector<Real> detadz_map_scaled;
1106  std::vector<Real> dzetadx_map_scaled;
1107  std::vector<Real> dzetady_map_scaled;
1108  std::vector<Real> dzetadz_map_scaled;
1109 
1110 
1111  // respectively weighted shape functions.
1112  //FIXME: these names are correct only in 3D
1113  // but avoid long and clumbsy names...
1114  std::vector<std::vector<Real>> phixr;
1115  std::vector<std::vector<RealGradient>> dphixr;
1116  std::vector<std::vector<RealGradient>> dphixr_sq;
1117 
1118  std::vector<Real> JxWxdecay;
1119  std::vector<Real> JxW;
1120 
1121  std::vector<Point> normals;
1122  std::vector<std::vector<Point>> tangents;
1123 
1124  // numbering scheme maps
1125 
1134  std::vector<unsigned int> _radial_node_index;
1135 
1144  std::vector<unsigned int> _base_node_index;
1145 
1154  std::vector<unsigned int> _radial_shape_index;
1155 
1164  std::vector<unsigned int> _base_shape_index;
1165 
1166  // some more protected members
1167 
1172  unsigned int _n_total_approx_sf;
1173 
1178  unsigned int _n_total_qp;
1179 
1184  std::vector<Real> _total_qrule_weights;
1185 
1190  std::unique_ptr<QBase> base_qrule;
1191 
1196  std::unique_ptr<QBase> radial_qrule;
1197 
1203  std::unique_ptr<const Elem> base_elem;
1204 
1211  std::unique_ptr<FEBase> base_fe;
1212 
1222 
1223 
1224 private:
1225 
1229  virtual bool shapes_need_reinit() const override;
1230 
1239 
1240 
1241 #ifdef DEBUG
1242 
1247  static bool _warned_for_shape;
1248  static bool _warned_for_dshape;
1249 
1250 #endif
1251 
1257  template <unsigned int friend_Dim, FEFamily friend_T_radial, InfMapType friend_T_map>
1258  friend class InfFE;
1259 
1260  friend class InfFEMap;
1261 };
1262 
1263 
1264 
1265 // InfFERadial class inline members
1266 // FIXME: decay in 3D is a/r
1267 // thus, this function fixes the mapping v<->r explicitly.
1268 // This is consistent with the current InfFEMap, which, however, is
1269 // written such that more general mappings can be implemented.
1270 inline
1271 Real InfFERadial::decay(const unsigned int dim, const Real v)
1272 {
1273  switch (dim)
1274  {
1275  case 3:
1276  return (1.-v)/2.;
1277 
1278  case 2:
1279  // according to P. Bettess "Infinite Elements",
1280  // the 2D case is rather tricky:
1281  // - the sqrt() requires special integration rules
1282  // if not both trial- and test function are involved
1283  // - the analytic solutions contain not only the sqrt, but
1284  // also a polynomial with rather many terms, so convergence
1285  // might be bad.
1286  return sqrt((1.-v)/2.);
1287 
1288  case 1:
1289  return 1.;
1290 
1291  default:
1292  libmesh_error_msg("Invalid dim = " << dim);
1293  }
1294 }
1295 
1296 inline
1297 Real InfFERadial::decay_deriv (const unsigned int dim, const Real v)
1298 {
1299  switch (dim)
1300  {
1301  case 3:
1302  return -0.5;
1303 
1304  case 2:
1305  // according to P. Bettess "Infinite Elements",
1306  // the 2D case is rather tricky:
1307  // - the sqrt() requires special integration rules
1308  // if not both trial- and test function are involved
1309  // - the analytic solutions contain not only the sqrt, but
1310  // also a polynomial with rather many terms, so convergence
1311  // might be bad.
1312 
1313  libmesh_assert (-1.-1.e-5 <= v && v < 1.);
1314  return -0.25/sqrt(1.-v);
1315 
1316  case 1:
1317  return 0.;
1318 
1319  default:
1320  libmesh_error_msg("Invalid dim = " << dim);
1321  }
1322 }
1323 
1324 } // namespace libMesh
1325 
1326 
1327 #endif // ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
1328 
1329 
1330 #endif // LIBMESH_INF_FE_H
std::vector< Real > weightxr_sq
Definition: inf_fe.h:1047
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:1093
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:182
std::vector< Real > dxidx_map
Definition: inf_fe.h:1088
bool calculate_map_scaled
Are we calculating scaled mapping functions?
Definition: inf_fe.h:969
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:326
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:640
ElemType
Defines an enum for geometric element types.
static ElemType get_elem_type(const ElemType type)
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:904
virtual const std::vector< std::vector< Point > > & get_tangents() const override
Definition: inf_fe.h:799
static Point map(const Elem *inf_elem, const Point &reference_point)
Definition: inf_fe.h:458
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:1070
static bool _warned_for_shape
Definition: inf_fe.h:1247
bool calculations_started
Have calculations with this object already been started? Then all get_* functions should already have...
Definition: fe_abstract.h:645
Order
defines an enum for polynomial orders.
Definition: enum_order.h:40
bool calculate_phi
Should we calculate shape functions?
Definition: fe_abstract.h:670
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:1103
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:974
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:1096
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:1054
static Real decay_deriv(const unsigned int dim, const Real)
Definition: inf_fe.h:1297
std::vector< Real > dxidx_map_scaled
Definition: inf_fe.h:1100
static constexpr Real TOLERANCE
std::vector< Real > dzetady_map
Definition: inf_fe.h:1095
unsigned int dim
std::vector< std::vector< Real > > phixr
Definition: inf_fe.h:1114
virtual const std::vector< RealGradient > & get_dxyzdeta() const override
Definition: inf_fe.h:648
std::vector< Real > dxidz_map
Definition: inf_fe.h:1090
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:1196
virtual const std::vector< Real > & get_detadz() const override
Definition: inf_fe.h:741
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:1144
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:479
virtual const std::vector< RealGradient > & get_d2xyzdxi2() const override
Definition: inf_fe.h:663
unsigned int _n_total_qp
The total number of quadrature points for the current configuration.
Definition: inf_fe.h:1178
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:1101
std::vector< Real > dxidy_map
Definition: inf_fe.h:1089
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:673
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:959
std::vector< Real > detadx_map
Definition: inf_fe.h:1091
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:696
Class that encapsulates mapping (i.e.
Definition: inf_fe_map.h:47
The Node constraint storage format.
Definition: dof_map.h:146
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:1104
ADRealEigenVector< T, D, asd > sqrt(const ADRealEigenVector< T, D, asd > &)
Definition: type_vector.h:53
static unsigned int n_shape_functions(const FEType &fet, const ElemType t)
Definition: inf_fe.h:381
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:435
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:1076
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:759
std::vector< Real > dzetady_map_scaled
Definition: inf_fe.h:1107
std::vector< RealGradient > dweightxr_sq
Definition: inf_fe.h:1056
std::vector< std::vector< RealGradient > > dphixr
Definition: inf_fe.h:1115
static Real D(const Real v)
Definition: inf_fe.h:82
static bool _warned_for_dshape
Definition: inf_fe.h:1248
~InfFE()=default
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:169
virtual const std::vector< Real > & get_dzetadz() const override
Definition: inf_fe.h:768
bool calculate_xyz
Are we calculating the positions of quadrature points?
Definition: inf_fe.h:985
std::vector< Real > detadz_map_scaled
Definition: inf_fe.h:1105
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:1118
bool calculate_jxw
Are we calculating the unscaled jacobian? We avoid it if not requested explicitly; this has the worst...
Definition: inf_fe.h:992
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:1258
virtual const std::vector< Real > & get_JxWxdecay_sq() const override
Definition: inf_fe.h:594
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:476
virtual const std::vector< Real > & get_Sobolev_weight() const override
Definition: inf_fe.h:780
static unsigned int n_dofs(const FEType &fet, const ElemType inf_elem_type)
Definition: inf_fe_static.C:66
virtual unsigned int n_quadrature_points() const override
Definition: inf_fe.h:561
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:1211
std::vector< Real > JxW
Definition: inf_fe.h:1119
virtual const std::vector< Real > & get_dxidz() const override
Definition: inf_fe.h:714
libmesh_assert(ctx)
static Real decay(const unsigned int dim, const Real v)
Definition: inf_fe.h:1271
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:1221
virtual const std::vector< Real > & get_detadx() const override
Definition: inf_fe.h:723
virtual const std::vector< RealGradient > & get_d2xyzdxideta() const override
Definition: inf_fe.h:678
bool calculate_dphi_scaled
Are we calculating scaled shape function gradients?
Definition: inf_fe.h:979
FEGenericBase< Real > FEBase
virtual const std::vector< Real > & get_dzetadx() const override
Definition: inf_fe.h:750
virtual bool shapes_need_reinit() const override
Definition: inf_fe.C:1148
unsigned int _n_total_approx_sf
The number of total approximation shape functions for the current configuration.
Definition: inf_fe.h:1172
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:415
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:1134
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:780
virtual const std::vector< Real > & get_dxidy() const override
Definition: inf_fe.h:705
std::vector< Real > som
the radial decay in local coordinates.
Definition: inf_fe.h:1065
static Point inverse_map(const Elem *elem, const Point &p, const Real tolerance=TOLERANCE, const bool secure=true)
Definition: inf_fe.h:466
virtual const std::vector< std::vector< OutputGradient > > & get_dphi_over_decay() const override
Definition: inf_fe.h:631
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
bool calculate_dphi
Should we calculate shape function gradients?
Definition: fe_abstract.h:675
virtual const std::vector< Real > & get_Sobolev_weightxR_sq() const override
Definition: inf_fe.h:823
bool calculate_map
Are we calculating mapping functions?
Definition: fe_abstract.h:665
virtual const std::vector< RealGradient > & get_d2xyzdeta2() const override
Definition: inf_fe.h:668
std::vector< Real > dzetadx_map
Definition: inf_fe.h:1094
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:1190
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)
Computes the indices of shape functions in the base base_shape and in radial direction radial_shape (...
std::vector< Real > dxidz_map_scaled
Definition: inf_fe.h:1102
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:1154
static Real D_deriv(const Real v)
Definition: inf_fe.h:90
std::vector< Point > xyz
Physical quadrature points.
Definition: inf_fe.h:1045
virtual const std::vector< RealGradient > & get_Sobolev_dweight() const override
Definition: inf_fe.h:790
std::vector< Point > normals
Definition: inf_fe.h:1121
InfFERadial()
Never use an object of this type.
Definition: inf_fe.h:60
std::vector< std::vector< RealGradient > > dphixr_sq
Definition: inf_fe.h:1116
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:1246
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:121
virtual const std::vector< std::vector< OutputShape > > & get_phi_over_decayxR() const override
Definition: inf_fe.h:609
virtual const std::vector< RealGradient > & get_d2xyzdetadzeta() const override
Definition: inf_fe.h:688
virtual const std::vector< Real > & get_detady() const override
Definition: inf_fe.h:732
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:579
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:1164
std::vector< Real > dzetadz_map_scaled
Definition: inf_fe.h:1108
virtual const std::vector< std::vector< OutputGradient > > & get_dphi_over_decayxR() const override
Definition: inf_fe.h:619
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:80
virtual const std::vector< RealGradient > & get_d2xyzdxidzeta() const override
Definition: inf_fe.h:683
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:527
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:1238
virtual const std::vector< Point > & get_normals() const override
Definition: inf_fe.h:806
std::vector< std::vector< Real > > dmodedv
the first local derivative of the radial approximation shapes.
Definition: inf_fe.h:1082
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:553
virtual bool is_hierarchic() const override
Definition: inf_fe.h:442
The constraint matrix storage format.
Definition: dof_map.h:98
std::vector< std::vector< Point > > tangents
Definition: inf_fe.h:1122
The QBase class provides the basic functionality from which various quadrature rules can be derived...
Definition: quadrature.h:53
std::vector< Real > detady_map
Definition: inf_fe.h:1092
virtual const std::vector< Point > & get_xyz() const override
Definition: inf_fe.h:569
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:1106
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:1184
virtual const std::vector< RealGradient > & get_Sobolev_dweightxR_sq() const override
Definition: inf_fe.h:834
virtual const std::vector< Real > & get_curvatures() const override
Definition: inf_fe.h:814
static unsigned int n_shape_functions(const FEType &fet, const Elem *inf_elem)
Definition: inf_fe.h:385
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:110
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:1203
virtual const std::vector< RealGradient > & get_dxyzdzeta() const override
Definition: inf_fe.h:656