libMesh
Public Types | Public Member Functions | Static Public Member Functions | Private Member Functions | Static Private Member Functions | List of all members
libMesh::FEInterface Class Reference

This class provides an encapsulated access to all static public member functions of finite element classes. More...

#include <fe_interface.h>

Public Types

typedef unsigned int(* n_dofs_at_node_ptr) (const ElemType, const Order, const unsigned int)
 

Public Member Functions

virtual ~FEInterface ()
 Destructor. More...
 

Static Public Member Functions

static unsigned int n_shape_functions (const unsigned int dim, const FEType &fe_t, const ElemType t)
 
static unsigned int 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 n_dofs_at_node_ptr n_dofs_at_node_function (const unsigned int dim, const FEType &fe_t)
 
static unsigned int n_dofs_per_elem (const unsigned int dim, const FEType &fe_t, const ElemType t)
 
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 Automatically decides which finite element class to use. More...
 
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 Automatically decides which finite element class to use. More...
 
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. More...
 
static Point map (unsigned int dim, const FEType &fe_t, const Elem *elem, const Point &p)
 
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)
 
static void inverse_map (const unsigned int dim, const FEType &fe_t, const Elem *elem, const std::vector< Point > &physical_points, std::vector< Point > &reference_points, const Real tolerance=TOLERANCE, const bool secure=true)
 
static bool on_reference_element (const Point &p, const ElemType t, const Real eps=TOLERANCE)
 
static Real shape (const unsigned int dim, const FEType &fe_t, const ElemType t, const unsigned int i, const Point &p)
 
static Real shape (const unsigned int dim, const FEType &fe_t, const Elem *elem, const unsigned int i, const Point &p)
 
template<typename OutputType >
static void shape (const unsigned int dim, const FEType &fe_t, const ElemType t, const unsigned int i, const Point &p, OutputType &phi)
 
template<typename OutputType >
static void shape (const unsigned int dim, const FEType &fe_t, const Elem *elem, const unsigned int i, const Point &p, OutputType &phi)
 
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, and sets the values in data. More...
 
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 variable number var_number. More...
 
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 variable number var_number. More...
 
static unsigned int max_order (const FEType &fe_t, const ElemType &el_t)
 
static bool extra_hanging_dofs (const FEType &fe_t)
 
static FEFieldType field_type (const FEType &fe_type)
 
static FEFieldType field_type (const FEFamily &fe_family)
 
static unsigned int n_vec_dim (const MeshBase &mesh, const FEType &fe_type)
 

Private Member Functions

 FEInterface ()
 Empty constructor. More...
 

Static Private Member Functions

static bool is_InfFE_elem (const ElemType et)
 
static unsigned int ifem_n_shape_functions (const unsigned int dim, const FEType &fe_t, const ElemType t)
 
static unsigned int ifem_n_dofs (const unsigned int dim, const FEType &fe_t, const ElemType t)
 
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 ifem_n_dofs_per_elem (const unsigned int dim, const FEType &fe_t, const ElemType t)
 
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)
 
static Point ifem_map (const unsigned int dim, const FEType &fe_t, const Elem *elem, const Point &p)
 
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)
 
static void ifem_inverse_map (const unsigned int dim, const FEType &fe_t, const Elem *elem, const std::vector< Point > &physical_points, std::vector< Point > &reference_points, const Real tolerance=TOLERANCE, const bool secure=true)
 
static bool ifem_on_reference_element (const Point &p, const ElemType t, const Real eps)
 
static Real ifem_shape (const unsigned int dim, const FEType &fe_t, const ElemType t, const unsigned int i, const Point &p)
 
static Real ifem_shape (const unsigned int dim, const FEType &fe_t, const Elem *elem, const unsigned int i, const Point &p)
 
static void ifem_compute_data (const unsigned int dim, const FEType &fe_t, const Elem *elem, FEComputeData &data)
 

Detailed Description

This class provides an encapsulated access to all static public member functions of finite element classes.

Using this class, one need not worry about the correct finite element class.

Author
Daniel Dreyer
Date
2002-2007 Interface class which provides access to FE functions.

Definition at line 64 of file fe_interface.h.

Member Typedef Documentation

