Go to the documentation of this file.
20 #ifndef LIBMESH_FE_INTERFACE_H
21 #define LIBMESH_FE_INTERFACE_H
24 #include "libmesh/libmesh_common.h"
25 #include "libmesh/vector_value.h"
50 #ifdef LIBMESH_ENABLE_PERIODIC
51 class PeriodicBoundaries;
52 class PointLocatorBase;
99 static unsigned int n_dofs(
const unsigned int dim,
110 static unsigned int n_dofs(
const unsigned int dim,
124 const unsigned int n);
157 const unsigned int dim,
160 std::vector<unsigned int> & di);
170 const unsigned int dim,
173 std::vector<unsigned int> & di);
189 const std::vector<Number> & elem_soln,
208 const bool secure =
true);
216 const std::vector<Point> & physical_points,
217 std::vector<Point> & reference_points,
219 const bool secure =
true);
244 const unsigned int i,
259 const unsigned int i,
271 template<
typename OutputType>
272 static void shape(
const unsigned int dim,
275 const unsigned int i,
288 template<
typename OutputType>
289 static void shape(
const unsigned int dim,
292 const unsigned int i,
298 const unsigned int i,
300 const bool add_p_level);
324 const unsigned int i,
325 const unsigned int j,
341 const unsigned int i,
342 const unsigned int j,
347 const unsigned int i,
348 const unsigned int j,
350 const bool add_p_level);
360 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
383 const unsigned int i,
384 const unsigned int j,
405 const unsigned int i,
406 const unsigned int j,
411 const unsigned int i,
412 const unsigned int j,
414 const bool add_p_level);
440 #ifdef LIBMESH_ENABLE_AMR
448 const unsigned int variable_number,
450 #endif // #ifdef LIBMESH_ENABLE_AMR
452 #ifdef LIBMESH_ENABLE_PERIODIC
463 const unsigned int variable_number,
465 #endif // #ifdef LIBMESH_ENABLE_PERIODIC
517 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
536 const unsigned int n);
545 const std::vector<Number> & elem_soln,
558 const bool secure =
true);
563 const std::vector<Point> & physical_points,
564 std::vector<Point> & reference_points,
566 const bool secure =
true);
576 const unsigned int i,
582 const unsigned int i,
589 const unsigned int i,
590 const unsigned int j,
596 const unsigned int i,
597 const unsigned int j,
612 #endif // LIBMESH_FE_INTERFACE_H
The constraint matrix storage format.
static shape_ptr shape_function(const unsigned int dim, const FEType &fe_t)
We're using a class instead of a typedef to allow forward declarations and future flexibility.
static bool on_reference_element(const Point &p, const ElemType t, const Real eps=TOLERANCE)
static bool extra_hanging_dofs(const FEType &fe_t)
class FEComputeData hides arbitrary data to be passed to and from children of FEBase through the FEIn...
static unsigned int n_shape_functions(const unsigned int dim, const FEType &fe_t, const ElemType t)
static FEContinuity get_continuity(const FEType &fe_type)
Returns the input FEType's FEContinuity based on the underlying FEFamily and potentially the Order,...
static unsigned int n_vec_dim(const MeshBase &mesh, const FEType &fe_type)
static Point inverse_map(const unsigned int dim, const FEType &fe_t, const Elem *elem, const Point &p, const Real tolerance=TOLERANCE, const bool secure=true)
This is now deprecated; use FEMap::inverse_map instead.
static void ifem_nodal_soln(const unsigned int dim, const FEType &fe_t, const Elem *elem, const std::vector< Number > &elem_soln, std::vector< Number > &nodal_soln)
This class provides an encapsulated access to all static public member functions of finite element cl...
static n_dofs_at_node_ptr n_dofs_at_node_function(const unsigned int dim, const FEType &fe_t)
static unsigned int ifem_n_dofs(const unsigned int dim, const FEType &fe_t, const ElemType t)
static unsigned int n_dofs_at_node(const unsigned int dim, const FEType &fe_t, const ElemType t, const unsigned int n)
static Real shape_second_deriv(const unsigned int dim, const FEType &fe_t, const ElemType t, const unsigned int i, const unsigned int j, const Point &p)
static unsigned int ifem_n_dofs_per_elem(const unsigned int dim, const FEType &fe_t, const ElemType t)
The libMesh namespace provides an interface to certain functionality in the library.
virtual ~FEInterface()
Destructor.
static shape_deriv_ptr shape_deriv_function(const unsigned int dim, const FEType &fe_t)
static const Real TOLERANCE
static bool is_InfFE_elem(const ElemType et)
Real(* shape_ptr)(const Elem *elem, const Order o, const unsigned int i, const Point &p, const bool add_p_level)
unsigned int(* n_dofs_at_node_ptr)(const ElemType, const Order, const unsigned int)
static void nodal_soln(const unsigned int dim, const FEType &fe_t, const Elem *elem, const std::vector< Number > &elem_soln, std::vector< Number > &nodal_soln)
Build the nodal soln from the element soln.
static void compute_periodic_constraints(DofConstraints &constraints, DofMap &dof_map, const PeriodicBoundaries &boundaries, const MeshBase &mesh, const PointLocatorBase *point_locator, const unsigned int variable_number, const Elem *elem)
Computes the constraint matrix contributions (for periodic boundary conditions) corresponding to vari...
static Real ifem_shape_deriv(const unsigned int dim, const FEType &fe_t, const ElemType t, const unsigned int i, const unsigned int j, const Point &p)
This is the MeshBase class.
static unsigned int n_dofs_per_elem(const unsigned int dim, const FEType &fe_t, const ElemType t)
A Point defines a location in LIBMESH_DIM dimensional Real space.
static void dofs_on_edge(const Elem *const elem, const unsigned int dim, const FEType &fe_t, unsigned int e, std::vector< unsigned int > &di)
Fills the vector di with the local degree of freedom indices associated with edge e of element elem A...
static unsigned int ifem_n_dofs_at_node(const unsigned int dim, const FEType &fe_t, const ElemType t, const unsigned int n)
static unsigned int n_dofs(const unsigned int dim, const FEType &fe_t, const ElemType t)
Real(* shape_deriv_ptr)(const Elem *elem, const Order o, const unsigned int i, const unsigned int j, const Point &p, const bool add_p_level)
static unsigned int ifem_n_shape_functions(const unsigned int dim, const FEType &fe_t, const ElemType t)
static Point map(unsigned int dim, const FEType &fe_t, const Elem *elem, const Point &p)
This is now deprecated; use FEMap::map instead.
static Real ifem_shape(const unsigned int dim, const FEType &fe_t, const ElemType t, const unsigned int i, const Point &p)
static void dofs_on_side(const Elem *const elem, const unsigned int dim, const FEType &fe_t, unsigned int s, std::vector< unsigned int > &di)
Fills the vector di with the local degree of freedom indices associated with side s of element elem A...
static Real shape_deriv(const unsigned int dim, const FEType &fe_t, const ElemType t, const unsigned int i, const unsigned int j, const Point &p)
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
static bool ifem_on_reference_element(const Point &p, const ElemType t, const Real eps)
This class handles the numbering of degrees of freedom on a mesh.
static FEFieldType field_type(const FEType &fe_type)
Real(* shape_second_deriv_ptr)(const Elem *elem, const Order o, const unsigned int i, const unsigned int j, const Point &p, const bool add_p_level)
This is the base class from which all geometric element types are derived.
static void compute_constraints(DofConstraints &constraints, DofMap &dof_map, const unsigned int variable_number, const Elem *elem)
Computes the constraint matrix contributions (for non-conforming adapted meshes) corresponding to var...
IterBase * data
Ideally this private member data should have protected access.
static Point ifem_inverse_map(const unsigned int dim, const FEType &fe_t, const Elem *elem, const Point &p, const Real tolerance=TOLERANCE, const bool secure=true)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
static shape_second_deriv_ptr shape_second_deriv_function(const unsigned int dim, const FEType &fe_t)
This is the base class for point locators.
static void ifem_compute_data(const unsigned int dim, const FEType &fe_t, const Elem *elem, FEComputeData &data)
static unsigned int max_order(const FEType &fe_t, const ElemType &el_t)
void ErrorVector unsigned int
ElemType
Defines an enum for geometric element types.
FEInterface()
Empty constructor.
static void compute_data(const unsigned int dim, const FEType &fe_t, const Elem *elem, FEComputeData &data)
Lets the appropriate child of FEBase compute the requested data for the input specified in data,...
static Real shape(const unsigned int dim, const FEType &fe_t, const ElemType t, const unsigned int i, const Point &p)
static Point ifem_map(const unsigned int dim, const FEType &fe_t, const Elem *elem, const Point &p)