Go to the documentation of this file.
20 #ifndef LIBMESH_INF_FE_H
21 #define LIBMESH_INF_FE_H
23 #include "libmesh/libmesh_config.h"
25 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
28 #include "libmesh/fe_base.h"
29 #include "libmesh/inf_fe_map.h"
82 static Real D (
const Real v) {
return (1.-v)*(1.-v)/4.; }
110 {
return static_cast<unsigned int>(o_radial)+1; }
122 const unsigned int n_onion);
133 {
return static_cast<unsigned int>(o_radial)+1; }
178 const Order base_mapping_order);
215 template <
unsigned int Dim, FEFamily T_radial, InfMapType T_map>
216 class InfFE :
public FEBase
243 InfFE(
const FEType & fet);
262 const unsigned int i,
279 const unsigned int i,
295 const Elem * inf_elem,
296 const unsigned int i,
297 const unsigned int j,
314 const unsigned int i,
315 const unsigned int j,
329 const Elem * inf_elem,
338 {
return n_dofs(fet, t); }
355 const unsigned int n);
386 const std::vector<Number> & elem_soln,
391 const Point & reference_point) {
400 const bool secure =
true) {
407 const std::vector<Point> & physical_points,
408 std::vector<Point> & reference_points,
410 const bool secure =
true) {
413 reference_points, tolerance, secure);
427 const std::vector<Point> *
const pts =
nullptr,
428 const std::vector<Real> *
const weights =
nullptr)
override;
436 const unsigned int s,
438 const std::vector<Point> *
const pts =
nullptr,
439 const std::vector<Real> *
const weights =
nullptr)
override;
447 const unsigned int edge,
449 const std::vector<Point> *
const pts =
nullptr,
450 const std::vector<Real> *
const weights =
nullptr)
override;
459 const std::vector<Point> & ,
460 std::vector<Point> & )
override
462 libmesh_not_implemented();
541 const Elem *)
override
542 { libmesh_not_implemented(); }
550 const std::vector<Point> * radial_pts =
nullptr);
559 const std::vector<Point> & base_qp,
560 const Elem * inf_elem);
567 const Elem * inf_side);
602 const unsigned int outer_node_index,
603 unsigned int & base_node,
604 unsigned int & radial_node);
615 const unsigned int outer_node_index,
616 unsigned int & base_node,
617 unsigned int & radial_node);
627 const unsigned int i,
628 unsigned int & base_shape,
629 unsigned int & radial_shape);
662 std::vector<std::vector<Real>>
mode;
840 template <
unsigned int friend_Dim, FEFamily friend_T_radial, InfMapType friend_T_map>
865 libmesh_error_msg(
"Invalid dim = " <<
dim);
872 #endif // ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
875 #endif // LIBMESH_INF_FE_H
The QBase class provides the basic functionality from which various quadrature rules can be derived.
virtual bool is_hierarchic() const override
static Real decay(const unsigned int dim, const Real v)
virtual unsigned int n_quadrature_points() const override
std::unique_ptr< Elem > base_elem
The base element associated with the current infinite element.
static Order mapping_order()
class FEComputeData hides arbitrary data to be passed to and from children of FEBase through the FEIn...
InfFEBase()
Never use an object of this type.
void update_base_elem(const Elem *inf_elem)
Updates the protected member base_elem to the appropriate base element for the given inf_elem.
friend class InfFE
Make all InfFE<Dim,T_radial,T_map> classes friends of each other, so that the protected eval() may be...
std::vector< Real > dphasedzeta
the third local derivative (for 3D, the derivative in radial direction) of the phase term in local co...
FEType current_fe_type
This FEType stores the characteristics for which the data structures phi, phi_map etc are currently i...
A specific instantiation of the FEBase class.
static Real eval_deriv(Real v, Order o_radial, unsigned int i)
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.
std::vector< Real > dphasedxi
the first local derivative (for 3D, the first in the base) of the phase term in local coordinates.
The libMesh namespace provides an interface to certain functionality in the library.
static unsigned int n_dofs(const FEType &fet, const ElemType inf_elem_type)
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.
void combine_base_radial(const Elem *inf_elem)
Combines the shape functions, which were formed in init_shape_functions(Elem *), with geometric data.
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...
std::unique_ptr< QBase > radial_qrule
The quadrature rule for the base element associated with the current infinite element.
static const Real TOLERANCE
static Point map(const Elem *inf_elem, const Point &reference_point)
std::vector< std::vector< Real > > mode
the radial approximation shapes in local coordinates Needed when setting up the overall shape functio...
virtual FEContinuity get_continuity() const override
static unsigned int n_shape_functions(const FEType &fet, const ElemType t)
static Real D_deriv(const Real v)
std::vector< Real > dsomdv
the first local derivative of the radial decay in local coordinates.
unsigned int _n_total_qp
The total number of quadrature points for the current configuration.
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>.
virtual bool shapes_need_reinit() const override
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,...
InfFERadial()
Never use an object of this type.
static Real decay_deriv(const Real)
Class that encapsulates mapping (i.e.
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)
static Point inverse_map(const unsigned int dim, const Elem *elem, const Point &p, const Real tolerance=TOLERANCE, const bool secure=true)
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...
FEGenericBase< Real > FEBase
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 bool _warned_for_dshape
std::vector< Real > dweightdv
the additional radial weight in local coordinates, over all quadrature points.
static unsigned int n_dofs_per_elem(const Order o_radial)
std::unique_ptr< QBase > base_qrule
The quadrature rule for the base element associated with the current infinite element.
A Point defines a location in LIBMESH_DIM dimensional Real space.
std::vector< Real > dphasedeta
the second local derivative (for 3D, the second in the base) of the phase term in local coordinates.
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.
std::vector< Real > som
the radial decay in local coordinates.
Infinite elements are in some sense directional, compared to conventional finite elements.
virtual unsigned int n_shape_functions() const override
static bool _warned_for_nodal_soln
static members that are used to issue warning messages only once.
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...
static unsigned int n_dofs(const Order o_radial)
std::unique_ptr< FEBase > base_fe
Have a FE<Dim-1,T_base> handy for base approximation.
std::vector< unsigned int > _base_shape_index
The internal structure of the InfFE – tensor product of base element shapes times radial shapes – has...
static Point inverse_map(const Elem *elem, const Point &p, const Real tolerance=TOLERANCE, const bool secure=true)
static Real shape_deriv(const FEType &fet, const Elem *inf_elem, const unsigned int i, const unsigned int j, const Point &p)
This nested class contains most of the static methods related to the base part of an infinite element...
static Point map(const unsigned int dim, const Elem *inf_elem, const Point &reference_point)
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
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...
static Elem * build_elem(const Elem *inf_elem)
Build the base element of an infinite element.
unsigned int _n_total_approx_sf
The number of total approximation shape functions for the current configuration.
static Real shape(const FEType &fet, const ElemType t, const unsigned int i, const Point &p)
std::vector< std::vector< Real > > radial_map
the radial mapping shapes in local coordinates
static unsigned int n_dofs_at_node(const Order o_radial, const unsigned int n_onion)
std::vector< Real > _total_qrule_weights
this vector contains the combined integration weights, so that FEAbstract::compute_map() can still be...
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...
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 (...
This is the base class from which all geometric element types are derived.
std::vector< std::vector< Real > > dmodedv
the first local derivative of the radial approximation shapes.
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...
IterBase * data
Ideally this private member data should have protected access.
static Real eval(Real v, Order o_radial, unsigned int i)
static unsigned int n_dofs_at_node(const FEType &fet, const ElemType inf_elem_type, const unsigned int n)
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.
static void compute_data(const FEType &fe_t, const Elem *inf_elem, FEComputeData &data)
Generalized version of shape(), takes an Elem *.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
std::vector< std::vector< Real > > dradialdv_map
the first local derivative of the radial mapping shapes
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.
static bool _warned_for_shape
std::vector< Real > dist
the radial distance of the base nodes from the origin
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::vector< unsigned int > _radial_shape_index
The internal structure of the InfFE – tensor product of base element shapes times radial shapes – has...
static Real D(const Real v)
static ElemType get_elem_type(const ElemType type)
ElemType
Defines an enum for geometric element types.
static unsigned int n_base_mapping_sf(const Elem &base_elem, const Order base_mapping_order)
static unsigned int n_dofs_per_elem(const FEType &fet, const ElemType inf_elem_type)