◆ n_dofs_at_node_ptr

typedef unsigned int(* libMesh::FEInterface::n_dofs_at_node_ptr) (const ElemType, const Order, const unsigned int)

Definition at line 114 of file fe_interface.h.

Constructor & Destructor Documentation

◆ FEInterface()

libMesh::FEInterface::FEInterface ( )
private

Empty constructor.

Do not create an object of this type.

◆ ~FEInterface()

virtual libMesh::FEInterface::~FEInterface ( )
virtual

Destructor.

Definition at line 78 of file fe_interface.h.

78 {}

Member Function Documentation

◆ compute_constraints()

static void libMesh::FEInterface::compute_constraints ( DofConstraints constraints,
DofMap dof_map,
const unsigned int  variable_number,
const Elem elem 
)
static

Computes the constraint matrix contributions (for non-conforming adapted meshes) corresponding to variable number var_number.

◆ compute_data()

static void libMesh::FEInterface::compute_data ( const unsigned int  dim,
const FEType fe_t,
const Elem elem,
FEComputeData data 
)
static

Lets the appropriate child of FEBase compute the requested data for the input specified in data, and sets the values in data.

See this as a generalization of shape(). With infinite elements disabled, computes values for all shape functions of elem evaluated at p.

Note
On a p-refined element, fe_t.order should be the base order of the element.

◆ compute_periodic_constraints()

static void libMesh::FEInterface::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 
)
static

Computes the constraint matrix contributions (for periodic boundary conditions) corresponding to variable number var_number.

◆ dofs_on_edge()

static void libMesh::FEInterface::dofs_on_edge ( const Elem *const  elem,
const unsigned int  dim,
const FEType fe_t,
unsigned int  e,
std::vector< unsigned int > &  di 
)
static

Fills the vector di with the local degree of freedom indices associated with edge e of element elem Automatically decides which finite element class to use.

On a p-refined element, fe_t.order should be the base order of the element.

Referenced by libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::operator()().

◆ dofs_on_side()

static void libMesh::FEInterface::dofs_on_side ( const Elem *const  elem,
const unsigned int  dim,
const FEType fe_t,
unsigned int  s,
std::vector< unsigned int > &  di 
)
static

Fills the vector di with the local degree of freedom indices associated with side s of element elem Automatically decides which finite element class to use.

On a p-refined element, fe_t.order should be the base order of the element.

Referenced by libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::operator()().

◆ extra_hanging_dofs()

static bool libMesh::FEInterface::extra_hanging_dofs ( const FEType fe_t)
static
Returns
true if separate degrees of freedom must be allocated for vertex DoFs and edge/face DoFs at a hanging node.

◆ field_type() [1/2]

static FEFieldType libMesh::FEInterface::field_type ( const FEType fe_type)
static
Returns
The number of components of a vector-valued element. Scalar-valued elements return 1.

◆ field_type() [2/2]

static FEFieldType libMesh::FEInterface::field_type ( const FEFamily fe_family)
static
Returns
The number of components of a vector-valued element. Scalar-valued elements return 1.

◆ ifem_compute_data()

static void libMesh::FEInterface::ifem_compute_data ( const unsigned int  dim,
const FEType fe_t,
const Elem elem,
FEComputeData data 
)
staticprivate

◆ ifem_inverse_map() [1/2]

static Point libMesh::FEInterface::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 
)
staticprivate

◆ ifem_inverse_map() [2/2]

static void libMesh::FEInterface::ifem_inverse_map ( const unsigned int  dim,
const FEType fe_t,
const Elem elem,
const std::vector< Point > &  physical_points,
std::vector< Point > &  reference_points,
const Real  tolerance = TOLERANCE,
const bool  secure = true 
)
staticprivate

◆ ifem_map()

static Point libMesh::FEInterface::ifem_map ( const unsigned int  dim,
const FEType fe_t,
const Elem elem,
const Point p 
)
staticprivate

◆ ifem_n_dofs()

static unsigned int libMesh::FEInterface::ifem_n_dofs ( const unsigned int  dim,
const FEType fe_t,
const ElemType  t 
)
staticprivate

◆ ifem_n_dofs_at_node()

