libMesh
inf_fe.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2019 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 
20 #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 
76  static Real decay_deriv (const Real) { return -.5; }
77 
82  static Real D (const Real v) { return (1.-v)*(1.-v)/4.; }
83 
87  static Real D_deriv (const Real v) { return (v-1.)/2.; }
88 
93  static Order mapping_order() { return FIRST; }
94 
109  static unsigned int n_dofs (const Order o_radial)
110  { return static_cast<unsigned int>(o_radial)+1; }
111 
121  static unsigned int n_dofs_at_node (const Order o_radial,
122  const unsigned int n_onion);
123 
132  static unsigned int n_dofs_per_elem (const Order o_radial)
133  { return static_cast<unsigned int>(o_radial)+1; }
134 
135 };
136 
137 
138 
148 {
149 private:
150 
155 
156 public:
157 
163  static Elem * build_elem (const Elem * inf_elem);
164 
170  static ElemType get_elem_type (const ElemType type);
171 
177  static unsigned int n_base_mapping_sf (const Elem & base_elem,
178  const Order base_mapping_order);
179 
180 };
181 
182 
183 
215 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map>
216 class InfFE : public FEBase
217 {
218 
219  /*
220  * Protect the nested class
221  */
222 protected:
223 
224 
225 public:
226 
227  // InfFE continued
228 
242  explicit
243  InfFE(const FEType & fet);
244  ~InfFE() {}
245 
246  // The static public members for access from FEInterface etc
247 
260  static Real shape(const FEType & fet,
261  const ElemType t,
262  const unsigned int i,
263  const Point & p);
264 
277  static Real shape(const FEType & fet,
278  const Elem * elem,
279  const unsigned int i,
280  const Point & p);
281 
294  static Real shape_deriv (const FEType & fet,
295  const Elem * inf_elem,
296  const unsigned int i,
297  const unsigned int j,
298  const Point & p);
299 
312  static Real shape_deriv (const FEType & fet,
313  const ElemType inf_elem_type,
314  const unsigned int i,
315  const unsigned int j,
316  const Point & p);
317 
328  static void compute_data(const FEType & fe_t,
329  const Elem * inf_elem,
330  FEComputeData & data);
331 
336  static unsigned int n_shape_functions (const FEType & fet,
337  const ElemType t)
338  { return n_dofs(fet, t); }
339 
346  static unsigned int n_dofs(const FEType & fet,
347  const ElemType inf_elem_type);
348 
353  static unsigned int n_dofs_at_node(const FEType & fet,
354  const ElemType inf_elem_type,
355  const unsigned int n);
356 
361  static unsigned int n_dofs_per_elem(const FEType & fet,
362  const ElemType inf_elem_type);
363 
367  virtual FEContinuity get_continuity() const override
368  { return C_ZERO; } // FIXME - is this true??
369 
374  virtual bool is_hierarchic() const override
375  { return false; } // FIXME - Inf FEs don't handle p elevation yet
376 
384  static void nodal_soln(const FEType & fet,
385  const Elem * elem,
386  const std::vector<Number> & elem_soln,
387  std::vector<Number> & nodal_soln);
388 
389 
390  static Point map (const Elem * inf_elem,
391  const Point & reference_point) {
392  // libmesh_deprecated(); // soon
393  return InfFEMap::map(Dim, inf_elem, reference_point);
394  }
395 
396 
397  static Point inverse_map (const Elem * elem,
398  const Point & p,
399  const Real tolerance = TOLERANCE,
400  const bool secure = true) {
401  // libmesh_deprecated(); // soon
402  return InfFEMap::inverse_map(Dim, elem, p, tolerance, secure);
403  }
404 
405 
406  static void inverse_map (const Elem * elem,
407  const std::vector<Point> & physical_points,
408  std::vector<Point> & reference_points,
409  const Real tolerance = TOLERANCE,
410  const bool secure = true) {
411  // libmesh_deprecated(); // soon
412  return InfFEMap::inverse_map(Dim, elem, physical_points,
413  reference_points, tolerance, secure);
414  }
415 
416 
417  // The workhorses of InfFE. These are often used during matrix assembly.
418 
426  virtual void reinit (const Elem * elem,
427  const std::vector<Point> * const pts = nullptr,
428  const std::vector<Real> * const weights = nullptr) override;
429 
435  virtual void reinit (const Elem * inf_elem,
436  const unsigned int s,
437  const Real tolerance = TOLERANCE,
438  const std::vector<Point> * const pts = nullptr,
439  const std::vector<Real> * const weights = nullptr) override;
440 
446  virtual void edge_reinit (const Elem * elem,
447  const unsigned int edge,
448  const Real tolerance = TOLERANCE,
449  const std::vector<Point> * const pts = nullptr,
450  const std::vector<Real> * const weights = nullptr) override;
451 
456  virtual void side_map (const Elem * /* elem */,
457  const Elem * /* side */,
458  const unsigned int /* s */,
459  const std::vector<Point> & /* reference_side_points */,
460  std::vector<Point> & /* reference_points */) override
461  {
462  libmesh_not_implemented();
463  }
464 
476  virtual void attach_quadrature_rule (QBase * q) override;
477 
482  virtual unsigned int n_shape_functions () const override
483  { return _n_total_approx_sf; }
484 
490  virtual unsigned int n_quadrature_points () const override
492 
493 
494 protected:
495 
496  // static members used by the workhorses
497 
514  static Real eval(Real v,
515  Order o_radial,
516  unsigned int i);
517 
523  static Real eval_deriv(Real v,
524  Order o_radial,
525  unsigned int i);
526 
527 
528 
529  // Non-static members used by the workhorses
530 
535  void update_base_elem (const Elem * inf_elem);
536 
540  virtual void init_base_shape_functions(const std::vector<Point> &,
541  const Elem *) override
542  { libmesh_not_implemented(); }
543 
549  void init_radial_shape_functions(const Elem * inf_elem,
550  const std::vector<Point> * radial_pts = nullptr);
551 
558  void init_shape_functions(const std::vector<Point> & radial_qp,
559  const std::vector<Point> & base_qp,
560  const Elem * inf_elem);
561 
566  void init_face_shape_functions (const std::vector<Point> &,
567  const Elem * inf_side);
568 
577  void combine_base_radial(const Elem * inf_elem);
578 
589  virtual void compute_shape_functions(const Elem *, const std::vector<Point> &) override;
590 
591 
592 
593  // Miscellaneous static members
594 
601  static void compute_node_indices (const ElemType inf_elem_type,
602  const unsigned int outer_node_index,
603  unsigned int & base_node,
604  unsigned int & radial_node);
605 
614  static void compute_node_indices_fast (const ElemType inf_elem_type,
615  const unsigned int outer_node_index,
616  unsigned int & base_node,
617  unsigned int & radial_node);
618 
625  static void compute_shape_indices (const FEType & fet,
626  const ElemType inf_elem_type,
627  const unsigned int i,
628  unsigned int & base_shape,
629  unsigned int & radial_shape);
630 
634  std::vector<Real> dist;
635 
642  std::vector<Real> dweightdv;
643 
651  std::vector<Real> som;
656  std::vector<Real> dsomdv;
657 
662  std::vector<std::vector<Real>> mode;
663 
668  std::vector<std::vector<Real>> dmodedv;
669 
673  std::vector<std::vector<Real>> radial_map;
674 
675 
679  std::vector<std::vector<Real>> dradialdv_map;
680 
686  std::vector<Real> dphasedxi;
687 
693  std::vector<Real> dphasedeta;
694 
700  std::vector<Real> dphasedzeta;
701 
702 
703 
704 
705  // numbering scheme maps
706 
715  std::vector<unsigned int> _radial_node_index;
716 
725  std::vector<unsigned int> _base_node_index;
726 
735  std::vector<unsigned int> _radial_shape_index;
736 
745  std::vector<unsigned int> _base_shape_index;
746 
747 
748 
749 
750  // some more protected members
751 
756  unsigned int _n_total_approx_sf;
757 
762  unsigned int _n_total_qp;
763 
768  std::vector<Real> _total_qrule_weights;
769 
774  std::unique_ptr<QBase> base_qrule;
775 
780  std::unique_ptr<QBase> radial_qrule;
781 
786  std::unique_ptr<Elem> base_elem;
787 
794  std::unique_ptr<FEBase> base_fe;
795 
805 
806 
807 private:
808 
812  virtual bool shapes_need_reinit() const override;
813 
822 
823 
824 #ifdef DEBUG
825 
830  static bool _warned_for_shape;
831  static bool _warned_for_dshape;
832 
833 #endif
834 
840  template <unsigned int friend_Dim, FEFamily friend_T_radial, InfMapType friend_T_map>
841  friend class InfFE;
842 
843  friend class InfFEMap;
844 };
845 
846 
847 
848 // InfFERadial class inline members
849 inline
850 Real InfFERadial::decay(const unsigned int dim, const Real v)
851 {
852  switch (dim)
853  //TODO:[DD] What decay do i have in 2D and 1D?
854  {
855  case 3:
856  return (1.-v)/2.;
857 
858  case 2:
859  return 0.;
860 
861  case 1:
862  return 0.;
863 
864  default:
865  libmesh_error_msg("Invalid dim = " << dim);
866  }
867 }
868 
869 } // namespace libMesh
870 
871 
872 #endif // ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
873 
874 
875 #endif // LIBMESH_INF_FE_H
libMesh::QBase
The QBase class provides the basic functionality from which various quadrature rules can be derived.
Definition: quadrature.h:61
libMesh::InfFE::is_hierarchic
virtual bool is_hierarchic() const override
Definition: inf_fe.h:374
libMesh::InfFERadial::decay
static Real decay(const unsigned int dim, const Real v)
Definition: inf_fe.h:850
libMesh::InfFE::n_quadrature_points
virtual unsigned int n_quadrature_points() const override
Definition: inf_fe.h:490
libMesh::InfFE::base_elem
std::unique_ptr< Elem > base_elem
The base element associated with the current infinite element.
Definition: inf_fe.h:786
libMesh::InfFERadial::mapping_order
static Order mapping_order()
Definition: inf_fe.h:93
libMesh::FEComputeData
class FEComputeData hides arbitrary data to be passed to and from children of FEBase through the FEIn...
Definition: fe_compute_data.h:51
libMesh::InfFEBase::InfFEBase
InfFEBase()
Never use an object of this type.
Definition: inf_fe.h:154
libMesh::InfFE::~InfFE
~InfFE()
Definition: inf_fe.h:244
libMesh::InfFE::update_base_elem
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:101
libMesh::InfFE::InfFE
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:841
libMesh::InfFE::dphasedzeta
std::vector< Real > dphasedzeta
the third local derivative (for 3D, the derivative in radial direction) of the phase term in local co...
Definition: inf_fe.h:700
libMesh::InfFE::current_fe_type
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:804
libMesh::InfFE
A specific instantiation of the FEBase class.
Definition: fe.h:40
libMesh::InfFE::eval_deriv
static Real eval_deriv(Real v, Order o_radial, unsigned int i)
libMesh::InfFE::init_face_shape_functions
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.
Definition: inf_fe_boundary.C:158
libMesh::InfFE::dphasedxi
std::vector< Real > dphasedxi
the first local derivative (for 3D, the first in the base) of the phase term in local coordinates.
Definition: inf_fe.h:686
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::InfFE::n_dofs
static unsigned int n_dofs(const FEType &fet, const ElemType inf_elem_type)
Definition: inf_fe_static.C:58
libMesh::InfFE::nodal_soln
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.
Definition: inf_fe_static.C:123
libMesh::Order
Order
Definition: enum_order.h:40
libMesh::InfFE::combine_base_radial
void combine_base_radial(const Elem *inf_elem)
Combines the shape functions, which were formed in init_shape_functions(Elem *), with geometric data.
Definition: inf_fe.C:751
libMesh::InfFE::_base_node_index
std::vector< unsigned int > _base_node_index
The internal structure of the InfFE – tensor product of base element times radial nodes – has to be d...
Definition: inf_fe.h:725
libMesh::InfFE::radial_qrule
std::unique_ptr< QBase > radial_qrule
The quadrature rule for the base element associated with the current infinite element.
Definition: inf_fe.h:780
libMesh::TOLERANCE
static const Real TOLERANCE
Definition: libmesh_common.h:128
libMesh::InfFE::map
static Point map(const Elem *inf_elem, const Point &reference_point)
Definition: inf_fe.h:390
libMesh::InfFE::mode
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:662
libMesh::InfFE::get_continuity
virtual FEContinuity get_continuity() const override
Definition: inf_fe.h:367
libMesh::InfFE::n_shape_functions
static unsigned int n_shape_functions(const FEType &fet, const ElemType t)
Definition: inf_fe.h:336
libMesh::InfFERadial::D_deriv
static Real D_deriv(const Real v)
Definition: inf_fe.h:87
libMesh::InfFE::dsomdv
std::vector< Real > dsomdv
the first local derivative of the radial decay in local coordinates.
Definition: inf_fe.h:656
libMesh::InfFE::_n_total_qp
unsigned int _n_total_qp
The total number of quadrature points for the current configuration.
Definition: inf_fe.h:762
libMesh::InfFE::init_base_shape_functions
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:540
dim
unsigned int dim
Definition: adaptivity_ex3.C:113
libMesh::InfFE::shapes_need_reinit
virtual bool shapes_need_reinit() const override
Definition: inf_fe.C:1001
libMesh::InfFE::init_shape_functions
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,...
Definition: inf_fe.C:422
libMesh::InfFERadial::InfFERadial
InfFERadial()
Never use an object of this type.
Definition: inf_fe.h:60
libMesh::InfFERadial::decay_deriv
static Real decay_deriv(const Real)
Definition: inf_fe.h:76
libMesh::InfFEMap
Class that encapsulates mapping (i.e.
Definition: inf_fe_map.h:47
libMesh::libmesh_assert
libmesh_assert(ctx)
libMesh::InfFE::inverse_map
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:406
libMesh::InfFEMap::inverse_map
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:88
libMesh::InfFE::compute_shape_functions
virtual void compute_shape_functions(const Elem *, const std::vector< Point > &) override
After having updated the jacobian and the transformation from local to global coordinates in FEAbstra...
Definition: inf_fe.C:899
libMesh::FEBase
FEGenericBase< Real > FEBase
Definition: exact_error_estimator.h:39
libMesh::InfFE::compute_node_indices
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...
Definition: inf_fe_static.C:523
libMesh::InfFE::_warned_for_dshape
static bool _warned_for_dshape
Definition: inf_fe.h:831
libMesh::InfFE::dweightdv
std::vector< Real > dweightdv
the additional radial weight in local coordinates, over all quadrature points.
Definition: inf_fe.h:642
libMesh::InfFERadial::n_dofs_per_elem
static unsigned int n_dofs_per_elem(const Order o_radial)
Definition: inf_fe.h:132
libMesh::InfFE::base_qrule
std::unique_ptr< QBase > base_qrule
The quadrature rule for the base element associated with the current infinite element.
Definition: inf_fe.h:774
libMesh::Point
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:38
libMesh::InfFE::dphasedeta
std::vector< Real > dphasedeta
the second local derivative (for 3D, the second in the base) of the phase term in local coordinates.
Definition: inf_fe.h:693
libMesh::InfFE::reinit
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:112
libMesh::C_ZERO
Definition: enum_fe_family.h:79
libMesh::InfFE::som
std::vector< Real > som
the radial decay in local coordinates.
Definition: inf_fe.h:651
libMesh::InfFERadial
Infinite elements are in some sense directional, compared to conventional finite elements.
Definition: inf_fe.h:53
libMesh::InfFE::n_shape_functions
virtual unsigned int n_shape_functions() const override
Definition: inf_fe.h:482
libMesh::InfFE::_warned_for_nodal_soln
static bool _warned_for_nodal_soln
static members that are used to issue warning messages only once.
Definition: inf_fe.h:829
libMesh::InfFE::side_map
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:456
libMesh::InfFERadial::n_dofs
static unsigned int n_dofs(const Order o_radial)
Definition: inf_fe.h:109
libMesh::InfFE::base_fe
std::unique_ptr< FEBase > base_fe
Have a FE<Dim-1,T_base> handy for base approximation.
Definition: inf_fe.h:794
libMesh::InfFE::_base_shape_index
std::vector< unsigned int > _base_shape_index
The internal structure of the InfFE – tensor product of base element shapes times radial shapes – has...
Definition: inf_fe.h:745
libMesh::InfFE::inverse_map
static Point inverse_map(const Elem *elem, const Point &p, const Real tolerance=TOLERANCE, const bool secure=true)
Definition: inf_fe.h:397
libMesh::InfFE::shape_deriv
static Real shape_deriv(const FEType &fet, const Elem *inf_elem, const unsigned int i, const unsigned int j, const Point &p)
Definition: inf_fe_static.C:284
libMesh::InfFEBase
This nested class contains most of the static methods related to the base part of an infinite element...
Definition: inf_fe.h:147
libMesh::InfFEMap::map
static Point map(const unsigned int dim, const Elem *inf_elem, const Point &reference_point)
Definition: inf_fe_map.C:40
libMesh::FEType
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:178
libMesh::InfFE::attach_quadrature_rule
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:71
libMesh::InfFEBase::build_elem
static Elem * build_elem(const Elem *inf_elem)
Build the base element of an infinite element.
Definition: inf_fe_base_radial.C:35
libMesh::InfFE::_n_total_approx_sf
unsigned int _n_total_approx_sf
The number of total approximation shape functions for the current configuration.
Definition: inf_fe.h:756
libMesh::InfFE::shape
static Real shape(const FEType &fet, const ElemType t, const unsigned int i, const Point &p)
Definition: inf_fe_static.C:158
libMesh::InfFE::radial_map
std::vector< std::vector< Real > > radial_map
the radial mapping shapes in local coordinates
Definition: inf_fe.h:673
libMesh::InfFERadial::n_dofs_at_node
static unsigned int n_dofs_at_node(const Order o_radial, const unsigned int n_onion)
Definition: inf_fe_base_radial.C:119
libMesh::InfFE::_total_qrule_weights
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:768
libMesh::InfFE::_compute_node_indices_fast_current_elem_type
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:821
libMesh::InfFE::compute_shape_indices
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 (...
Definition: inf_fe_static.C:839
libMesh::Elem
This is the base class from which all geometric element types are derived.
Definition: elem.h:100
libMesh::InfFE::dmodedv
std::vector< std::vector< Real > > dmodedv
the first local derivative of the radial approximation shapes.
Definition: inf_fe.h:668
libMesh::InfFE::_radial_node_index
std::vector< unsigned int > _radial_node_index
The internal structure of the InfFE – tensor product of base element times radial nodes – has to be d...
Definition: inf_fe.h:715
data
IterBase * data
Ideally this private member data should have protected access.
Definition: variant_filter_iterator.h:337
libMesh::InfFE::eval
static Real eval(Real v, Order o_radial, unsigned int i)
libMesh::InfFE::n_dofs_at_node
static unsigned int n_dofs_at_node(const FEType &fet, const ElemType inf_elem_type, const unsigned int n)
Definition: inf_fe_static.C:76
libMesh::InfFE::init_radial_shape_functions
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:336
libMesh::InfFE::compute_data
static void compute_data(const FEType &fe_t, const Elem *inf_elem, FEComputeData &data)
Generalized version of shape(), takes an Elem *.
Definition: inf_fe_static.C:337
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::InfFE::dradialdv_map
std::vector< std::vector< Real > > dradialdv_map
the first local derivative of the radial mapping shapes
Definition: inf_fe.h:679
libMesh::InfFE::compute_node_indices_fast
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.
Definition: inf_fe_static.C:730
libMesh::InfFE::_warned_for_shape
static bool _warned_for_shape
Definition: inf_fe.h:830
libMesh::FEContinuity
FEContinuity
Definition: enum_fe_family.h:77
libMesh::FIRST
Definition: enum_order.h:42
libMesh::InfFE::dist
std::vector< Real > dist
the radial distance of the base nodes from the origin
Definition: inf_fe.h:634
libMesh::InfFE::edge_reinit
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.
Definition: inf_fe_boundary.C:140
libMesh::InfFE::_radial_shape_index
std::vector< unsigned int > _radial_shape_index
The internal structure of the InfFE – tensor product of base element shapes times radial shapes – has...
Definition: inf_fe.h:735
libMesh::InfFERadial::D
static Real D(const Real v)
Definition: inf_fe.h:82
libMesh::InfFEBase::get_elem_type
static ElemType get_elem_type(const ElemType type)
Definition: inf_fe_base_radial.C:49
libMesh::ElemType
ElemType
Defines an enum for geometric element types.
Definition: enum_elem_type.h:33
libMesh::InfFEBase::n_base_mapping_sf
static unsigned int n_base_mapping_sf(const Elem &base_elem, const Order base_mapping_order)
Definition: inf_fe_base_radial.C:95
libMesh::InfFE::n_dofs_per_elem
static unsigned int n_dofs_per_elem(const FEType &fet, const ElemType inf_elem_type)
Definition: inf_fe_static.C:105