Interface class which provides access to FE functions. 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 () |
template<> | |
void | shape (const unsigned int dim, const FEType &fe_t, const ElemType t, const unsigned int i, const Point &p, Real &phi) |
template<> | |
void | shape (const unsigned int dim, const FEType &fe_t, const Elem *elem, const unsigned int i, const Point &p, Real &phi) |
template<> | |
void | shape (const unsigned int dim, const FEType &fe_t, const ElemType t, const unsigned int i, const Point &p, RealGradient &phi) |
template<> | |
void | shape (const unsigned int dim, const FEType &fe_t, const Elem *elem, const unsigned int i, const Point &p, RealGradient &phi) |
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) |
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) |
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) |
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) |
static void | compute_constraints (DofConstraints &constraints, DofMap &dof_map, const unsigned int variable_number, const Elem *elem) |
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) |
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 () | |
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) |
Interface class which provides access to FE functions.
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.
Definition at line 64 of file fe_interface.h.
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.
|
private |
Empty constructor. Do not create an object of this type.
Definition at line 35 of file fe_interface.C.
|
inlinevirtual |
|
static |
Computes the constraint matrix contributions (for non-conforming adapted meshes) corresponding to variable number var_number
.
Definition at line 873 of file fe_interface.C.
References libMesh::BERNSTEIN, libMesh::CLOUGH, libMesh::FE< Dim, T >::compute_constraints(), libMesh::Elem::dim(), libMesh::FEType::family, libMesh::HERMITE, libMesh::HIERARCHIC, libMesh::L2_HIERARCHIC, libMesh::LAGRANGE, libMesh::LAGRANGE_VEC, libMesh::SZABAB, and libMesh::DofMap::variable_type().
|
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
.
fe_t.order
should be the base order of the element. Definition at line 837 of file fe_interface.C.
References data, ifem_compute_data(), is_InfFE_elem(), n_dofs(), libMesh::FEType::order, libMesh::Elem::p_level(), shape(), and libMesh::Elem::type().
Referenced by libMesh::MeshFunction::discontinuous_value(), libMesh::MeshFunction::operator()(), and libMesh::System::point_value().
|
static |
Computes the constraint matrix contributions (for periodic boundary conditions) corresponding to variable number var_number
.
Definition at line 1017 of file fe_interface.C.
References libMesh::FEGenericBase< OutputType >::compute_periodic_constraints(), and mesh.
|
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.
Definition at line 536 of file fe_interface.C.
References libMesh::FEType::order.
Referenced by libMesh::FEGenericBase< FEOutputType< T >::type >::coarsened_dof_values(), and libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::operator()().
|
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.
Definition at line 523 of file fe_interface.C.
References libMesh::FEType::order.
Referenced by libMesh::FEGenericBase< FEOutputType< T >::type >::coarsened_dof_values(), libMesh::FEGenericBase< FEOutputType< T >::type >::compute_periodic_constraints(), libMesh::FEGenericBase< FEOutputType< T >::type >::compute_proj_constraints(), and libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::operator()().
|
static |
true
if separate degrees of freedom must be allocated for vertex DoFs and edge/face DoFs at a hanging node. Definition at line 1412 of file fe_interface.C.
References libMesh::BERNSTEIN, libMesh::CLOUGH, libMesh::FEType::family, libMesh::HERMITE, libMesh::HIERARCHIC, libMesh::L2_HIERARCHIC, libMesh::L2_LAGRANGE, libMesh::LAGRANGE, libMesh::LAGRANGE_VEC, libMesh::MONOMIAL, libMesh::NEDELEC_ONE, libMesh::SUBDIVISION, libMesh::SZABAB, and libMesh::XYZ.
Referenced by libMesh::DofMap::_dof_indices(), libMesh::DofMap::old_dof_indices(), and libMesh::DofMap::reinit().
|
static |
Definition at line 1437 of file fe_interface.C.
References libMesh::FEType::family.
Referenced by libMesh::ExactSolution::_compute_error(), libMesh::EquationSystems::build_parallel_solution_vector(), libMesh::EquationSystems::build_variable_names(), libMesh::ExactSolution::compute_error(), libMesh::FEGenericBase< FEOutputType< T >::type >::determine_calculations(), libMesh::ExactSolution::error_norm(), and libMesh::FEMSystem::init_context().
|
static |
Definition at line 1442 of file fe_interface.C.
References libMesh::LAGRANGE_VEC, libMesh::NEDELEC_ONE, libMesh::TYPE_SCALAR, and libMesh::TYPE_VECTOR.
|
staticprivate |
Definition at line 808 of file fe_interface_inf_fe.C.
References libMesh::InfFE< Dim, T_radial, T_map >::compute_data(), data, libMesh::INFINITE_MAP, libMesh::JACOBI_20_00, libMesh::JACOBI_30_00, libMesh::LAGRANGE, libMesh::LEGENDRE, and libMesh::FEType::radial_family.
Referenced by compute_data().
|
staticprivate |
Definition at line 447 of file fe_interface_inf_fe.C.
References libMesh::CARTESIAN, libMesh::ELLIPSOIDAL, libMesh::FEType::inf_map, libMesh::InfFE< Dim, T_radial, T_map >::inverse_map(), and libMesh::SPHERICAL.
Referenced by inverse_map().
|
staticprivate |
Definition at line 540 of file fe_interface_inf_fe.C.
References libMesh::CARTESIAN, libMesh::FEType::inf_map, and libMesh::InfFE< Dim, T_radial, T_map >::inverse_map().
|
staticprivate |
Definition at line 416 of file fe_interface_inf_fe.C.
References libMesh::CARTESIAN, libMesh::ELLIPSOIDAL, libMesh::FEType::inf_map, libMesh::InfFE< Dim, T_radial, T_map >::map(), and libMesh::SPHERICAL.
Referenced by map().
|
staticprivate |
Definition at line 71 of file fe_interface_inf_fe.C.
References libMesh::InfFE< Dim, T_radial, T_map >::n_dofs().
Referenced by n_dofs().
|
staticprivate |
Definition at line 102 of file fe_interface_inf_fe.C.
References libMesh::InfFE< Dim, T_radial, T_map >::n_dofs_at_node().
Referenced by n_dofs_at_node().
|
staticprivate |
Definition at line 135 of file fe_interface_inf_fe.C.
References libMesh::InfFE< Dim, T_radial, T_map >::n_dofs_per_elem().
Referenced by n_dofs_per_elem().
|
staticprivate |
Definition at line 39 of file fe_interface_inf_fe.C.
References libMesh::InfFE< Dim, T_radial, T_map >::n_shape_functions().
Referenced by n_shape_functions().
|
staticprivate |
Definition at line 166 of file fe_interface_inf_fe.C.
References libMesh::CARTESIAN, libMesh::FEType::inf_map, libMesh::INFINITE_MAP, libMesh::JACOBI_20_00, libMesh::JACOBI_30_00, libMesh::LAGRANGE, libMesh::LEGENDRE, nodal_soln(), libMesh::InfFE< Dim, T_radial, T_map >::nodal_soln(), and libMesh::FEType::radial_family.
Referenced by nodal_soln().
|
staticprivate |
Definition at line 602 of file fe_interface_inf_fe.C.
References libMesh::FEAbstract::on_reference_element().
|
staticprivate |
Definition at line 612 of file fe_interface_inf_fe.C.
References libMesh::INFINITE_MAP, libMesh::JACOBI_20_00, libMesh::JACOBI_30_00, libMesh::LAGRANGE, libMesh::LEGENDRE, libMesh::FEType::radial_family, and libMesh::InfFE< Dim, T_radial, T_map >::shape().
Referenced by shape().
|
staticprivate |
Definition at line 710 of file fe_interface_inf_fe.C.
References libMesh::INFINITE_MAP, libMesh::JACOBI_20_00, libMesh::JACOBI_30_00, libMesh::LAGRANGE, libMesh::LEGENDRE, libMesh::FEType::radial_family, and libMesh::InfFE< Dim, T_radial, T_map >::shape().
|
static |
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 , and the iteration is terminated when Definition at line 590 of file fe_interface.C.
References ifem_inverse_map(), is_InfFE_elem(), and libMesh::Elem::type().
Referenced by libMesh::HPCoarsenTest::add_projection(), libMesh::FEMContext::build_new_fe(), libMesh::FEGenericBase< FEOutputType< T >::type >::coarsened_dof_values(), libMesh::FEAbstract::compute_node_constraints(), libMesh::FEGenericBase< FEOutputType< T >::type >::compute_periodic_constraints(), libMesh::FEAbstract::compute_periodic_node_constraints(), libMesh::FEGenericBase< FEOutputType< T >::type >::compute_proj_constraints(), libMesh::InfQuad4::contains_point(), libMesh::InfPrism::contains_point(), libMesh::InfHex::contains_point(), libMesh::MeshFunction::discontinuous_gradient(), libMesh::MeshFunction::discontinuous_value(), libMesh::MeshFunction::gradient(), libMesh::MeshFunction::hessian(), inverse_map(), libMesh::DGFEMContext::neighbor_side_fe_reinit(), libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::operator()(), libMesh::MeshFunction::operator()(), libMesh::System::point_gradient(), libMesh::System::point_hessian(), libMesh::Elem::point_test(), libMesh::System::point_value(), libMesh::JumpErrorEstimator::reinit_sides(), and libMesh::HPCoarsenTest::select_refinement().
|
static |
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 , and the iteration is terminated when Definition at line 610 of file fe_interface.C.
References libMesh::err, ifem_inverse_map(), inverse_map(), is_InfFE_elem(), libMesh::TypeVector< T >::size(), and libMesh::Elem::type().
|
staticprivate |
true
if et
is an element to be processed by class InfFE
, false otherwise or when !LIBMESH_ENABLE_INFINITE_ELEMENTS. Definition at line 45 of file fe_interface.C.
Referenced by compute_data(), inverse_map(), map(), n_dofs(), n_dofs_at_node(), n_dofs_per_elem(), n_shape_functions(), nodal_soln(), and shape().
|
static |
p
which is passed in. Definition at line 574 of file fe_interface.C.
References ifem_map(), is_InfFE_elem(), and libMesh::Elem::type().
Referenced by libMesh::Elem::point_test().
Definition at line 1039 of file fe_interface.C.
References libMesh::BERNSTEIN, libMesh::CLOUGH, libMesh::EDGE2, libMesh::EDGE3, libMesh::EDGE4, libMesh::FEType::family, libMesh::HERMITE, libMesh::HEX20, libMesh::HEX27, libMesh::HEX8, libMesh::HIERARCHIC, libMesh::L2_HIERARCHIC, libMesh::L2_LAGRANGE, libMesh::LAGRANGE, libMesh::LAGRANGE_VEC, libMesh::MONOMIAL, libMesh::NEDELEC_ONE, libMesh::PRISM15, libMesh::PRISM18, libMesh::PRISM6, libMesh::PYRAMID13, libMesh::PYRAMID14, libMesh::PYRAMID5, libMesh::QUAD4, libMesh::QUAD8, libMesh::QUAD9, libMesh::QUADSHELL4, libMesh::QUADSHELL8, libMesh::SUBDIVISION, libMesh::SZABAB, libMesh::TET10, libMesh::TET4, libMesh::TRI3, libMesh::TRI3SUBDIVISION, libMesh::TRI6, libMesh::TRISHELL3, and libMesh::XYZ.
Referenced by libMesh::DofMap::reinit().
|
static |
On a p-refined element, fe_t.order
should be the total order of the element.
Definition at line 454 of file fe_interface.C.
References ifem_n_dofs(), is_InfFE_elem(), and libMesh::FEType::order.
Referenced by libMesh::DofMap::_dof_indices(), libMesh::FEGenericBase< FEOutputType< T >::type >::coarsened_dof_values(), compute_data(), libMesh::FEAbstract::compute_node_constraints(), libMesh::FEAbstract::compute_periodic_node_constraints(), libMesh::InfFE< Dim, T_radial, T_map >::n_dofs(), and libMesh::HPCoarsenTest::select_refinement().
|
static |
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.
Definition at line 473 of file fe_interface.C.
References ifem_n_dofs_at_node(), is_InfFE_elem(), and libMesh::FEType::order.
Referenced by libMesh::DofMap::_dof_indices(), libMesh::FEGenericBase< FEOutputType< T >::type >::coarsened_dof_values(), libMesh::InfFE< Dim, T_radial, T_map >::compute_shape_indices(), libMesh::InfFE< Dim, T_radial, T_map >::n_dofs_at_node(), n_dofs_at_node_function(), libMesh::DofMap::old_dof_indices(), libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::operator()(), and libMesh::DofMap::reinit().
|
static |
Definition at line 493 of file fe_interface.C.
References n_dofs_at_node().
Referenced by libMesh::DofMap::_dof_indices(), and libMesh::DofMap::old_dof_indices().
|
static |
On a p-refined element, fe_t.order
should be the total order of the element.
Definition at line 504 of file fe_interface.C.
References ifem_n_dofs_per_elem(), is_InfFE_elem(), and libMesh::FEType::order.
Referenced by libMesh::DofMap::_dof_indices(), libMesh::InfFE< Dim, T_radial, T_map >::compute_shape_indices(), libMesh::InfFE< Dim, T_radial, T_map >::n_dofs_per_elem(), libMesh::DofMap::old_dof_indices(), and libMesh::DofMap::reinit().
|
static |
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.
Definition at line 428 of file fe_interface.C.
References ifem_n_shape_functions(), is_InfFE_elem(), and libMesh::FEType::order.
|
static |
Definition at line 1454 of file fe_interface.C.
References libMesh::FEType::family, libMesh::LAGRANGE_VEC, mesh, and libMesh::NEDELEC_ONE.
Referenced by libMesh::ExactSolution::_compute_error(), libMesh::EquationSystems::build_parallel_solution_vector(), and libMesh::EquationSystems::build_variable_names().
|
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.
Definition at line 550 of file fe_interface.C.
References ifem_nodal_soln(), is_InfFE_elem(), libMesh::FEType::order, and libMesh::Elem::type().
Referenced by libMesh::EquationSystems::build_discontinuous_solution_vector(), libMesh::EquationSystems::build_parallel_solution_vector(), ifem_nodal_soln(), libMesh::EnsightIO::write_scalar_ascii(), and libMesh::EnsightIO::write_vector_ascii().
|
static |
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, becomes .
Definition at line 647 of file fe_interface.C.
References libMesh::FEAbstract::on_reference_element().
Referenced by libMesh::InfQuad4::contains_point(), libMesh::InfPrism::contains_point(), libMesh::InfHex::contains_point(), and libMesh::Elem::point_test().
|
static |
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.fe_t.order
should be the total order of the element. Definition at line 657 of file fe_interface.C.
References ifem_shape(), is_InfFE_elem(), and libMesh::FEType::order.
Referenced by libMesh::InfFE< Dim, T_radial, T_map >::compute_data(), compute_data(), libMesh::FEAbstract::compute_node_constraints(), libMesh::FEAbstract::compute_periodic_node_constraints(), shape(), and libMesh::InfFE< Dim, T_radial, T_map >::shape().
|
static |
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.fe_t.order
should be the base order of the element. Definition at line 675 of file fe_interface.C.
References ifem_shape(), is_InfFE_elem(), libMesh::FEType::order, shape(), and libMesh::Elem::type().
|
static |
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.fe_t.order
should be the total order of the element.
|
static |
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.fe_t.order
should be the total order of the element. void libMesh::FEInterface::shape | ( | const unsigned int | dim, |
const FEType & | fe_t, | ||
const ElemType | t, | ||
const unsigned int | i, | ||
const Point & | p, | ||
Real & | phi | ||
) |
Definition at line 694 of file fe_interface.C.
void libMesh::FEInterface::shape | ( | const unsigned int | dim, |
const FEType & | fe_t, | ||
const Elem * | elem, | ||
const unsigned int | i, | ||
const Point & | p, | ||
Real & | phi | ||
) |
Definition at line 735 of file fe_interface.C.
void libMesh::FEInterface::shape | ( | const unsigned int | dim, |
const FEType & | fe_t, | ||
const ElemType | t, | ||
const unsigned int | i, | ||
const Point & | p, | ||
RealGradient & | phi | ||
) |
Definition at line 775 of file fe_interface.C.
void libMesh::FEInterface::shape | ( | const unsigned int | dim, |
const FEType & | fe_t, | ||
const Elem * | elem, | ||
const unsigned int | i, | ||
const Point & | p, | ||
RealGradient & | phi | ||
) |
Definition at line 807 of file fe_interface.C.