static unsigned int libMesh::FEInterface::ifem_n_dofs_at_node ( const unsigned int  dim,
const FEType fe_t,
const ElemType  t,
const unsigned int  n 
)
staticprivate

◆ ifem_n_dofs_per_elem()

static unsigned int libMesh::FEInterface::ifem_n_dofs_per_elem ( const unsigned int  dim,
const FEType fe_t,
const ElemType  t 
)
staticprivate

◆ ifem_n_shape_functions()

static unsigned int libMesh::FEInterface::ifem_n_shape_functions ( const unsigned int  dim,
const FEType fe_t,
const ElemType  t 
)
staticprivate

◆ ifem_nodal_soln()

static void libMesh::FEInterface::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 
)
staticprivate

◆ ifem_on_reference_element()

static bool libMesh::FEInterface::ifem_on_reference_element ( const Point p,
const ElemType  t,
const Real  eps 
)
staticprivate

◆ ifem_shape() [1/2]

static Real libMesh::FEInterface::ifem_shape ( const unsigned int  dim,
const FEType fe_t,
const ElemType  t,
const unsigned int  i,
const Point p 
)
staticprivate

◆ ifem_shape() [2/2]

static Real libMesh::FEInterface::ifem_shape ( const unsigned int  dim,
const FEType fe_t,
const Elem elem,
const unsigned int  i,
const Point p 
)
staticprivate

◆ inverse_map() [1/2]

static Point libMesh::FEInterface::inverse_map ( const unsigned int  dim,
const FEType fe_t,
const Elem elem,
const Point p,
const Real  tolerance = TOLERANCE,
const bool  secure = true 
)
static
Returns
The location (on the reference element) of the point p located in physical space. This function requires inverting the (probably nonlinear) transformation map, so it is not trivial. The optional parameter tolerance defines how close is "good enough." The map inversion iteration computes the sequence $ \{ p_n \} $, and the iteration is terminated when $ \|p - p_n\| < \mbox{\texttt{tolerance}} $

Referenced by libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::operator()(), FETest< order, family, elem_type >::testGradU(), FETest< order, family, elem_type >::testGradUComp(), and FETest< order, family, elem_type >::testU().

◆ inverse_map() [2/2]

static void libMesh::FEInterface::inverse_map ( const unsigned int  dim,
const FEType fe_t,
const Elem elem,
const std::vector< Point > &  physical_points,
std::vector< Point > &  reference_points,
const Real  tolerance = TOLERANCE,
const bool  secure = true 
)
static
Returns
The location (on the reference element) of the points physical_points located in physical space. This function requires inverting the (probably nonlinear) transformation map, so it is not trivial. The location of each point on the reference element is returned in the vector reference_points. The optional parameter tolerance defines how close is "good enough." The map inversion iteration computes the sequence $ \{ p_n \} $, and the iteration is terminated when $ \|p - p_n\| < \mbox{\texttt{tolerance}} $

◆ is_InfFE_elem()

static bool libMesh::FEInterface::is_InfFE_elem ( const ElemType  et)
staticprivate
Returns
true if et is an element to be processed by class InfFE, false otherwise or when !LIBMESH_ENABLE_INFINITE_ELEMENTS.

◆ map()

static Point libMesh::FEInterface::map ( unsigned int  dim,
const FEType fe_t,
const Elem elem,
const Point p 
)
static
Returns
The point in physical space corresponding to the reference point p which is passed in.

◆ max_order()

static unsigned int libMesh::FEInterface::max_order ( const FEType fe_t,
const ElemType el_t 
)
static
Returns
The maximum polynomial degree that the given finite element family can support on the given geometric element.

◆ n_dofs()

static unsigned int libMesh::FEInterface::n_dofs ( const unsigned int  dim,
const FEType fe_t,
const ElemType  t 
)
static
Returns
The number of shape functions associated with this finite element. Automatically decides which finite element class to use.

On a p-refined element, fe_t.order should be the total order of the element.

◆ n_dofs_at_node()

static unsigned int libMesh::FEInterface::n_dofs_at_node ( const unsigned int  dim,
const FEType fe_t,
const ElemType  t,
const unsigned int  n 
)
static
Returns
The number of dofs at node n for a finite element of type fe_t. Automatically decides which finite element class to use.

