Go to the documentation of this file.
20 #ifndef LIBMESH_FE_ABSTRACT_H
21 #define LIBMESH_FE_ABSTRACT_H
24 #include "libmesh/reference_counted_object.h"
25 #include "libmesh/point.h"
26 #include "libmesh/vector_value.h"
27 #include "libmesh/fe_type.h"
28 #include "libmesh/fe_map.h"
29 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
30 #include "libmesh/tensor_value.h"
33 #ifdef LIBMESH_FORWARD_DECLARE_ENUMS
39 #include "libmesh/enum_elem_type.h"
52 template <
typename T>
class DenseMatrix;
53 template <
typename T>
class DenseVector;
59 template <
typename T>
class NumericVector;
62 #ifdef LIBMESH_ENABLE_NODE_CONSTRAINTS
63 class NodeConstraints;
66 #ifdef LIBMESH_ENABLE_PERIODIC
67 class PeriodicBoundaries;
68 class PointLocatorBase;
71 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
72 template <
unsigned int Dim, FEFamily T_radial, InfMapType T_map>
125 static std::unique_ptr<FEAbstract>
build (
const unsigned int dim,
143 const std::vector<Point> *
const pts =
nullptr,
144 const std::vector<Real> *
const weights =
nullptr) = 0;
155 const unsigned int side,
157 const std::vector<Point> *
const pts =
nullptr,
158 const std::vector<Real> *
const weights =
nullptr) = 0;
169 const unsigned int edge,
171 const std::vector<Point> * pts =
nullptr,
172 const std::vector<Real> * weights =
nullptr) = 0;
180 const unsigned int s,
181 const std::vector<Point> & reference_side_points,
182 std::vector<Point> & reference_points) = 0;
199 std::vector<Point> & nodes);
202 #ifdef LIBMESH_ENABLE_NODE_CONSTRAINTS
209 #endif // LIBMESH_ENABLE_NODE_CONSTRAINTS
211 #ifdef LIBMESH_ENABLE_PERIODIC
213 #ifdef LIBMESH_ENABLE_NODE_CONSTRAINTS
223 #endif // LIBMESH_ENABLE_NODE_CONSTRAINTS
225 #endif // LIBMESH_ENABLE_PERIODIC
266 {
return _fe_map->get_dxyzdzeta(); }
268 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
383 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
471 virtual void print_phi(std::ostream & os)
const =0;
478 virtual void print_dphi(std::ostream & os)
const =0;
480 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
487 virtual void print_d2phi(std::ostream & os)
const =0;
553 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
622 #endif // LIBMESH_FE_ABSTRACT_H
The QBase class provides the basic functionality from which various quadrature rules can be derived.
virtual unsigned int n_shape_functions() const =0
We're using a class instead of a typedef to allow forward declarations and future flexibility.
const std::vector< RealGradient > & get_d2xyzdxideta() const
virtual void side_map(const Elem *elem, const Elem *side, const unsigned int s, const std::vector< Point > &reference_side_points, std::vector< Point > &reference_points)=0
Computes the reference space quadrature points on the side of an element based on the side quadrature...
const std::vector< Point > & get_xyz() const
bool calculate_dphi
Should we calculate shape function gradients?
FEFamily family
The type of finite element.
virtual bool shapes_need_reinit() const =0
void print_xyz(std::ostream &os) const
Prints the spatial location of each quadrature point (on the physical element).
const std::vector< Point > & get_normals() const
FEType fe_type
The finite element type for this object.
This class implements reference counting.
The libMesh namespace provides an interface to certain functionality in the library.
const std::vector< std::vector< Point > > & get_tangents() const
const std::vector< RealGradient > & get_d2xyzdzeta2() const
const std::vector< Real > & get_detadx() const
This class forms the foundation from which generic finite elements may be derived.
virtual void print_dphi(std::ostream &os) const =0
Prints the value of each shape function's derivative at each quadrature point.
bool calculate_dphiref
Should we calculate reference shape function gradients?
const std::vector< Real > & get_JxW() const
bool calculations_started
Have calculations with this object already been started? Then all get_* functions should already have...
const std::vector< Real > & get_dzetadx() const
static const Real TOLERANCE
const std::vector< Real > & get_dxidy() const
bool calculate_phi
Should we calculate shape functions?
bool calculate_map
Are we calculating mapping functions?
FEFamily get_family() const
bool calculate_curl_phi
Should we calculate shape function curls?
const std::vector< RealGradient > & get_d2xyzdeta2() const
const std::vector< RealGradient > & get_dxyzdzeta() const
bool shapes_on_quadrature
A flag indicating if current data structures correspond to quadrature rule points.
void set_fe_order(int new_order)
Sets the base FE order of the finite element.
friend std::ostream & operator<<(std::ostream &os, const FEAbstract &fe)
Same as above, but allows you to print to a stream.
const std::vector< Real > & get_detady() const
const FEMap & get_fe_map() const
const std::vector< RealGradient > & get_d2xyzdetadzeta() const
virtual FEContinuity get_continuity() const =0
void print_info(std::ostream &os) const
Prints all the relevant information about the current element.
This is the MeshBase class.
unsigned int get_p_level() const
unsigned int get_dim() const
virtual void print_phi(std::ostream &os) const =0
Prints the value of each shape function at each quadrature point.
bool calculate_div_phi
Should we calculate shape function divergences?
static std::unique_ptr< FEAbstract > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type.
const std::vector< Real > & get_dzetady() const
QBase * qrule
A pointer to the quadrature rule employed.
ElemType get_type() const
A Point defines a location in LIBMESH_DIM dimensional Real space.
std::unique_ptr< FEMap > _fe_map
virtual void edge_reinit(const Elem *elem, const unsigned int edge, const Real tolerance=TOLERANCE, const std::vector< Point > *pts=nullptr, const std::vector< Real > *weights=nullptr)=0
Reinitializes all the physical element-dependent data based on the edge of the element elem.
const std::vector< Real > & get_curvatures() const
const std::vector< Real > & get_dxidz() const
virtual unsigned int n_quadrature_points() const =0
const std::vector< Real > & get_dzetadz() const
const unsigned int dim
The dimensionality of the object.
const std::vector< Real > & get_dxidx() const
static void compute_node_constraints(NodeConstraints &constraints, const Elem *elem)
Computes the nodal constraint contributions (for non-conforming adapted meshes), using Lagrange geome...
FEType get_fe_type() const
const std::vector< Real > & get_detadz() const
static void compute_periodic_node_constraints(NodeConstraints &constraints, const PeriodicBoundaries &boundaries, const MeshBase &mesh, const PointLocatorBase *point_locator, const Elem *elem)
Computes the node position constraint equation contributions (for meshes with periodic boundary condi...
const std::vector< RealGradient > & get_dxyzdeta() const
static bool on_reference_element(const Point &p, const ElemType t, const Real eps=TOLERANCE)
The Node constraint storage format.
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
OrderWrapper order
The approximation order of the element.
virtual ~FEAbstract()
Destructor.
unsigned int _p_level
The p refinement level the current data structures are set up for.
This is the base class from which all geometric element types are derived.
virtual void print_d2phi(std::ostream &os) const =0
Prints the value of each shape function's second derivatives at each quadrature point.
ElemType elem_type
The element type the current data structures are set up for.
static void get_refspace_nodes(const ElemType t, std::vector< Point > &nodes)
virtual void compute_shape_functions(const Elem *, const std::vector< Point > &)=0
After having updated the jacobian and the transformation from local to global coordinates in FEMap::c...
Class contained in FE that encapsulates mapping (i.e.
const std::vector< RealGradient > & get_dxyzdxi() const
const std::vector< RealGradient > & get_d2xyzdxi2() const
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual void reinit(const Elem *elem, const std::vector< Point > *const pts=nullptr, const std::vector< Real > *const weights=nullptr)=0
This is at the core of this class.
virtual void attach_quadrature_rule(QBase *q)=0
Provides the class with the quadrature rule.
FEAbstract(const unsigned int dim, const FEType &fet)
Constructor.
This is the base class for point locators.
virtual bool is_hierarchic() const =0
void print_JxW(std::ostream &os) const
Prints the Jacobian times the weight for each quadrature point.
bool calculate_d2phi
Should we calculate shape function hessians?
void ErrorVector unsigned int
ElemType
Defines an enum for geometric element types.
const std::vector< RealGradient > & get_d2xyzdxidzeta() const