$extrastylesheet
#include <fe_interface.h>
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 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_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) |
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.
| libMesh::FEInterface::FEInterface | ( | ) | [private] |
Empty constructor. Do not create an object of this type.
Definition at line 32 of file fe_interface.C.
{
libmesh_error_msg("ERROR: Do not define an object of this type.");
}
| virtual libMesh::FEInterface::~FEInterface | ( | ) | [inline, virtual] |
| 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.
Definition at line 840 of file fe_interface.C.
References libMesh::BERNSTEIN, libMesh::CLOUGH, libMesh::Elem::dim(), libMesh::FEType::family, libMesh::HERMITE, libMesh::HIERARCHIC, libMesh::L2_HIERARCHIC, libMesh::LAGRANGE, libMesh::LAGRANGE_VEC, libMesh::libmesh_assert(), libMesh::SZABAB, and libMesh::DofMap::variable_type().
{
libmesh_assert(elem);
const FEType& fe_t = dof_map.variable_type(variable_number);
switch (elem->dim())
{
case 0:
case 1:
{
// No constraints in 0D/1D.
return;
}
case 2:
{
switch (fe_t.family)
{
case CLOUGH:
FE<2,CLOUGH>::compute_constraints (constraints,
dof_map,
variable_number,
elem); return;
case HERMITE:
FE<2,HERMITE>::compute_constraints (constraints,
dof_map,
variable_number,
elem); return;
case LAGRANGE:
FE<2,LAGRANGE>::compute_constraints (constraints,
dof_map,
variable_number,
elem); return;
case HIERARCHIC:
FE<2,HIERARCHIC>::compute_constraints (constraints,
dof_map,
variable_number,
elem); return;
case L2_HIERARCHIC:
FE<2,L2_HIERARCHIC>::compute_constraints (constraints,
dof_map,
variable_number,
elem); return;
case LAGRANGE_VEC:
FE<2,LAGRANGE_VEC>::compute_constraints (constraints,
dof_map,
variable_number,
elem); return;
#ifdef LIBMESH_ENABLE_HIGHER_ORDER_SHAPES
case SZABAB:
FE<2,SZABAB>::compute_constraints (constraints,
dof_map,
variable_number,
elem); return;
case BERNSTEIN:
FE<2,BERNSTEIN>::compute_constraints (constraints,
dof_map,
variable_number,
elem); return;
#endif
default:
return;
}
}
case 3:
{
switch (fe_t.family)
{
case HERMITE:
FE<3,HERMITE>::compute_constraints (constraints,
dof_map,
variable_number,
elem); return;
case LAGRANGE:
FE<3,LAGRANGE>::compute_constraints (constraints,
dof_map,
variable_number,
elem); return;
case HIERARCHIC:
FE<3,HIERARCHIC>::compute_constraints (constraints,
dof_map,
variable_number,
elem); return;
case L2_HIERARCHIC:
FE<3,L2_HIERARCHIC>::compute_constraints (constraints,
dof_map,
variable_number,
elem); return;
case LAGRANGE_VEC:
FE<3,LAGRANGE_VEC>::compute_constraints (constraints,
dof_map,
variable_number,
elem); return;
#ifdef LIBMESH_ENABLE_HIGHER_ORDER_SHAPES
case SZABAB:
FE<3,SZABAB>::compute_constraints (constraints,
dof_map,
variable_number,
elem); return;
case BERNSTEIN:
FE<3,BERNSTEIN>::compute_constraints (constraints,
dof_map,
variable_number,
elem); return;
#endif
default:
return;
}
}
default:
libmesh_error_msg("Invalid dimension = " << elem->dim());
}
}
| 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 returns the values also through data. See this as a generalization of shape(). Currently, with disabled infinite elements, returns a vector of all shape functions of elem evaluated ap p.
On a p-refined element, fe_t.order should be the base order of the element.
Definition at line 804 of file fe_interface.C.
References ifem_compute_data(), libMesh::FEComputeData::init(), is_InfFE_elem(), n_dofs(), libMesh::FEType::order, libMesh::FEComputeData::p, libMesh::Elem::p_level(), libMesh::FEComputeData::shape, shape(), and libMesh::Elem::type().
Referenced by ifem_compute_data(), and libMesh::MeshFunction::operator()().
{
#ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
if ( elem && is_InfFE_elem(elem->type()) )
{
data.init();
ifem_compute_data(dim, fe_t, elem, data);
return;
}
#endif
FEType p_refined = fe_t;
p_refined.order = static_cast<Order>(p_refined.order + elem->p_level());
const unsigned int n_dof = n_dofs (dim, p_refined, elem->type());
const Point& p = data.p;
data.shape.resize(n_dof);
// set default values for all the output fields
data.init();
for (unsigned int n=0; n<n_dof; n++)
data.shape[n] = shape(dim, p_refined, elem, n, p);
return;
}
| 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.
Definition at line 984 of file fe_interface.C.
{
// No element-specific optimizations currently exist
FEBase::compute_periodic_constraints (constraints,
dof_map,
boundaries,
mesh,
point_locator,
variable_number,
elem);
}
| 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.
Definition at line 497 of file fe_interface.C.
References libMesh::FEType::order.
Referenced by libMesh::FEGenericBase< OutputType >::coarsened_dof_values(), libMesh::ProjectSolution::operator()(), libMesh::ProjectFEMSolution::operator()(), and libMesh::BoundaryProjectSolution::operator()().
{
const Order o = fe_t.order;
void_fe_with_vec_switch(dofs_on_edge(elem, o, e, di));
libmesh_error_msg("We'll never get here!");
}
| 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.
Definition at line 482 of file fe_interface.C.
References libMesh::FEType::order.
Referenced by libMesh::FEGenericBase< OutputType >::coarsened_dof_values(), libMesh::FEGenericBase< OutputType >::compute_periodic_constraints(), libMesh::FEGenericBase< OutputType >::compute_proj_constraints(), libMesh::ProjectSolution::operator()(), libMesh::ProjectFEMSolution::operator()(), and libMesh::BoundaryProjectSolution::operator()().
{
const Order o = fe_t.order;
void_fe_with_vec_switch(dofs_on_side(elem, o, s, di));
libmesh_error_msg("We'll never get here!");
}
| 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.
Definition at line 1351 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().
{
switch (fe_t.family)
{
case LAGRANGE:
case L2_LAGRANGE:
case MONOMIAL:
case L2_HIERARCHIC:
case XYZ:
case SUBDIVISION:
case LAGRANGE_VEC:
case NEDELEC_ONE:
return false;
case CLOUGH:
case HERMITE:
case HIERARCHIC:
#ifdef LIBMESH_ENABLE_HIGHER_ORDER_SHAPES
case BERNSTEIN:
case SZABAB:
#endif
default:
return true;
}
}
| 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.
Definition at line 1376 of file fe_interface.C.
References libMesh::FEType::family.
Referenced by libMesh::ExactSolution::_compute_error(), libMesh::EquationSystems::build_solution_vector(), libMesh::EquationSystems::build_variable_names(), libMesh::ExactSolution::compute_error(), libMesh::ExactSolution::error_norm(), libMesh::FEMSystem::init_context(), and libMesh::FE< Dim, T >::init_shape_functions().
{
return FEInterface::field_type( fe_type.family );
}
| 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.
Definition at line 1381 of file fe_interface.C.
References libMesh::LAGRANGE_VEC, libMesh::NEDELEC_ONE, libMesh::TYPE_SCALAR, and libMesh::TYPE_VECTOR.
{
switch (fe_family)
{
case LAGRANGE_VEC:
case NEDELEC_ONE:
return TYPE_VECTOR;
default:
return TYPE_SCALAR;
}
}
| void libMesh::FEInterface::ifem_compute_data | ( | const unsigned int | dim, |
| const FEType & | fe_t, | ||
| const Elem * | elem, | ||
| FEComputeData & | data | ||
| ) | [static, private] |
Definition at line 802 of file fe_interface_inf_fe.C.
References compute_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().
{
switch (dim)
{
// 1D
case 1:
{
switch (fe_t.radial_family)
{
/*
* For no derivatives (and local coordinates, as
* given in \p p) the infinite element shapes
* are independent of mapping type
*/
case INFINITE_MAP:
InfFE<1,INFINITE_MAP,CARTESIAN>::compute_data(fe_t, elem, data);
break;
case JACOBI_20_00:
InfFE<1,JACOBI_20_00,CARTESIAN>::compute_data(fe_t, elem, data);
break;
case JACOBI_30_00:
InfFE<1,JACOBI_30_00,CARTESIAN>::compute_data(fe_t, elem, data);
break;
case LEGENDRE:
InfFE<1,LEGENDRE,CARTESIAN>::compute_data(fe_t, elem, data);
break;
case LAGRANGE:
InfFE<1,LAGRANGE,CARTESIAN>::compute_data(fe_t, elem, data);
break;
default:
libmesh_error_msg("Invalid radial family = " << fe_t.radial_family);
}
break;
}
// 2D
case 2:
{
switch (fe_t.radial_family)
{
case INFINITE_MAP:
InfFE<2,INFINITE_MAP,CARTESIAN>::compute_data(fe_t, elem, data);
break;
case JACOBI_20_00:
InfFE<2,JACOBI_20_00,CARTESIAN>::compute_data(fe_t, elem, data);
break;
case JACOBI_30_00:
InfFE<2,JACOBI_30_00,CARTESIAN>::compute_data(fe_t, elem, data);
break;
case LEGENDRE:
InfFE<2,LEGENDRE,CARTESIAN>::compute_data(fe_t, elem, data);
break;
case LAGRANGE:
InfFE<2,LAGRANGE,CARTESIAN>::compute_data(fe_t, elem, data);
break;
default:
libmesh_error_msg("Invalid radial family = " << fe_t.radial_family);
}
break;
}
// 3D
case 3:
{
switch (fe_t.radial_family)
{
case INFINITE_MAP:
InfFE<3,INFINITE_MAP,CARTESIAN>::compute_data(fe_t, elem, data);
break;
case JACOBI_20_00:
InfFE<3,JACOBI_20_00,CARTESIAN>::compute_data(fe_t, elem, data);
break;
case JACOBI_30_00:
InfFE<3,JACOBI_30_00,CARTESIAN>::compute_data(fe_t, elem, data);
break;
case LEGENDRE:
InfFE<3,LEGENDRE,CARTESIAN>::compute_data(fe_t, elem, data);
break;
case LAGRANGE:
InfFE<3,LAGRANGE,CARTESIAN>::compute_data(fe_t, elem, data);
break;
default:
libmesh_error_msg("Invalid radial family = " << fe_t.radial_family);
}
break;
}
default:
libmesh_error_msg("Invalid dim = " << dim);
break;
}
}
| 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 |
||
| ) | [static, private] |
Definition at line 431 of file fe_interface_inf_fe.C.
References libMesh::CARTESIAN, libMesh::ELLIPSOIDAL, libMesh::FEType::inf_map, inverse_map(), and libMesh::SPHERICAL.
Referenced by inverse_map().
{
switch (dim)
{
// 1D
case 1:
{
switch (fe_t.inf_map)
{
case CARTESIAN:
return InfFE<1,JACOBI_20_00,CARTESIAN>::inverse_map(elem, p, tolerance, secure);
case SPHERICAL:
case ELLIPSOIDAL:
libmesh_not_implemented_msg("ERROR: Spherical and Ellipsoidal IFEMs not (yet) implemented.");
/*
case SPHERICAL:
return InfFE<1,JACOBI_20_00,SPHERICAL>::inverse_map(elem, p, tolerance);
case ELLIPSOIDAL:
return InfFE<1,JACOBI_20_00,ELLIPSOIDAL>::inverse_map(elem, p, tolerance);
*/
default:
libmesh_error_msg("Invalid map = " << fe_t.inf_map);
}
}
// 2D
case 2:
{
switch (fe_t.inf_map)
{
case CARTESIAN:
return InfFE<2,JACOBI_20_00,CARTESIAN>::inverse_map(elem, p, tolerance, secure);
case SPHERICAL:
case ELLIPSOIDAL:
libmesh_not_implemented_msg("ERROR: Spherical and Ellipsoidal IFEMs not (yet) implemented.");
/*
case SPHERICAL:
return InfFE<2,JACOBI_20_00,SPHERICAL>::inverse_map(elem, p, tolerance);
case ELLIPSOIDAL:
return InfFE<2,JACOBI_20_00,ELLIPSOIDAL>::inverse_map(elem, p, tolerance);
*/
default:
libmesh_error_msg("Invalid map = " << fe_t.inf_map);
}
}
// 3D
case 3:
{
switch (fe_t.inf_map)
{
case CARTESIAN:
return InfFE<3,JACOBI_20_00,CARTESIAN>::inverse_map(elem, p, tolerance, secure);
case SPHERICAL:
case ELLIPSOIDAL:
libmesh_not_implemented_msg("ERROR: Spherical and Ellipsoidal IFEMs not (yet) implemented.");
/*
case SPHERICAL:
return InfFE<3,JACOBI_20_00,SPHERICAL>::inverse_map(elem, p, tolerance);
case ELLIPSOIDAL:
return InfFE<3,JACOBI_20_00,ELLIPSOIDAL>::inverse_map(elem, p, tolerance);
*/
default:
libmesh_error_msg("Invalid map = " << fe_t.inf_map);
}
}
default:
libmesh_error_msg("Invalid dim = " << dim);
}
libmesh_error_msg("We'll never get here!");
Point pt;
return pt;
}
| 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 |
||
| ) | [static, private] |
Definition at line 528 of file fe_interface_inf_fe.C.
References libMesh::CARTESIAN, libMesh::FEType::inf_map, and inverse_map().
{
switch (dim)
{
// 1D
case 1:
{
switch (fe_t.inf_map)
{
case CARTESIAN:
InfFE<1,JACOBI_20_00,CARTESIAN>::inverse_map(elem, physical_points, reference_points, tolerance, secure);
return;
default:
libmesh_error_msg("Invalid map = " << fe_t.inf_map);
}
}
// 2D
case 2:
{
switch (fe_t.inf_map)
{
case CARTESIAN:
InfFE<2,JACOBI_20_00,CARTESIAN>::inverse_map(elem, physical_points, reference_points, tolerance, secure);
return;
default:
libmesh_error_msg("Invalid map = " << fe_t.inf_map);
}
}
// 3D
case 3:
{
switch (fe_t.inf_map)
{
case CARTESIAN:
InfFE<3,JACOBI_20_00,CARTESIAN>::inverse_map(elem, physical_points, reference_points, tolerance, secure);
return;
default:
libmesh_error_msg("Invalid map = " << fe_t.inf_map);
}
}
default:
libmesh_error_msg("Invalid dim = " << dim);
}
}
| unsigned int libMesh::FEInterface::ifem_n_dofs | ( | const unsigned int | dim, |
| const FEType & | fe_t, | ||
| const ElemType | t | ||
| ) | [static, private] |
Definition at line 74 of file fe_interface_inf_fe.C.
References n_dofs().
Referenced by n_dofs().
{
switch (dim)
{
// 1D
case 1:
/*
* Since InfFE<Dim,T_radial,T_map>::n_dofs(...)
* is actually independent of T_radial and T_map, we can use
* just any T_radial and T_map
*/
return InfFE<1,JACOBI_20_00,CARTESIAN>::n_dofs(fe_t, t);
// 2D
case 2:
return InfFE<2,JACOBI_20_00,CARTESIAN>::n_dofs(fe_t, t);
// 3D
case 3:
return InfFE<3,JACOBI_20_00,CARTESIAN>::n_dofs(fe_t, t);
default:
libmesh_error_msg("Unsupported dim = " << dim);
}
libmesh_error_msg("We'll never get here!");
return 0;
}
| unsigned int libMesh::FEInterface::ifem_n_dofs_at_node | ( | const unsigned int | dim, |
| const FEType & | fe_t, | ||
| const ElemType | t, | ||
| const unsigned int | n | ||
| ) | [static, private] |
Definition at line 108 of file fe_interface_inf_fe.C.
References n_dofs_at_node().
Referenced by n_dofs_at_node().
{
switch (dim)
{
// 1D
case 1:
/*
* Since InfFE<Dim,T_radial,T_map>::n_dofs_at_node(...)
* is actually independent of T_radial and T_map, we can use
* just any T_radial and T_map
*/
return InfFE<1,JACOBI_20_00,CARTESIAN>::n_dofs_at_node(fe_t, t, n);
// 2D
case 2:
return InfFE<2,JACOBI_20_00,CARTESIAN>::n_dofs_at_node(fe_t, t, n);
// 3D
case 3:
return InfFE<3,JACOBI_20_00,CARTESIAN>::n_dofs_at_node(fe_t, t, n);
default:
libmesh_error_msg("Unsupported dim = " << dim);
}
libmesh_error_msg("We'll never get here!");
return 0;
}
| unsigned int libMesh::FEInterface::ifem_n_dofs_per_elem | ( | const unsigned int | dim, |
| const FEType & | fe_t, | ||
| const ElemType | t | ||
| ) | [static, private] |
Definition at line 144 of file fe_interface_inf_fe.C.
References n_dofs_per_elem().
Referenced by n_dofs_per_elem().
{
switch (dim)
{
// 1D
case 1:
/*
* Since InfFE<Dim,T_radial,T_map>::n_dofs(...)
* is actually independent of T_radial and T_map, we can use
* just any T_radial and T_map
*/
return InfFE<1,JACOBI_20_00,CARTESIAN>::n_dofs_per_elem(fe_t, t);
// 2D
case 2:
return InfFE<2,JACOBI_20_00,CARTESIAN>::n_dofs_per_elem(fe_t, t);
// 3D
case 3:
return InfFE<3,JACOBI_20_00,CARTESIAN>::n_dofs_per_elem(fe_t, t);
default:
libmesh_error_msg("Unsupported dim = " << dim);
}
libmesh_error_msg("We'll never get here!");
return 0;
}
| unsigned int libMesh::FEInterface::ifem_n_shape_functions | ( | const unsigned int | dim, |
| const FEType & | fe_t, | ||
| const ElemType | t | ||
| ) | [static, private] |
Definition at line 39 of file fe_interface_inf_fe.C.
References n_shape_functions().
Referenced by n_shape_functions().
{
switch (dim)
{
// 1D
case 1:
/*
* Since InfFE<Dim,T_radial,T_map>::n_shape_functions(...)
* is actually independent of T_radial and T_map, we can use
* just any T_radial and T_map
*/
return InfFE<1,JACOBI_20_00,CARTESIAN>::n_shape_functions(fe_t, t);
// 2D
case 2:
return InfFE<2,JACOBI_20_00,CARTESIAN>::n_shape_functions(fe_t, t);
// 3D
case 3:
return InfFE<3,JACOBI_20_00,CARTESIAN>::n_shape_functions(fe_t, t);
default:
libmesh_error_msg("Unsupported dim = " << dim);
}
libmesh_error_msg("We'll never get here!");
return 0;
}
| 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 | ||
| ) | [static, private] |
Definition at line 178 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(), and libMesh::FEType::radial_family.
Referenced by nodal_soln().
{
switch (dim)
{
// 1D
case 1:
{
switch (fe_t.radial_family)
{
case INFINITE_MAP:
libmesh_error_msg("ERROR: INFINTE_MAP is not a valid shape family for radial approximation.");
case JACOBI_20_00:
{
switch (fe_t.inf_map)
{
case CARTESIAN:
{
InfFE<1,JACOBI_20_00,CARTESIAN>::nodal_soln(fe_t, elem, elem_soln, nodal_soln);
break;
}
default:
libmesh_error_msg("ERROR: Spherical & Ellipsoidal IFEMs not implemented.");
}
break;
}
case JACOBI_30_00:
{
switch (fe_t.inf_map)
{
case CARTESIAN:
{
InfFE<1,JACOBI_30_00,CARTESIAN>::nodal_soln(fe_t, elem, elem_soln, nodal_soln);
break;
}
default:
libmesh_error_msg("ERROR: Spherical & Ellipsoidal IFEMs not implemented.");
}
break;
}
case LEGENDRE:
{
switch (fe_t.inf_map)
{
case CARTESIAN:
{
InfFE<1,LEGENDRE,CARTESIAN>::nodal_soln(fe_t, elem, elem_soln, nodal_soln);
break;
}
default:
libmesh_error_msg("ERROR: Spherical & Ellipsoidal IFEMs not implemented.");
}
break;
}
case LAGRANGE:
{
switch (fe_t.inf_map)
{
case CARTESIAN:
{
InfFE<1,LAGRANGE,CARTESIAN>::nodal_soln(fe_t, elem, elem_soln, nodal_soln);
break;
}
default:
libmesh_error_msg("ERROR: Spherical & Ellipsoidal IFEMs not implemented.");
}
break;
}
default:
libmesh_error_msg("ERROR: Bad FEType.radial_family= " << fe_t.radial_family);
}
break;
}
// 2D
case 2:
{
switch (fe_t.radial_family)
{
case INFINITE_MAP:
libmesh_error_msg("ERROR: INFINTE_MAP is not a valid shape family for radial approximation.");
case JACOBI_20_00:
{
switch (fe_t.inf_map)
{
case CARTESIAN:
{
InfFE<2,JACOBI_20_00,CARTESIAN>::nodal_soln(fe_t, elem, elem_soln, nodal_soln);
break;
}
default:
libmesh_error_msg("ERROR: Spherical & Ellipsoidal IFEMs not implemented.");
}
break;
}
case JACOBI_30_00:
{
switch (fe_t.inf_map)
{
case CARTESIAN:
{
InfFE<2,JACOBI_30_00,CARTESIAN>::nodal_soln(fe_t, elem, elem_soln, nodal_soln);
break;
}
default:
libmesh_error_msg("ERROR: Spherical & Ellipsoidal IFEMs not implemented.");
}
break;
}
case LEGENDRE:
{
switch (fe_t.inf_map)
{
case CARTESIAN:
{
InfFE<2,LEGENDRE,CARTESIAN>::nodal_soln(fe_t, elem, elem_soln, nodal_soln);
break;
}
default:
libmesh_error_msg("ERROR: Spherical & Ellipsoidal IFEMs not implemented.");
}
break;
}
case LAGRANGE:
{
switch (fe_t.inf_map)
{
case CARTESIAN:
{
InfFE<2,LAGRANGE,CARTESIAN>::nodal_soln(fe_t, elem, elem_soln, nodal_soln);
break;
}
default:
libmesh_error_msg("ERROR: Spherical & Ellipsoidal IFEMs not implemented.");
}
break;
}
default:
libmesh_error_msg("ERROR: Bad FEType.radial_family= " << fe_t.radial_family);
}
break;
}
// 3D
case 3:
{
switch (fe_t.radial_family)
{
case INFINITE_MAP:
libmesh_error_msg("ERROR: INFINTE_MAP is not a valid shape family for radial approximation.");
case JACOBI_20_00:
{
switch (fe_t.inf_map)
{
case CARTESIAN:
{
InfFE<3,JACOBI_20_00,CARTESIAN>::nodal_soln(fe_t, elem, elem_soln, nodal_soln);
break;
}
default:
libmesh_error_msg("ERROR: Spherical & Ellipsoidal IFEMs not implemented.");
}
break;
}
case JACOBI_30_00:
{
switch (fe_t.inf_map)
{
case CARTESIAN:
{
InfFE<3,JACOBI_30_00,CARTESIAN>::nodal_soln(fe_t, elem, elem_soln, nodal_soln);
break;
}
default:
libmesh_error_msg("ERROR: Spherical & Ellipsoidal IFEMs not implemented.");
}
break;
}
case LEGENDRE:
{
switch (fe_t.inf_map)
{
case CARTESIAN:
{
InfFE<3,LEGENDRE,CARTESIAN>::nodal_soln(fe_t, elem, elem_soln, nodal_soln);
break;
}
default:
libmesh_error_msg("ERROR: Spherical & Ellipsoidal IFEMs not implemented.");
}
break;
}
case LAGRANGE:
{
switch (fe_t.inf_map)
{
case CARTESIAN:
{
InfFE<3,LAGRANGE,CARTESIAN>::nodal_soln(fe_t, elem, elem_soln, nodal_soln);
break;
}
default:
libmesh_error_msg("ERROR: Spherical & Ellipsoidal IFEMs not implemented.");
}
break;
}
default:
libmesh_error_msg("ERROR: Bad FEType.radial_family= " << fe_t.radial_family);
}
break;
}
default:
libmesh_error_msg("Invalid dim = " << dim);
}
}
| bool libMesh::FEInterface::ifem_on_reference_element | ( | const Point & | p, |
| const ElemType | t, | ||
| const Real | eps | ||
| ) | [static, private] |
Definition at line 590 of file fe_interface_inf_fe.C.
References on_reference_element().
{
return FEBase::on_reference_element(p,t,eps);
}
| Real libMesh::FEInterface::ifem_shape | ( | const unsigned int | dim, |
| const FEType & | fe_t, | ||
| const ElemType | t, | ||
| const unsigned int | i, | ||
| const Point & | p | ||
| ) | [static, private] |
Definition at line 600 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 shape().
Referenced by shape().
{
switch (dim)
{
// 1D
case 1:
{
switch (fe_t.radial_family)
{
/*
* For no derivatives (and local coordinates, as
* given in \p p) the infinite element shapes
* are independent of mapping type
*/
case INFINITE_MAP:
return InfFE<1,INFINITE_MAP,CARTESIAN>::shape(fe_t, t, i, p);
case JACOBI_20_00:
return InfFE<1,JACOBI_20_00,CARTESIAN>::shape(fe_t, t, i, p);
case JACOBI_30_00:
return InfFE<1,JACOBI_30_00,CARTESIAN>::shape(fe_t, t, i, p);
case LEGENDRE:
return InfFE<1,LEGENDRE,CARTESIAN>::shape(fe_t, t, i, p);
case LAGRANGE:
return InfFE<1,LAGRANGE,CARTESIAN>::shape(fe_t, t, i, p);
default:
libmesh_error_msg("Invalid radial family = " << fe_t.radial_family);
}
}
// 2D
case 2:
{
switch (fe_t.radial_family)
{
case INFINITE_MAP:
return InfFE<2,INFINITE_MAP,CARTESIAN>::shape(fe_t, t, i, p);
case JACOBI_20_00:
return InfFE<2,JACOBI_20_00,CARTESIAN>::shape(fe_t, t, i, p);
case JACOBI_30_00:
return InfFE<2,JACOBI_30_00,CARTESIAN>::shape(fe_t, t, i, p);
case LEGENDRE:
return InfFE<2,LEGENDRE,CARTESIAN>::shape(fe_t, t, i, p);
case LAGRANGE:
return InfFE<2,LAGRANGE,CARTESIAN>::shape(fe_t, t, i, p);
default:
libmesh_error_msg("Invalid radial family = " << fe_t.radial_family);
}
}
// 3D
case 3:
{
switch (fe_t.radial_family)
{
case INFINITE_MAP:
return InfFE<3,INFINITE_MAP,CARTESIAN>::shape(fe_t, t, i, p);
case JACOBI_20_00:
return InfFE<3,JACOBI_20_00,CARTESIAN>::shape(fe_t, t, i, p);
case JACOBI_30_00:
return InfFE<3,JACOBI_30_00,CARTESIAN>::shape(fe_t, t, i, p);
case LEGENDRE:
return InfFE<3,LEGENDRE,CARTESIAN>::shape(fe_t, t, i, p);
case LAGRANGE:
return InfFE<3,LAGRANGE,CARTESIAN>::shape(fe_t, t, i, p);
default:
libmesh_error_msg("Invalid radial family = " << fe_t.radial_family);
}
}
default:
libmesh_error_msg("Invalid dim = " << dim);
}
libmesh_error_msg("We'll never get here!");
return 0.;
}
| Real libMesh::FEInterface::ifem_shape | ( | const unsigned int | dim, |
| const FEType & | fe_t, | ||
| const Elem * | elem, | ||
| const unsigned int | i, | ||
| const Point & | p | ||
| ) | [static, private] |
Definition at line 701 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 shape().
{
switch (dim)
{
// 1D
case 1:
{
switch (fe_t.radial_family)
{
/*
* For no derivatives (and local coordinates, as
* given in \p p) the infinite element shapes
* are independent of mapping type
*/
case INFINITE_MAP:
return InfFE<1,INFINITE_MAP,CARTESIAN>::shape(fe_t, elem, i, p);
case JACOBI_20_00:
return InfFE<1,JACOBI_20_00,CARTESIAN>::shape(fe_t, elem, i, p);
case JACOBI_30_00:
return InfFE<1,JACOBI_30_00,CARTESIAN>::shape(fe_t, elem, i, p);
case LEGENDRE:
return InfFE<1,LEGENDRE,CARTESIAN>::shape(fe_t, elem, i, p);
case LAGRANGE:
return InfFE<1,LAGRANGE,CARTESIAN>::shape(fe_t, elem, i, p);
default:
libmesh_error_msg("Invalid radial family = " << fe_t.radial_family);
}
}
// 2D
case 2:
{
switch (fe_t.radial_family)
{
case INFINITE_MAP:
return InfFE<2,INFINITE_MAP,CARTESIAN>::shape(fe_t, elem, i, p);
case JACOBI_20_00:
return InfFE<2,JACOBI_20_00,CARTESIAN>::shape(fe_t, elem, i, p);
case JACOBI_30_00:
return InfFE<2,JACOBI_30_00,CARTESIAN>::shape(fe_t, elem, i, p);
case LEGENDRE:
return InfFE<2,LEGENDRE,CARTESIAN>::shape(fe_t, elem, i, p);
case LAGRANGE:
return InfFE<2,LAGRANGE,CARTESIAN>::shape(fe_t, elem, i, p);
default:
libmesh_error_msg("Invalid radial family = " << fe_t.radial_family);
}
}
// 3D
case 3:
{
switch (fe_t.radial_family)
{
case INFINITE_MAP:
return InfFE<3,INFINITE_MAP,CARTESIAN>::shape(fe_t, elem, i, p);
case JACOBI_20_00:
return InfFE<3,JACOBI_20_00,CARTESIAN>::shape(fe_t, elem, i, p);
case JACOBI_30_00:
return InfFE<3,JACOBI_30_00,CARTESIAN>::shape(fe_t, elem, i, p);
case LEGENDRE:
return InfFE<3,LEGENDRE,CARTESIAN>::shape(fe_t, elem, i, p);
case LAGRANGE:
return InfFE<3,LAGRANGE,CARTESIAN>::shape(fe_t, elem, i, p);
default:
libmesh_error_msg("Invalid radial family = " << fe_t.radial_family);
}
}
default:
libmesh_error_msg("Invalid dim = " << dim);
}
libmesh_error_msg("We'll never get here!");
return 0.;
}
| 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] |
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 552 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< OutputType >::coarsened_dof_values(), libMesh::FEAbstract::compute_node_constraints(), libMesh::FEGenericBase< OutputType >::compute_periodic_constraints(), libMesh::FEAbstract::compute_periodic_node_constraints(), libMesh::FEGenericBase< OutputType >::compute_proj_constraints(), libMesh::InfQuad4::contains_point(), libMesh::InfPrism6::contains_point(), libMesh::InfHex8::contains_point(), libMesh::MeshFunction::gradient(), libMesh::MeshFunction::hessian(), ifem_inverse_map(), inverse_map(), libMesh::DGFEMContext::neighbor_side_fe_reinit(), 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().
{
#ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
if ( is_InfFE_elem(elem->type()) )
return ifem_inverse_map(dim, fe_t, elem, p,tolerance, secure);
#endif
fe_with_vec_switch(inverse_map(elem, p, tolerance, secure));
libmesh_error_msg("We'll never get here!");
return Point();
}
| 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] |
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 575 of file fe_interface.C.
References libMesh::err, ifem_inverse_map(), inverse_map(), is_InfFE_elem(), libMesh::TypeVector< T >::size(), and libMesh::Elem::type().
{
const std::size_t n_pts = physical_points.size();
// Resize the vector
reference_points.resize(n_pts);
if (n_pts == 0)
{
libMesh::err << "WARNING: empty vector physical_points!"
<< std::endl;
libmesh_here();
return;
}
#ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
if ( is_InfFE_elem(elem->type()) )
{
ifem_inverse_map(dim, fe_t, elem, physical_points, reference_points, tolerance, secure);
return;
// libmesh_not_implemented();
}
#endif
void_fe_with_vec_switch(inverse_map(elem, physical_points, reference_points, tolerance, secure));
libmesh_error_msg("We'll never get here!");
}
| bool libMesh::FEInterface::is_InfFE_elem | ( | const ElemType | et | ) | [inline, static, private] |
et is an element to be processed by class InfFE. Otherwise, it returns false. For compatibility with disabled infinite elements it always returns false. Definition at line 449 of file fe_interface.h.
Referenced by compute_data(), inverse_map(), n_dofs(), n_dofs_at_node(), n_dofs_per_elem(), n_shape_functions(), nodal_soln(), and shape().
{
return false;
}
| Point libMesh::FEInterface::map | ( | unsigned int | dim, |
| const FEType & | fe_t, | ||
| const Elem * | elem, | ||
| const Point & | p | ||
| ) | [static] |
Returns the point in physical space of the reference point refpoint which is passed in.
Definition at line 537 of file fe_interface.C.
Referenced by libMesh::Elem::point_test().
{
fe_with_vec_switch(map(elem, p));
libmesh_error_msg("We'll never get here!");
return Point();
}
| 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.
Definition at line 1006 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::SUBDIVISION, libMesh::SZABAB, libMesh::TET10, libMesh::TET4, libMesh::TRI3, libMesh::TRI3SUBDIVISION, libMesh::TRI6, and libMesh::XYZ.
Referenced by libMesh::DofMap::reinit().
{
// Yeah, I know, infinity is much larger than 11, but our
// solvers don't seem to like high degree polynomials, and our
// quadrature rules and number_lookups tables
// need to go up higher.
const unsigned int unlimited = 11;
// If we used 0 as a default, then elements missing from this
// table (e.g. infinite elements) would be considered broken.
const unsigned int unknown = unlimited;
switch (fe_t.family)
{
case LAGRANGE:
case L2_LAGRANGE: // TODO: L2_LAGRANGE can have higher "max_order" than LAGRANGE
case LAGRANGE_VEC:
switch (el_t)
{
case EDGE2:
case EDGE3:
case EDGE4:
return 3;
case TRI3:
return 1;
case TRI6:
return 2;
case QUAD4:
return 1;
case QUAD8:
case QUAD9:
return 2;
case TET4:
return 1;
case TET10:
return 2;
case HEX8:
return 1;
case HEX20:
case HEX27:
return 2;
case PRISM6:
return 1;
case PRISM15:
case PRISM18:
return 2;
case PYRAMID5:
return 1;
case PYRAMID13:
case PYRAMID14:
return 2;
default:
return unknown;
}
break;
case MONOMIAL:
switch (el_t)
{
case EDGE2:
case EDGE3:
case EDGE4:
case TRI3:
case TRI6:
case QUAD4:
case QUAD8:
case QUAD9:
case TET4:
case TET10:
case HEX8:
case HEX20:
case HEX27:
case PRISM6:
case PRISM15:
case PRISM18:
case PYRAMID5:
case PYRAMID13:
case PYRAMID14:
return unlimited;
default:
return unknown;
}
break;
#ifdef LIBMESH_ENABLE_HIGHER_ORDER_SHAPES
case BERNSTEIN:
switch (el_t)
{
case EDGE2:
case EDGE3:
case EDGE4:
return unlimited;
case TRI3:
return 0;
case TRI6:
return 6;
case QUAD4:
return 0;
case QUAD8:
case QUAD9:
return unlimited;
case TET4:
return 1;
case TET10:
return 2;
case HEX8:
return 0;
case HEX20:
return 2;
case HEX27:
return 4;
case PRISM6:
case PRISM15:
case PRISM18:
case PYRAMID5:
case PYRAMID13:
case PYRAMID14:
return 0;
default:
return unknown;
}
break;
case SZABAB:
switch (el_t)
{
case EDGE2:
case EDGE3:
case EDGE4:
return 7;
case TRI3:
return 0;
case TRI6:
return 7;
case QUAD4:
return 0;
case QUAD8:
case QUAD9:
return 7;
case TET4:
case TET10:
case HEX8:
case HEX20:
case HEX27:
case PRISM6:
case PRISM15:
case PRISM18:
case PYRAMID5:
case PYRAMID13:
case PYRAMID14:
return 0;
default:
return unknown;
}
break;
#endif
case XYZ:
switch (el_t)
{
case EDGE2:
case EDGE3:
case EDGE4:
case TRI3:
case TRI6:
case QUAD4:
case QUAD8:
case QUAD9:
case TET4:
case TET10:
case HEX8:
case HEX20:
case HEX27:
case PRISM6:
case PRISM15:
case PRISM18:
case PYRAMID5:
case PYRAMID13:
case PYRAMID14:
return unlimited;
default:
return unknown;
}
break;
case CLOUGH:
switch (el_t)
{
case EDGE2:
case EDGE3:
return 3;
case EDGE4:
case TRI3:
return 0;
case TRI6:
return 3;
case QUAD4:
case QUAD8:
case QUAD9:
case TET4:
case TET10:
case HEX8:
case HEX20:
case HEX27:
case PRISM6:
case PRISM15:
case PRISM18:
case PYRAMID5:
case PYRAMID13:
case PYRAMID14:
return 0;
default:
return unknown;
}
break;
case HERMITE:
switch (el_t)
{
case EDGE2:
case EDGE3:
return unlimited;
case EDGE4:
case TRI3:
case TRI6:
return 0;
case QUAD4:
return 3;
case QUAD8:
case QUAD9:
return unlimited;
case TET4:
case TET10:
return 0;
case HEX8:
return 3;
case HEX20:
case HEX27:
return unlimited;
case PRISM6:
case PRISM15:
case PRISM18:
case PYRAMID5:
case PYRAMID13:
case PYRAMID14:
return 0;
default:
return unknown;
}
break;
case HIERARCHIC:
switch (el_t)
{
case EDGE2:
case EDGE3:
case EDGE4:
return unlimited;
case TRI3:
return 1;
case TRI6:
return unlimited;
case QUAD4:
return 1;
case QUAD8:
case QUAD9:
return unlimited;
case TET4:
case TET10:
return 0;
case HEX8:
case HEX20:
return 1;
case HEX27:
return unlimited;
case PRISM6:
case PRISM15:
case PRISM18:
case PYRAMID5:
case PYRAMID13:
case PYRAMID14:
return 0;
default:
return unknown;
}
break;
case L2_HIERARCHIC:
switch (el_t)
{
case EDGE2:
case EDGE3:
case EDGE4:
return unlimited;
case TRI3:
return 1;
case TRI6:
return unlimited;
case QUAD4:
return 1;
case QUAD8:
case QUAD9:
return unlimited;
case TET4:
case TET10:
return 0;
case HEX8:
case HEX20:
return 1;
case HEX27:
return unlimited;
case PRISM6:
case PRISM15:
case PRISM18:
case PYRAMID5:
case PYRAMID13:
case PYRAMID14:
return 0;
default:
return unknown;
}
break;
case SUBDIVISION:
switch (el_t)
{
case TRI3SUBDIVISION:
return unlimited;
default:
return unknown;
}
break;
case NEDELEC_ONE:
switch (el_t)
{
case TRI6:
case QUAD8:
case QUAD9:
case HEX20:
case HEX27:
return 1;
default:
return 0;
}
break;
default:
return 0;
break;
}
}
| unsigned int libMesh::FEInterface::n_dofs | ( | const unsigned int | dim, |
| const FEType & | fe_t, | ||
| const ElemType | t | ||
| ) | [static] |
On a p-refined element, fe_t.order should be the total order of the element.
Definition at line 414 of file fe_interface.C.
References ifem_n_dofs(), is_InfFE_elem(), and libMesh::FEType::order.
Referenced by libMesh::FEGenericBase< OutputType >::coarsened_dof_values(), compute_data(), libMesh::FEAbstract::compute_node_constraints(), libMesh::FEAbstract::compute_periodic_node_constraints(), ifem_n_dofs(), libMesh::InfFE< Dim, T_radial, T_map >::n_dofs(), and libMesh::HPCoarsenTest::select_refinement().
{
#ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
if ( is_InfFE_elem(t) )
return ifem_n_dofs(dim, fe_t, t);
#endif
const Order o = fe_t.order;
fe_with_vec_switch(n_dofs(t, o));
libmesh_error_msg("We'll never get here!");
return 0;
}
| unsigned int libMesh::FEInterface::n_dofs_at_node | ( | const unsigned int | dim, |
| const FEType & | fe_t, | ||
| const ElemType | t, | ||
| const unsigned int | n | ||
| ) | [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 ifem_n_dofs_at_node(), is_InfFE_elem(), and libMesh::FEType::order.
Referenced by libMesh::DofMap::_dof_indices(), libMesh::FEGenericBase< OutputType >::coarsened_dof_values(), libMesh::InfFE< Dim, T_radial, T_map >::compute_shape_indices(), ifem_n_dofs_at_node(), libMesh::InfFE< Dim, T_radial, T_map >::n_dofs_at_node(), libMesh::DofMap::old_dof_indices(), libMesh::ProjectSolution::operator()(), libMesh::ProjectFEMSolution::operator()(), libMesh::BoundaryProjectSolution::operator()(), and libMesh::DofMap::reinit().
{
#ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
if ( is_InfFE_elem(t) )
return ifem_n_dofs_at_node(dim, fe_t, t, n);
#endif
const Order o = fe_t.order;
fe_with_vec_switch(n_dofs_at_node(t, o, n));
libmesh_error_msg("We'll never get here!");
return 0;
}
| unsigned int libMesh::FEInterface::n_dofs_per_elem | ( | const unsigned int | dim, |
| const FEType & | fe_t, | ||
| const ElemType | t | ||
| ) | [static] |
On a p-refined element, fe_t.order should be the total order of the element.
Definition at line 460 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(), ifem_n_dofs_per_elem(), libMesh::InfFE< Dim, T_radial, T_map >::n_dofs_per_elem(), libMesh::DofMap::old_dof_indices(), and libMesh::DofMap::reinit().
{
#ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
if ( is_InfFE_elem(t) )
return ifem_n_dofs_per_elem(dim, fe_t, t);
#endif
const Order o = fe_t.order;
fe_with_vec_switch(n_dofs_per_elem(t, o));
libmesh_error_msg("We'll never get here!");
return 0;
}
| unsigned int libMesh::FEInterface::n_shape_functions | ( | const unsigned int | dim, |
| const FEType & | fe_t, | ||
| const ElemType | t | ||
| ) | [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 385 of file fe_interface.C.
References ifem_n_shape_functions(), is_InfFE_elem(), and libMesh::FEType::order.
Referenced by ifem_n_shape_functions().
{
#ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
/*
* Since the FEType, stored in DofMap/(some System child), has to
* be the _same_ for InfFE and FE, we have to catch calls
* to infinite elements through the element type.
*/
if ( is_InfFE_elem(t) )
return ifem_n_shape_functions(dim, fe_t, t);
#endif
const Order o = fe_t.order;
fe_with_vec_switch(n_shape_functions(t, o));
libmesh_error_msg("We'll never get here!");
return 0;
}
| 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.
Definition at line 1393 of file fe_interface.C.
References libMesh::FEType::family, libMesh::LAGRANGE_VEC, libMesh::MeshBase::mesh_dimension(), and libMesh::NEDELEC_ONE.
Referenced by libMesh::ExactSolution::_compute_error(), libMesh::EquationSystems::build_solution_vector(), and libMesh::EquationSystems::build_variable_names().
{
switch (fe_type.family)
{
//FIXME: We currently assume that the number of vector components is tied
// to the mesh dimension. This will break for mixed-dimension meshes.
case LAGRANGE_VEC:
case NEDELEC_ONE:
return mesh.mesh_dimension();
default:
return 1;
}
}
| 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.
Definition at line 513 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_solution_vector(), ifem_nodal_soln(), libMesh::EnsightIO::write_scalar_ascii(), and libMesh::EnsightIO::write_vector_ascii().
{
#ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
if ( is_InfFE_elem(elem->type()) )
{
ifem_nodal_soln(dim, fe_t, elem, elem_soln, nodal_soln);
return;
}
#endif
const Order order = fe_t.order;
void_fe_with_vec_switch(nodal_soln(elem, order, elem_soln, nodal_soln));
}
| bool libMesh::FEInterface::on_reference_element | ( | const Point & | p, |
| const ElemType | t, | ||
| const Real | eps = TOLERANCE |
||
| ) | [static] |
Since we are doing floating point comparisons here the parameter eps can be specified to indicate a tolerance. For example,
becomes
.
Definition at line 614 of file fe_interface.C.
Referenced by libMesh::InfQuad4::contains_point(), libMesh::InfPrism6::contains_point(), libMesh::InfHex8::contains_point(), ifem_on_reference_element(), and libMesh::Elem::point_test().
{
return FEBase::on_reference_element(p,t,eps);
}
| Real libMesh::FEInterface::shape | ( | const unsigned int | dim, |
| const FEType & | fe_t, | ||
| const ElemType | t, | ||
| const unsigned int | i, | ||
| const Point & | p | ||
| ) | [static] |
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.On a p-refined element, fe_t.order should be the total order of the element.
Definition at line 624 of file fe_interface.C.
References ifem_shape(), is_InfFE_elem(), and libMesh::FEType::order.
Referenced by compute_data(), libMesh::InfFE< Dim, T_radial, T_map >::compute_data(), libMesh::FEAbstract::compute_node_constraints(), libMesh::FEAbstract::compute_periodic_node_constraints(), ifem_shape(), shape(), and libMesh::InfFE< Dim, T_radial, T_map >::shape().
{
#ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
if ( is_InfFE_elem(t) )
return ifem_shape(dim, fe_t, t, i, p);
#endif
const Order o = fe_t.order;
fe_switch(shape(t,o,i,p));
libmesh_error_msg("We'll never get here!");
return 0.;
}
| Real libMesh::FEInterface::shape | ( | const unsigned int | dim, |
| const FEType & | fe_t, | ||
| const Elem * | elem, | ||
| const unsigned int | i, | ||
| const Point & | p | ||
| ) | [static] |
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.On a p-refined element, fe_t.order should be the base order of the element.
Definition at line 645 of file fe_interface.C.
References ifem_shape(), is_InfFE_elem(), libMesh::FEType::order, shape(), and libMesh::Elem::type().
{
#ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
if ( elem && is_InfFE_elem(elem->type()) )
return ifem_shape(dim, fe_t, elem, i, p);
#endif
const Order o = fe_t.order;
fe_switch(shape(elem,o,i,p));
libmesh_error_msg("We'll never get here!");
return 0.;
}
| 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] |
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.On a p-refined element, fe_t.order should be the total order of the element.
| 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] |
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.On a p-refined element, 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 667 of file fe_interface.C.
{
#ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
if ( is_InfFE_elem(t) )
phi = ifem_shape(dim, fe_t, t, i, p);
#endif
const Order o = fe_t.order;
switch(dim)
{
case 0:
fe_scalar_vec_error_switch(0, shape(t,o,i,p), phi = , ; break;);
break;
case 1:
fe_scalar_vec_error_switch(1, shape(t,o,i,p), phi = , ; break;);
break;
case 2:
fe_scalar_vec_error_switch(2, shape(t,o,i,p), phi = , ; break;);
break;
case 3:
fe_scalar_vec_error_switch(3, shape(t,o,i,p), phi = , ; break;);
break;
default:
libmesh_error_msg("Invalid dimension = " << dim);
}
return;
}
| 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 705 of file fe_interface.C.
{
#ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
if ( elem && is_InfFE_elem(elem->type()) )
phi = ifem_shape(dim, fe_t, elem, i, p);
#endif
const Order o = fe_t.order;
switch(dim)
{
case 0:
fe_scalar_vec_error_switch(0, shape(elem,o,i,p), phi = , ; break;);
break;
case 1:
fe_scalar_vec_error_switch(1, shape(elem,o,i,p), phi = , ; break;);
break;
case 2:
fe_scalar_vec_error_switch(2, shape(elem,o,i,p), phi = , ; break;);
break;
case 3:
fe_scalar_vec_error_switch(3, shape(elem,o,i,p), phi = , ; break;);
break;
default:
libmesh_error_msg("Invalid dimension = " << dim);
}
return;
}
| 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 743 of file fe_interface.C.
{
const Order o = fe_t.order;
switch(dim)
{
case 0:
fe_vector_scalar_error_switch(0, shape(t,o,i,p), phi = , ; break;);
break;
case 1:
fe_vector_scalar_error_switch(1, shape(t,o,i,p), phi = , ; break;);
break;
case 2:
fe_vector_scalar_error_switch(2, shape(t,o,i,p), phi = , ; break;);
break;
case 3:
fe_vector_scalar_error_switch(3, shape(t,o,i,p), phi = , ; break;);
break;
default:
libmesh_error_msg("Invalid dimension = " << dim);
}
return;
}
| 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 774 of file fe_interface.C.
{
const Order o = fe_t.order;
switch(dim)
{
case 0:
fe_vector_scalar_error_switch(0, shape(elem,o,i,p), phi = , ; break;);
break;
case 1:
fe_vector_scalar_error_switch(1, shape(elem,o,i,p), phi = , ; break;);
break;
case 2:
fe_vector_scalar_error_switch(2, shape(elem,o,i,p), phi = , ; break;);
break;
case 3:
fe_vector_scalar_error_switch(3, shape(elem,o,i,p), phi = , ; break;);
break;
default:
libmesh_error_msg("Invalid dimension = " << dim);
}
return;
}