On a p-refined element, fe_t.order should be the total order of the element.

Referenced by libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::operator()().

◆ n_dofs_at_node_function()

static n_dofs_at_node_ptr libMesh::FEInterface::n_dofs_at_node_function ( const unsigned int  dim,
const FEType fe_t 
)
static
Returns
A function which evaluates n_dofs_at_node for the requested FE type and dimension.

◆ n_dofs_per_elem()

static unsigned int libMesh::FEInterface::n_dofs_per_elem ( const unsigned int  dim,
const FEType fe_t,
const ElemType  t 
)
static
Returns
The number of dofs interior to the element, not associated with any interior nodes. Automatically decides which finite element class to use.

On a p-refined element, fe_t.order should be the total order of the element.

◆ n_shape_functions()

static unsigned int libMesh::FEInterface::n_shape_functions ( const unsigned int  dim,
const FEType fe_t,
const ElemType  t 
)
static
Returns
The number of shape functions associated with this finite element of type fe_t. Automatically decides which finite element class to use.

On a p-refined element, fe_t.order should be the total order of the element.

◆ n_vec_dim()

static unsigned int libMesh::FEInterface::n_vec_dim ( const MeshBase mesh,
const FEType fe_type 
)
static
Returns
The number of components of a vector-valued element. Scalar-valued elements return 1.

◆ nodal_soln()

static void libMesh::FEInterface::nodal_soln ( const unsigned int  dim,
const FEType fe_t,
const Elem elem,
const std::vector< Number > &  elem_soln,
std::vector< Number > &  nodal_soln 
)
static

Build the nodal soln from the element soln.

This is the solution that will be plotted. Automatically passes the request to the appropriate finite element class member. To indicate that results from this specific implementation of nodal_soln should not be used, the vector nodal_soln is returned empty.

On a p-refined element, fe_t.order should be the base order of the element.

◆ on_reference_element()

static bool libMesh::FEInterface::on_reference_element ( const Point p,
const ElemType  t,
const Real  eps = TOLERANCE 
)
static
Returns
true if the point p is located on the reference element for element type t, false otherwise.

Since we are doing floating point comparisons, the parameter eps can be specified to indicate a tolerance. For example, $ \xi \le 1 $ becomes $ \xi \le 1 + \epsilon $.

◆ shape() [1/4]

static Real libMesh::FEInterface::shape ( const unsigned int  dim,
const FEType fe_t,
const ElemType  t,
const unsigned int  i,
const Point p 
)
static
Returns
The value of the $ i^{th} $ shape function at point p. This method allows you to specify the dimension, element type, and order directly. Automatically passes the request to the appropriate finite element class member.
Note
On a p-refined element, fe_t.order should be the total order of the element.

◆ shape() [2/4]

static Real libMesh::FEInterface::shape ( const unsigned int  dim,
const FEType fe_t,
const Elem elem,
const unsigned int  i,
const Point p 
)
static
Returns
The value of the $ i^{th} $ shape function at point p. This method allows you to specify the dimension, element type, and order directly. Automatically passes the request to the appropriate finite element class member.
Note
On a p-refined element, fe_t.order should be the base order of the element.

◆ shape() [3/4]

template<typename OutputType >
static void libMesh::FEInterface::shape ( const unsigned int  dim,
const FEType fe_t,
const ElemType  t,
const unsigned int  i,
const Point p,
OutputType &  phi 
)
static
Returns
The value of the $ i^{th} $ shape function at point p. This method allows you to specify the dimension, element type, and order directly. Automatically passes the request to the appropriate scalar finite element class member.
Note
On a p-refined element, fe_t.order should be the total order of the element.

◆ shape() [4/4]

template<typename OutputType >
static void libMesh::FEInterface::shape ( const unsigned int  dim,
const FEType fe_t,
const Elem elem,
const unsigned int  i,
const Point p,
OutputType &  phi 
)
static
Returns
The value of the $ i^{th} $ shape function at point p. This method allows you to specify the dimension, element type, and order directly. Automatically passes the request to the appropriate scalar finite element class member.
Note
On a p-refined element, fe_t.order should be the total order of the element.

The documentation for this class was generated from the following file: