libMesh
|
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) |
typedef Real(* | shape_ptr) (const FEType fe_t, const Elem *elem, const unsigned int i, const Point &p, const bool add_p_level) |
Typedef for pointer to a function that returns FE shape function values. More... | |
typedef Real(* | shape_deriv_ptr) (const FEType fet, const Elem *elem, const unsigned int i, const unsigned int j, const Point &p, const bool add_p_level) |
Typedef for pointer to a function that returns FE shape function derivative values. More... | |
typedef Real(* | shape_second_deriv_ptr) (const FEType fet, const Elem *elem, const unsigned int i, const unsigned int j, const Point &p, const bool add_p_level) |
Typedef for pointer to a function that returns FE shape function second derivative values. More... | |
Public Member Functions | |
virtual | ~FEInterface ()=default |
Destructor. More... | |
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 FEType &fe_t, const Elem *elem, const unsigned int i, const Point &p, Real &phi) |
template<> | |
void | shape (const FEType &fe_t, int extra_order, const Elem *elem, const unsigned int i, const Point &p, Real &phi) |
template<> | |
void | shapes (const unsigned int dim, const FEType &fe_t, const Elem *elem, const unsigned int i, const std::vector< Point > &p, std::vector< Real > &phi, const bool add_p_level) |
template<> | |
void | all_shapes (const unsigned int dim, const FEType &fe_t, const Elem *elem, const std::vector< Point > &p, std::vector< std::vector< Real >> &phi, const bool add_p_level) |
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 FEType &fe_t, const Elem *elem, const unsigned int i, const Point &p, RealGradient &phi) |
template<> | |
void | shape (const FEType &fe_t, int extra_order, const Elem *elem, const unsigned int i, const Point &p, RealGradient &phi) |
template<> | |
void | shapes (const unsigned int dim, const FEType &fe_t, const Elem *elem, const unsigned int i, const std::vector< Point > &p, std::vector< RealGradient > &phi, const bool add_p_level) |
template<> | |
void | all_shapes (const unsigned int dim, const FEType &fe_t, const Elem *elem, const std::vector< Point > &p, std::vector< std::vector< RealGradient >> &phi, const bool add_p_level) |
template<> | |
void | shape_derivs (const FEType &fe_t, const Elem *elem, const unsigned int i, const unsigned int j, const std::vector< Point > &p, std::vector< Real > &dphi, const bool add_p_level) |
template<> | |
void | all_shape_derivs (const unsigned int dim, const FEType &fe_t, const Elem *elem, const std::vector< Point > &p, std::vector< std::vector< Real >> *comps[3], const bool add_p_level) |
template<> | |
void | shape_derivs (const FEType &fe_t, const Elem *elem, const unsigned int i, const unsigned int j, const std::vector< Point > &p, std::vector< RealGradient > &dphi, const bool add_p_level) |
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_shape_functions (const FEType &fe_t, const Elem *elem, const bool add_p_level=true) |
static unsigned int | n_shape_functions (const FEType &fe_t, const int extra_order, const Elem *elem) |
Same as above, but ignores the elem->p_level() and uses the specified extra_order instead. More... | |
static unsigned int | n_dofs (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 Elem *elem) |
Similar to the function above but takes an Elem * and accounts for p-refinement internally, if any. More... | |
static unsigned int | n_dofs (const FEType &fe_t, const Elem *elem, const bool add_p_level=true) |
static unsigned int | n_dofs (const FEType &fe_t, int extra_order, const Elem *elem) |
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 n_dofs_at_node_ptr | n_dofs_at_node_function (const FEType &fe_t, const Elem *elem) |
static unsigned int | n_dofs_at_node (const FEType &fe_t, const Elem *elem, const unsigned int n, const bool add_p_level=true) |
static unsigned int | n_dofs_at_node (const FEType &fe_t, const int extra_order, const Elem *elem, const unsigned int n) |
static unsigned int | n_dofs_per_elem (const unsigned int dim, const FEType &fe_t, const ElemType t) |
static unsigned int | n_dofs_per_elem (const FEType &fe_t, const Elem *elem, const bool add_p_level=true) |
static unsigned int | n_dofs_per_elem (const FEType &fe_t, const int extra_order, const Elem *elem) |
Same thing but internally elem->p_level() is ignored and extra_order is used instead. More... | |
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, const bool add_p_level=true) |
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, const bool add_p_level=true) |
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, const bool add_p_level=true, const unsigned int vdim=1) |
Build the nodal soln from the element soln. More... | |
static void | side_nodal_soln (const FEType &fe_t, const Elem *elem, const unsigned int side, const std::vector< Number > &elem_soln, std::vector< Number > &nodal_soln, const bool add_p_level=true, const unsigned int vdim=1) |
Build the nodal soln on one side from the (full) element soln. More... | |
static Point | map (unsigned int dim, const FEType &fe_t, const Elem *elem, const Point &p) |
This is now deprecated; use FEMap::map instead. More... | |
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. More... | |
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) |
This is now deprecated; use FEMap::inverse_map instead. More... | |
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) |
static Real | shape (const FEType &fe_t, const Elem *elem, const unsigned int i, const Point &p, const bool add_p_level=true) |
static Real | shape (const FEType &fe_t, int extra_order, 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) |
template<typename OutputType > | |
static void | shape (const FEType &fe_t, const Elem *elem, const unsigned int i, const Point &p, OutputType &phi) |
template<typename OutputType > | |
static void | shape (const FEType &fe_t, int extra_order, const Elem *elem, const unsigned int i, const Point &p, OutputType &phi) |
template<typename OutputType > | |
static void | shapes (const unsigned int dim, const FEType &fe_t, const Elem *elem, const unsigned int i, const std::vector< Point > &p, std::vector< OutputType > &phi, const bool add_p_level=true) |
Fills phi with the values of the \( i^{th} \) shape function at point p . More... | |
template<typename OutputType > | |
static void | all_shapes (const unsigned int dim, const FEType &fe_t, const Elem *elem, const std::vector< Point > &p, std::vector< std::vector< OutputType >> &phi, const bool add_p_level=true) |
static shape_ptr | shape_function (const unsigned int dim, const FEType &fe_t, const ElemType t) |
static shape_ptr | shape_function (const FEType &fe_t, const Elem *elem) |
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) |
static Real | shape_deriv (const unsigned int dim, const FEType &fe_t, const Elem *elem, const unsigned int i, const unsigned int j, const Point &p) |
static Real | shape_deriv (const FEType &fe_t, const Elem *elem, const unsigned int i, const unsigned int j, const Point &p) |
static Real | shape_deriv (const FEType &fe_t, int extra_order, const Elem *elem, const unsigned int i, const unsigned int j, const Point &p) |
Non-deprecated version of function above. More... | |
template<typename OutputType > | |
static void | shape_derivs (const FEType &fe_t, const Elem *elem, const unsigned int i, const unsigned int j, const std::vector< Point > &p, std::vector< OutputType > &dphi, const bool add_p_level=true) |
Fills dphi with the derivatives of the \( i^{th} \) shape function at point p in direction j . More... | |
template<typename OutputType > | |
static void | all_shape_derivs (const unsigned int dim, const FEType &fe_t, const Elem *elem, const std::vector< Point > &p, std::vector< std::vector< OutputType >> *comps[3], const bool add_p_level=true) |
static shape_deriv_ptr | shape_deriv_function (const unsigned int dim, const FEType &fe_t, const ElemType t) |
static shape_deriv_ptr | shape_deriv_function (const FEType &fe_t, const Elem *elem) |
Non-deprecated version of the function above. More... | |
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 Real | shape_second_deriv (const unsigned int dim, const FEType &fe_t, const Elem *elem, const unsigned int i, const unsigned int j, const Point &p) |
static Real | shape_second_deriv (const FEType &fe_t, const Elem *elem, const unsigned int i, const unsigned int j, const Point &p) |
static Real | shape_second_deriv (const FEType &fe_t, int extra_order, const Elem *elem, const unsigned int i, const unsigned int j, const Point &p) |
Non-deprecated version of function above taking an extra_order parameter. More... | |
static shape_second_deriv_ptr | shape_second_deriv_function (const unsigned int dim, const FEType &fe_t, const ElemType t) |
static shape_second_deriv_ptr | shape_second_deriv_function (const FEType &fe_t, const Elem *elem) |
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 bool | orientation_dependent (const FEFamily &fe_family) |
static unsigned int | n_vec_dim (const MeshBase &mesh, const FEType &fe_type) |
static FEContinuity | get_continuity (const FEType &fe_type) |
Returns the input FEType's FEContinuity based on the underlying FEFamily and potentially the Order, although we do not currently support FEs with order-dependent continuity. More... | |
static bool | is_hierarchic (const FEType &fe_type) |
Returns whether or not the input FEType's higher-order shape functions are always hierarchic. More... | |
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_shape_functions (const FEType &fe_t, const Elem *elem) |
static unsigned int | ifem_n_dofs (const unsigned int dim, const FEType &fe_t, const ElemType t) |
static unsigned int | ifem_n_dofs (const FEType &fe_t, const Elem *elem) |
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_at_node (const FEType &fe_t, const Elem *elem, const unsigned int n) |
static unsigned int | ifem_n_dofs_per_elem (const unsigned int dim, const FEType &fe_t, const ElemType t) |
static unsigned int | ifem_n_dofs_per_elem (const FEType &fe_t, const Elem *elem) |
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 Real | ifem_shape (const FEType &fe_t, const Elem *t, const unsigned int i, const Point &p) |
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) |
static Real | ifem_shape_deriv (const unsigned int dim, const FEType &fe_t, const Elem *elem, const unsigned int i, const unsigned int j, const Point &p) |
static Real | ifem_shape_deriv (const FEType &fe_t, const Elem *elem, const unsigned int i, const unsigned int j, const Point &p) |
static void | ifem_compute_data (const unsigned int dim, const FEType &fe_t, const Elem *elem, FEComputeData &data) |
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 65 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 151 of file fe_interface.h.
typedef Real(* libMesh::FEInterface::shape_deriv_ptr) (const FEType fet, const Elem *elem, const unsigned int i, const unsigned int j, const Point &p, const bool add_p_level) |
Typedef for pointer to a function that returns FE shape function derivative values.
The p_level()
of the passed-in elem
is accounted for internally when the add_p_level
flag is set to true. For more information, see fe.h.
Definition at line 648 of file fe_interface.h.
typedef Real(* libMesh::FEInterface::shape_ptr) (const FEType fe_t, const Elem *elem, const unsigned int i, const Point &p, const bool add_p_level) |
Typedef for pointer to a function that returns FE shape function values.
The p_level()
of the passed-in elem
is accounted for internally when the add_p_level
flag is set to true. For more information, see fe.h.
Definition at line 541 of file fe_interface.h.
typedef Real(* libMesh::FEInterface::shape_second_deriv_ptr) (const FEType fet, const Elem *elem, const unsigned int i, const unsigned int j, const Point &p, const bool add_p_level) |
Typedef for pointer to a function that returns FE shape function second derivative values.
The p_level()
of the passed-in elem
is accounted for internally when the add_p_level
flag is set to true. For more information, see fe.h.
Definition at line 753 of file fe_interface.h.
|
private |
Empty constructor.
Do not create an object of this type.
Definition at line 42 of file fe_interface.C.
|
virtualdefault |
Destructor.
|
static |
Referenced by libMesh::FEMap::init_reference_to_physical_map().
void libMesh::FEInterface::all_shape_derivs | ( | const unsigned int | dim, |
const FEType & | fe_t, | ||
const Elem * | elem, | ||
const std::vector< Point > & | p, | ||
std::vector< std::vector< Real >> * | comps[3], | ||
const bool | add_p_level | ||
) |
Definition at line 1552 of file fe_interface.C.
References dim, libMesh::index_range(), and libMesh::make_range().
|
static |
Referenced by libMesh::FEMap::init_reference_to_physical_map().
void libMesh::FEInterface::all_shapes | ( | const unsigned int | dim, |
const FEType & | fe_t, | ||
const Elem * | elem, | ||
const std::vector< Point > & | p, | ||
std::vector< std::vector< Real >> & | phi, | ||
const bool | add_p_level | ||
) |
Definition at line 1092 of file fe_interface.C.
References dim, and libMesh::index_range().
void libMesh::FEInterface::all_shapes | ( | const unsigned int | dim, |
const FEType & | fe_t, | ||
const Elem * | elem, | ||
const std::vector< Point > & | p, | ||
std::vector< std::vector< RealGradient >> & | phi, | ||
const bool | add_p_level | ||
) |
Definition at line 1308 of file fe_interface.C.
References dim.
|
static |
Computes the constraint matrix contributions (for non-conforming adapted meshes) corresponding to variable number var_number
.
Definition at line 1963 of file fe_interface.C.
References libMesh::CLOUGH, libMesh::FE< Dim, T >::compute_constraints(), libMesh::Elem::dim(), libMesh::FEType::family, libMesh::HERMITE, libMesh::HIERARCHIC, libMesh::HIERARCHIC_VEC, libMesh::LAGRANGE, libMesh::LAGRANGE_VEC, libMesh::libmesh_assert(), libMesh::SIDE_HIERARCHIC, 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.dim
argument. Definition at line 1889 of file fe_interface.C.
References libMesh::FEGenericBase< OutputType >::build(), dim, libMesh::FEComputeData::dshape, ifem_compute_data(), libMesh::FEComputeData::init(), is_InfFE_elem(), libMesh::FEComputeData::local_transform, n_dofs(), libMesh::FEComputeData::need_derivative(), libMesh::FEComputeData::p, libMesh::FEComputeData::shape, shape(), shape_deriv(), and libMesh::Elem::type().
Referenced by libMesh::MeshFunction::_gradient_on_elem(), libMesh::MeshFunction::discontinuous_value(), libMesh::DTKEvaluator::evaluate(), ifem_compute_data(), libMesh::MeshFunction::operator()(), libMesh::System::point_gradient(), and libMesh::System::point_value().
|
static |
Computes the constraint matrix contributions (for periodic boundary conditions) corresponding to variable number var_number
.
Definition at line 2099 of file fe_interface.C.
References libMesh::FEGenericBase< OutputType >::compute_periodic_constraints(), field_type(), mesh, libMesh::TYPE_SCALAR, libMesh::TYPE_VECTOR, and libMesh::DofMap::variable_type().
|
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.
dim
argument. Definition at line 611 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 >::ProjectEdges::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.
dim
argument. Definition at line 597 of file fe_interface.C.
References libMesh::FEType::order.
Referenced by alternative_fe_assembly(), 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 >::ProjectSides::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 2613 of file fe_interface.C.
References libMesh::CLOUGH, libMesh::FEType::family, libMesh::HERMITE, libMesh::HIERARCHIC, libMesh::HIERARCHIC_VEC, libMesh::L2_HIERARCHIC, libMesh::L2_HIERARCHIC_VEC, libMesh::L2_LAGRANGE, libMesh::L2_LAGRANGE_VEC, libMesh::L2_RAVIART_THOMAS, libMesh::LAGRANGE, libMesh::LAGRANGE_VEC, libMesh::MONOMIAL, libMesh::MONOMIAL_VEC, libMesh::NEDELEC_ONE, libMesh::RAVIART_THOMAS, libMesh::SIDE_HIERARCHIC, libMesh::SUBDIVISION, and libMesh::XYZ.
Referenced by libMesh::DofMap::_dof_indices(), libMesh::DofMap::_node_dof_indices(), libMesh::DofMap::old_dof_indices(), libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::SortAndCopy::operator()(), libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::ProjectVertices::operator()(), and libMesh::DofMap::reinit().
|
static |
Definition at line 2642 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::System::calculate_norm(), libMesh::ExactSolution::compute_error(), compute_periodic_constraints(), libMesh::FEGenericBase< FEOutputType< T >::type >::determine_calculations(), libMesh::InfFE< Dim, T_radial, T_map >::determine_calculations(), libMesh::ExactSolution::error_norm(), libMesh::EquationSystems::find_variable_numbers(), libMesh::FEMSystem::init_context(), main(), libMesh::Variable::n_components(), n_vec_dim(), and libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::SubFunctor::SubFunctor().
|
static |
Definition at line 2647 of file fe_interface.C.
References libMesh::TYPE_SCALAR, and libMesh::TYPE_VECTOR.
|
static |
Returns the input FEType's FEContinuity based on the underlying FEFamily and potentially the Order, although we do not currently support FEs with order-dependent continuity.
These should exactly match the FEBase::get_continuity() specializations/overrides for the different FE types.
Definition at line 2725 of file fe_interface.C.
References libMesh::BERNSTEIN, libMesh::C_ONE, libMesh::C_ZERO, libMesh::CLOUGH, libMesh::DISCONTINUOUS, libMesh::Utility::enum_to_string(), libMesh::FEType::family, libMesh::H_CURL, libMesh::H_DIV, libMesh::HERMITE, libMesh::HIERARCHIC, libMesh::HIERARCHIC_VEC, libMesh::INFINITE_MAP, libMesh::JACOBI_20_00, libMesh::JACOBI_30_00, libMesh::L2_HIERARCHIC, libMesh::L2_HIERARCHIC_VEC, libMesh::L2_LAGRANGE, libMesh::L2_LAGRANGE_VEC, libMesh::L2_RAVIART_THOMAS, libMesh::LAGRANGE, libMesh::LAGRANGE_VEC, libMesh::LEGENDRE, libMesh::MONOMIAL, libMesh::MONOMIAL_VEC, libMesh::NEDELEC_ONE, libMesh::RATIONAL_BERNSTEIN, libMesh::RAVIART_THOMAS, libMesh::SCALAR, libMesh::SIDE_DISCONTINUOUS, libMesh::SIDE_HIERARCHIC, libMesh::SUBDIVISION, libMesh::SZABAB, and libMesh::XYZ.
Referenced by libMesh::EquationSystems::build_discontinuous_solution_vector(), FETest< order, family, elem_type >::testFEInterface(), and libMesh::DofMap::use_coupled_neighbor_dofs().
|
staticprivate |
Definition at line 825 of file fe_interface_inf_fe.C.
References compute_data(), and dim.
Referenced by compute_data().
|
staticprivate |
Definition at line 578 of file fe_interface_inf_fe.C.
References libMesh::CARTESIAN, dim, libMesh::ELLIPSOIDAL, libMesh::Utility::enum_to_string(), libMesh::FEType::inf_map, libMesh::InfFE< Dim, T_radial, T_map >::inverse_map(), and libMesh::SPHERICAL.
Referenced by inverse_map().
|
staticprivate |
Definition at line 671 of file fe_interface_inf_fe.C.
References libMesh::CARTESIAN, dim, libMesh::Utility::enum_to_string(), libMesh::FEType::inf_map, and libMesh::InfFE< Dim, T_radial, T_map >::inverse_map().
|
staticprivate |
Definition at line 547 of file fe_interface_inf_fe.C.
References libMesh::CARTESIAN, dim, libMesh::ELLIPSOIDAL, libMesh::Utility::enum_to_string(), libMesh::FEType::inf_map, libMesh::InfFE< Dim, T_radial, T_map >::map(), and libMesh::SPHERICAL.
Referenced by map().
|
staticprivate |
Definition at line 103 of file fe_interface_inf_fe.C.
References dim, and libMesh::InfFE< Dim, T_radial, T_map >::n_dofs().
Referenced by n_dofs().
|
staticprivate |
Definition at line 137 of file fe_interface_inf_fe.C.
References libMesh::Elem::dim(), and libMesh::InfFE< Dim, T_radial, T_map >::n_dofs().
|
staticprivate |
Definition at line 167 of file fe_interface_inf_fe.C.
References dim, and libMesh::InfFE< Dim, T_radial, T_map >::n_dofs_at_node().
Referenced by n_dofs_at_node().
|
staticprivate |
Definition at line 201 of file fe_interface_inf_fe.C.
References libMesh::Elem::dim(), and libMesh::InfFE< Dim, T_radial, T_map >::n_dofs_at_node().
|
staticprivate |
Definition at line 234 of file fe_interface_inf_fe.C.
References dim, and libMesh::InfFE< Dim, T_radial, T_map >::n_dofs_per_elem().
Referenced by n_dofs_per_elem().
|
staticprivate |
Definition at line 267 of file fe_interface_inf_fe.C.
References libMesh::Elem::dim(), and libMesh::InfFE< Dim, T_radial, T_map >::n_dofs_per_elem().
|
staticprivate |
Definition at line 40 of file fe_interface_inf_fe.C.
References dim, and libMesh::InfFE< Dim, T_radial, T_map >::n_shape_functions().
Referenced by n_shape_functions().
|
staticprivate |
Definition at line 73 of file fe_interface_inf_fe.C.
References libMesh::Elem::dim(), and libMesh::InfFE< Dim, T_radial, T_map >::n_shape_functions().
|
staticprivate |
Definition at line 297 of file fe_interface_inf_fe.C.
References libMesh::CARTESIAN, dim, libMesh::Utility::enum_to_string(), 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 733 of file fe_interface_inf_fe.C.
References libMesh::FEAbstract::on_reference_element().
|
staticprivate |
Definition at line 742 of file fe_interface_inf_fe.C.
References shape().
Referenced by shape().
|
staticprivate |
Definition at line 755 of file fe_interface_inf_fe.C.
References shape().
|
staticprivate |
Definition at line 769 of file fe_interface_inf_fe.C.
References dim, libMesh::Elem::dim(), and shape().
|
staticprivate |
Definition at line 797 of file fe_interface_inf_fe.C.
References shape_deriv().
Referenced by shape_deriv().
|
staticprivate |
Definition at line 783 of file fe_interface_inf_fe.C.
References shape_deriv().
|
staticprivate |
Definition at line 812 of file fe_interface_inf_fe.C.
References dim, libMesh::Elem::dim(), and shape_deriv().
|
static |
This is now deprecated; use FEMap::inverse_map instead.
Definition at line 692 of file fe_interface.C.
References dim, ifem_inverse_map(), is_InfFE_elem(), and libMesh::Elem::type().
Referenced by inverse_map().
|
static |
This is now deprecated; use FEMap::inverse_map instead.
Definition at line 712 of file fe_interface.C.
References dim, libMesh::err, ifem_inverse_map(), inverse_map(), is_InfFE_elem(), and libMesh::Elem::type().
|
static |
Returns whether or not the input FEType's higher-order shape functions are always hierarchic.
Definition at line 2688 of file fe_interface.C.
References libMesh::BERNSTEIN, libMesh::CLOUGH, libMesh::Utility::enum_to_string(), libMesh::FEType::family, libMesh::HERMITE, libMesh::HIERARCHIC, libMesh::HIERARCHIC_VEC, libMesh::L2_HIERARCHIC, libMesh::L2_HIERARCHIC_VEC, libMesh::L2_LAGRANGE, libMesh::L2_LAGRANGE_VEC, libMesh::L2_RAVIART_THOMAS, libMesh::LAGRANGE, libMesh::LAGRANGE_VEC, libMesh::MONOMIAL, libMesh::MONOMIAL_VEC, libMesh::NEDELEC_ONE, libMesh::RATIONAL_BERNSTEIN, libMesh::RAVIART_THOMAS, libMesh::SCALAR, libMesh::SIDE_HIERARCHIC, libMesh::SUBDIVISION, libMesh::SZABAB, and libMesh::XYZ.
Referenced by FETest< order, family, elem_type >::testFEInterface().
|
staticprivate |
true
if et
is an element to be processed by class InfFE
, false otherwise or when !LIBMESH_ENABLE_INFINITE_ELEMENTS. Definition at line 52 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(), shape(), shape_deriv(), shape_deriv_function(), shape_function(), shape_second_deriv(), shape_second_deriv_function(), and side_nodal_soln().
|
static |
This is now deprecated; use FEMap::map instead.
Definition at line 677 of file fe_interface.C.
References dim, ifem_map(), is_InfFE_elem(), and libMesh::Elem::type().
Definition at line 2131 of file fe_interface.C.
References libMesh::BERNSTEIN, libMesh::C0POLYGON, libMesh::C0POLYHEDRON, libMesh::CLOUGH, libMesh::EDGE2, libMesh::EDGE3, libMesh::EDGE4, libMesh::FEType::family, libMesh::HERMITE, libMesh::HEX20, libMesh::HEX27, libMesh::HEX8, libMesh::HIERARCHIC, libMesh::HIERARCHIC_VEC, libMesh::L2_HIERARCHIC, libMesh::L2_HIERARCHIC_VEC, libMesh::L2_LAGRANGE, libMesh::L2_LAGRANGE_VEC, libMesh::L2_RAVIART_THOMAS, libMesh::LAGRANGE, libMesh::LAGRANGE_VEC, libMesh::MONOMIAL, libMesh::MONOMIAL_VEC, libMesh::NEDELEC_ONE, libMesh::PRISM15, libMesh::PRISM18, libMesh::PRISM20, libMesh::PRISM21, libMesh::PRISM6, libMesh::PYRAMID13, libMesh::PYRAMID14, libMesh::PYRAMID18, libMesh::PYRAMID5, libMesh::QUAD4, libMesh::QUAD8, libMesh::QUAD9, libMesh::QUADSHELL4, libMesh::QUADSHELL8, libMesh::QUADSHELL9, libMesh::RATIONAL_BERNSTEIN, libMesh::RAVIART_THOMAS, libMesh::SIDE_HIERARCHIC, libMesh::SUBDIVISION, libMesh::SZABAB, libMesh::TET10, libMesh::TET14, libMesh::TET4, libMesh::TRI3, libMesh::TRI3SUBDIVISION, libMesh::TRI6, libMesh::TRI7, libMesh::TRISHELL3, and libMesh::XYZ.
Referenced by libMesh::FEAbstract::compute_node_constraints(), and libMesh::DofMap::reinit().
|
static |
Definition at line 355 of file fe_interface.C.
References dim, ifem_n_dofs(), is_InfFE_elem(), and libMesh::FEType::order.
Referenced by libMesh::DofMap::_dof_indices(), assemble_func(), assemble_SchroedingerEquation(), assemble_wave(), 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 >::compute_shape_functions(), libMesh::InfFE< Dim, T_radial, T_map >::inf_compute_constraints(), libMesh::InfFE< Dim, T_radial, T_map >::init_shape_functions(), main(), n_dofs(), libMesh::InfFE< Dim, T_radial, T_map >::n_dofs(), n_shape_functions(), and libMesh::HPCoarsenTest::select_refinement().
|
static |
Similar to the function above but takes an Elem * and accounts for p-refinement internally, if any.
This function is designed to prevent users from needing to trick FEInterface::n_dofs() into giving them the right number of dofs when working with p-refined elements. See, e.g. FEInterface::compute_data().
Definition at line 377 of file fe_interface.C.
References n_dofs(), libMesh::FEType::order, and libMesh::Elem::p_level().
|
static |
elem
for finite element type fe_t
The p_level() of elem
is accounted for internally by increasing the Order of the passed-in FEType if add_p_level
is true.
Definition at line 390 of file fe_interface.C.
References dim, libMesh::Elem::dim(), ifem_n_dofs(), is_InfFE_elem(), n_dofs(), libMesh::FEType::order, libMesh::Elem::p_level(), and libMesh::Elem::type().
|
static |
elem
for finite element type fe_t
Definition at line 412 of file fe_interface.C.
References dim, libMesh::Elem::dim(), ifem_n_dofs(), is_InfFE_elem(), n_dofs(), libMesh::FEType::order, and libMesh::Elem::type().
|
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 436 of file fe_interface.C.
References dim, ifem_n_dofs_at_node(), is_InfFE_elem(), and libMesh::FEType::order.
Referenced by libMesh::DofMap::_dof_indices(), libMesh::DofMap::_node_dof_indices(), libMesh::FEGenericBase< FEOutputType< T >::type >::coarsened_dof_values(), libMesh::InfFE< Dim, T_radial, T_map >::compute_shape_indices(), libMesh::DofMap::constrain_p_dofs(), n_dofs_at_node(), 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 >::SortAndCopy::operator()(), libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::ProjectVertices::operator()(), and libMesh::DofMap::reinit().
|
static |
fe_t
. Accounts for Elem::p_level() internally if add_p_level
is true. Definition at line 482 of file fe_interface.C.
References dim, libMesh::Elem::dim(), ifem_n_dofs_at_node(), is_InfFE_elem(), n_dofs_at_node(), libMesh::FEType::order, libMesh::Elem::p_level(), and libMesh::Elem::type().
|
static |
fe_t
. Ignores Elem::p_level() and computes a total Order given by fe_t.order + extra_order when determining the number of DOFs. Definition at line 506 of file fe_interface.C.
References dim, libMesh::Elem::dim(), ifem_n_dofs_at_node(), is_InfFE_elem(), n_dofs_at_node(), libMesh::FEType::order, and libMesh::Elem::type().
|
static |
The behavior is otherwise exactly the same, since this function does not depend on the Elem::p_level().
Definition at line 458 of file fe_interface.C.
References n_dofs_at_node().
Referenced by libMesh::DofMap::_dof_indices(), and libMesh::DofMap::old_dof_indices().
|
static |
Definition at line 470 of file fe_interface.C.
References dim, libMesh::Elem::dim(), and n_dofs_at_node().
|
static |
Definition at line 530 of file fe_interface.C.
References dim, 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::OldSolutionValue< Output, point_output >::eval_old_dofs(), n_dofs_per_elem(), libMesh::InfFE< Dim, T_radial, T_map >::n_dofs_per_elem(), libMesh::DofMap::old_dof_indices(), libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::SortAndCopy::operator()(), and libMesh::DofMap::reinit().
|
static |
On a p-refined element, fe_t.order
should be the total order of the element.
Definition at line 552 of file fe_interface.C.
References dim, libMesh::Elem::dim(), ifem_n_dofs_per_elem(), is_InfFE_elem(), n_dofs_per_elem(), libMesh::FEType::order, libMesh::Elem::p_level(), and libMesh::Elem::type().
|
static |
Same thing but internally elem->p_level() is ignored and extra_order is used instead.
Definition at line 575 of file fe_interface.C.
References dim, libMesh::Elem::dim(), ifem_n_dofs_per_elem(), is_InfFE_elem(), n_dofs_per_elem(), libMesh::FEType::order, and libMesh::Elem::type().
|
static |
Definition at line 272 of file fe_interface.C.
References dim, ifem_n_shape_functions(), is_InfFE_elem(), and libMesh::FEType::order.
Referenced by libMesh::FEMap::compute_face_map(), libMesh::FEMap::init_edge_shape_functions(), libMesh::FEMap::init_face_shape_functions(), libMesh::FEMap::init_reference_to_physical_map(), libMesh::LIBMESH_DEFAULT_VECTORIZED_FE(), libMesh::FEMap::map(), libMesh::FEMap::map_deriv(), and FETest< order, family, elem_type >::testFEInterface().
|
static |
elem
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.
Definition at line 299 of file fe_interface.C.
References dim, libMesh::Elem::dim(), ifem_n_shape_functions(), is_InfFE_elem(), n_dofs(), libMesh::FEType::order, libMesh::Elem::p_level(), and libMesh::Elem::type().
|
static |
Same as above, but ignores the elem->p_level() and uses the specified extra_order instead.
Definition at line 327 of file fe_interface.C.
References dim, libMesh::Elem::dim(), ifem_n_shape_functions(), is_InfFE_elem(), n_dofs(), libMesh::FEType::order, and libMesh::Elem::type().
|
static |
Definition at line 2677 of file fe_interface.C.
References libMesh::FEType::family, field_type(), mesh, and libMesh::TYPE_VECTOR.
Referenced by libMesh::ExactSolution::_compute_error(), libMesh::EquationSystems::build_discontinuous_solution_vector(), libMesh::EquationSystems::build_parallel_solution_vector(), libMesh::EquationSystems::build_variable_names(), and libMesh::Variable::n_components().
|
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.
fe_t.order
should be the base order of the element. The Elem::p_level(), if any, is accounted for internally by this routine.dim
argument. Definition at line 625 of file fe_interface.C.
References dim, 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(), side_nodal_soln(), libMesh::Nemesis_IO_Helper::write_nodal_solution(), 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, \( \xi \le 1 \) becomes \( \xi \le 1 + \epsilon \).
Elem::on_reference_element()
instead. Definition at line 751 of file fe_interface.C.
References libMesh::FEAbstract::on_reference_element().
|
static |
Definition at line 2658 of file fe_interface.C.
References libMesh::BERNSTEIN, libMesh::HIERARCHIC, libMesh::HIERARCHIC_VEC, libMesh::L2_HIERARCHIC, libMesh::L2_HIERARCHIC_VEC, libMesh::NEDELEC_ONE, libMesh::RATIONAL_BERNSTEIN, libMesh::RAVIART_THOMAS, and libMesh::SZABAB.
|
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 760 of file fe_interface.C.
References dim, 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(), ifem_shape(), libMesh::InfFE< Dim, T_radial, T_map >::inf_compute_constraints(), libMesh::LIBMESH_DEFAULT_VECTORIZED_FE(), main(), libMesh::HCurlFETransformation< OutputShape >::map_phi(), libMesh::HDivFETransformation< OutputShape >::map_phi(), libMesh::InfFE< Dim, T_radial, T_map >::shape(), shape(), libMesh::InfFE< Dim, T_radial, T_map >::shape_deriv(), shape_function(), and NavierSystem::side_constraint().
|
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 780 of file fe_interface.C.
References ifem_shape(), is_InfFE_elem(), libMesh::FEType::order, shape(), and libMesh::Elem::type().
|
static |
p
.Non-deprecated version of the shape() function. The Elem::p_level() is accounted for internally if add_p_level
Definition at line 804 of file fe_interface.C.
References dim, libMesh::Elem::dim(), ifem_shape(), is_InfFE_elem(), libMesh::FEType::order, shape(), and libMesh::Elem::type().
|
static |
p
.Non-deprecated version of the shape() function. The Elem::p_level() is ignored and extra_order is used instead.
Definition at line 831 of file fe_interface.C.
References dim, libMesh::Elem::dim(), 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.
|
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.Non-deprecated version of templated shape() function that accounts for Elem::p_level() internally.
|
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.Non-deprecated version of templated shape() function that ignores Elem::p_level() and instead uses extra_order internally.
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 863 of file fe_interface.C.
References dim.
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 908 of file fe_interface.C.
References dim.
void libMesh::FEInterface::shape | ( | const FEType & | fe_t, |
const Elem * | elem, | ||
const unsigned int | i, | ||
const Point & | p, | ||
Real & | phi | ||
) |
Definition at line 953 of file fe_interface.C.
References dim.
void libMesh::FEInterface::shape | ( | const FEType & | fe_t, |
int | extra_order, | ||
const Elem * | elem, | ||
const unsigned int | i, | ||
const Point & | p, | ||
Real & | phi | ||
) |
Definition at line 996 of file fe_interface.C.
References dim.
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 1136 of file fe_interface.C.
References dim.
void libMesh::FEInterface::shape | ( | const FEType & | fe_t, |
const Elem * | elem, | ||
const unsigned int | i, | ||
const Point & | p, | ||
RealGradient & | phi | ||
) |
Definition at line 1177 of file fe_interface.C.
References dim.
void libMesh::FEInterface::shape | ( | const FEType & | fe_t, |
int | extra_order, | ||
const Elem * | elem, | ||
const unsigned int | i, | ||
const Point & | p, | ||
RealGradient & | phi | ||
) |
Definition at line 1222 of file fe_interface.C.
References dim.
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 1851 of file fe_interface.C.
References dim.
|
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.Definition at line 1380 of file fe_interface.C.
References dim, ifem_shape_deriv(), is_InfFE_elem(), and libMesh::FEType::order.
Referenced by libMesh::InfFE< Dim, T_radial, T_map >::compute_data(), compute_data(), ifem_shape_deriv(), libMesh::InfFE< Dim, T_radial, T_map >::shape_deriv(), shape_deriv(), and shape_deriv_function().
|
static |
Definition at line 1404 of file fe_interface.C.
References dim, ifem_shape_deriv(), libMesh::Elem::infinite(), libMesh::FEType::order, and shape_deriv().
|
static |
p
. This method allows you to specify the dimension, element, 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. Definition at line 1443 of file fe_interface.C.
References dim, libMesh::Elem::dim(), ifem_shape_deriv(), is_InfFE_elem(), libMesh::FEType::order, shape_deriv(), and libMesh::Elem::type().
|
static |
Non-deprecated version of function above.
Definition at line 1471 of file fe_interface.C.
References dim, libMesh::Elem::dim(), ifem_shape_deriv(), libMesh::Elem::infinite(), libMesh::FEType::order, and shape_deriv().
|
static |
Referenced by libMesh::FEMap::init_edge_shape_functions(), libMesh::FEMap::init_face_shape_functions(), libMesh::FEMap::init_reference_to_physical_map(), and libMesh::FEMap::map_deriv().
|
static |
Non-deprecated version of the function above.
Definition at line 1650 of file fe_interface.C.
References dim, libMesh::Elem::dim(), is_InfFE_elem(), shape_deriv(), and libMesh::Elem::type().
|
static |
Fills dphi
with the derivatives of the \( i^{th} \) shape function at point p
in direction j
.
void libMesh::FEInterface::shape_derivs | ( | const FEType & | fe_t, |
const Elem * | elem, | ||
const unsigned int | i, | ||
const unsigned int | j, | ||
const std::vector< Point > & | p, | ||
std::vector< Real > & | dphi, | ||
const bool | add_p_level | ||
) |
Definition at line 1516 of file fe_interface.C.
void libMesh::FEInterface::shape_derivs | ( | const FEType & | fe_t, |
const Elem * | elem, | ||
const unsigned int | i, | ||
const unsigned int | j, | ||
const std::vector< Point > & | p, | ||
std::vector< RealGradient > & | dphi, | ||
const bool | add_p_level | ||
) |
Definition at line 1596 of file fe_interface.C.
|
static |
Referenced by libMesh::FEMap::init_edge_shape_functions(), libMesh::FEMap::init_face_shape_functions(), and libMesh::FEMap::map().
|
static |
Definition at line 1364 of file fe_interface.C.
References dim, libMesh::Elem::dim(), is_InfFE_elem(), 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.dim
and takes an Elem * instead. Definition at line 1668 of file fe_interface.C.
References dim, is_InfFE_elem(), and libMesh::FEType::order.
Referenced by shape_second_deriv(), and shape_second_deriv_function().
|
static |
dim
and takes an Elem * instead.Definition at line 1702 of file fe_interface.C.
References dim, libMesh::Elem::infinite(), libMesh::FEType::order, and shape_second_deriv().
|
static |
p
.fe_t.order
should be the total order of the element. Definition at line 1737 of file fe_interface.C.
References dim, libMesh::Elem::dim(), is_InfFE_elem(), libMesh::FEType::order, shape_second_deriv(), and libMesh::Elem::type().
|
static |
Non-deprecated version of function above taking an extra_order
parameter.
Definition at line 1774 of file fe_interface.C.
References dim, libMesh::Elem::dim(), libMesh::Elem::infinite(), libMesh::FEType::order, and shape_second_deriv().
|
static |
Referenced by libMesh::FEMap::init_edge_shape_functions(), libMesh::FEMap::init_face_shape_functions(), and libMesh::FEMap::init_reference_to_physical_map().
|
static |
Definition at line 1834 of file fe_interface.C.
References dim, libMesh::Elem::dim(), is_InfFE_elem(), shape_second_deriv(), and libMesh::Elem::type().
|
static |
Fills phi
with the values of the \( i^{th} \) shape function at point p
.
This method allows you to specify the dimension, element type, and order directly.
true
for add_p_level
if you want the Elem::p_level() to be accounted for internally, pass false if you want fe_t.order to be used instead.dim
as a parameter. This is a relatively large changeset with little benefit if we go the deprecation route, so it would probably be cleaner to just break backwards compatibility... these functions seem to mainly be used internally by the library and changing them is unlikely to break application codes. void libMesh::FEInterface::shapes | ( | const unsigned int | dim, |
const FEType & | fe_t, | ||
const Elem * | elem, | ||
const unsigned int | i, | ||
const std::vector< Point > & | p, | ||
std::vector< Real > & | phi, | ||
const bool | add_p_level | ||
) |
Definition at line 1047 of file fe_interface.C.
References dim, libMesh::index_range(), and libMesh::FEType::order.
void libMesh::FEInterface::shapes | ( | const unsigned int | dim, |
const FEType & | fe_t, | ||
const Elem * | elem, | ||
const unsigned int | i, | ||
const std::vector< Point > & | p, | ||
std::vector< RealGradient > & | phi, | ||
const bool | add_p_level | ||
) |
Definition at line 1270 of file fe_interface.C.
References dim.
|
static |
Build the nodal soln on one side from the (full) element soln.
This is the solution that will be plotted on side-elements.
fe_t.order
should be the base order of the element. The Elem::p_level(), if any, is accounted for internally by this routine. Definition at line 650 of file fe_interface.C.
References dim, libMesh::Elem::dim(), is_InfFE_elem(), nodal_soln(), libMesh::FEType::order, and libMesh::Elem::type().
Referenced by libMesh::EquationSystems::build_discontinuous_solution_vector(), and libMesh::EquationSystems::build_parallel_solution_vector().