libMesh::FEType Class Reference

Manages the family, order, etc. parameters for a given FE. More...

#include <fe_type.h>

Public Member Functions

 FEType (const int o=1, const FEFamily f=LAGRANGE)
 
 FEType (const int o=1, const FEFamily f=LAGRANGE, const int ro=THIRD, const FEFamily rf=JACOBI_20_00, const InfMapType im=CARTESIAN)
 
bool operator== (const FEType &f2) const
 
bool operator!= (const FEType &f2) const
 
bool operator< (const FEType &f2) const
 
Order default_quadrature_order () const
 
UniquePtr< QBasedefault_quadrature_rule (const unsigned int dim, const int extraorder=0) const
 

Public Attributes

OrderWrapper order
 
FEFamily family
 
OrderWrapper radial_order
 
FEFamily radial_family
 
InfMapType inf_map
 

Detailed Description

Manages the family, order, etc. parameters for a given FE.

class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialized finite element families.

Author
Benjamin S. Kirk
Date
2002

Definition at line 178 of file fe_type.h.

Constructor & Destructor Documentation

libMesh::FEType::FEType ( const int  o = 1,
const FEFamily  f = LAGRANGE 
)
inline

Constructor. Optionally takes the approximation Order and the finite element family FEFamily

Definition at line 188 of file fe_type.h.

189  :
190  order(o),
191  family(f)
192  {}
FEFamily family
Definition: fe_type.h:206
OrderWrapper order
Definition: fe_type.h:200
libMesh::FEType::FEType ( const int  o = 1,
const FEFamily  f = LAGRANGE,
const int  ro = THIRD,
const FEFamily  rf = JACOBI_20_00,
const InfMapType  im = CARTESIAN 
)
inline

Constructor. Optionally takes the approximation Order and the finite element family FEFamily. Note that for non-infinite elements the order and base order are the same, as with the family and base_family. It must be so, otherwise what we switch on would change when infinite elements are not compiled in.

Definition at line 218 of file fe_type.h.

222  :
223  order(o),
224  radial_order(ro),
225  family(f),
226  radial_family(rf),
227  inf_map(im)
228  {}
FEFamily family
Definition: fe_type.h:206
OrderWrapper radial_order
Definition: fe_type.h:238
OrderWrapper order
Definition: fe_type.h:200
InfMapType inf_map
Definition: fe_type.h:259
FEFamily radial_family
Definition: fe_type.h:251

Member Function Documentation

Order libMesh::FEType::default_quadrature_order ( ) const
inline
Returns
the default quadrature order for this FEType. The default quadrature order is calculated assuming a polynomial of degree order and is based on integrating the mass matrix for such an element exactly.

Definition at line 334 of file fe_type.h.

Referenced by libMesh::FEGenericBase< OutputType >::compute_periodic_constraints(), libMesh::FEGenericBase< OutputType >::compute_proj_constraints(), default_quadrature_rule(), and libMesh::Elem::volume().

335 {
336  return static_cast<Order>(2*static_cast<unsigned int>(order.get_order()) + 1);
337 }
int get_order() const
Definition: fe_type.h:77
OrderWrapper order
Definition: fe_type.h:200
UniquePtr< QBase > libMesh::FEType::default_quadrature_rule ( const unsigned int  dim,
const int  extraorder = 0 
) const
Returns
a quadrature rule of appropriate type and order for this FEType. The default quadrature rule is based on integrating the mass matrix for such an element exactly. Higher or lower degree rules can be chosen by changing the extraorder parameter.

Definition at line 31 of file fe_type.C.

References libMesh::CLOUGH, default_quadrature_order(), family, std::max(), and libMesh::SUBDIVISION.

Referenced by libMesh::ExactSolution::_compute_error(), libMesh::UniformRefinementEstimator::_estimate_error(), libMesh::System::calculate_norm(), libMesh::FEGenericBase< OutputType >::coarsened_dof_values(), libMesh::ExactErrorEstimator::estimate_error(), libMesh::FEMContext::FEMContext(), libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::operator()(), libMesh::WeightedPatchRecoveryErrorEstimator::EstimateError::operator()(), libMesh::PatchRecoveryErrorEstimator::EstimateError::operator()(), libMesh::BoundaryProjectSolution::operator()(), and libMesh::HPCoarsenTest::select_refinement().

33 {
34 
35  // Clough elements have at least piecewise cubic functions
36  if (family == CLOUGH)
37  {
38  // this seems ridiculous but for some reason gcc 3.3.5 wants
39  // this when using complex numbers (spetersen 04/20/06)
40  const unsigned int seven = 7;
41 
42  return UniquePtr<QBase>
43  (new QClough(dim,
44  static_cast<Order>
45  (std::max(static_cast<unsigned int>
46  (this->default_quadrature_order()), seven + extraorder))));
47  }
48 
49  if (family == SUBDIVISION)
50  return UniquePtr<QBase>
51  (new QGauss(dim, static_cast<Order>(1 + extraorder)));
52 
53  return UniquePtr<QBase>
54  (new QGauss(dim, static_cast<Order>(this->default_quadrature_order()
55  + extraorder)));
56 }
FEFamily family
Definition: fe_type.h:206
long double max(long double a, double b)
Order default_quadrature_order() const
Definition: fe_type.h:334
bool libMesh::FEType::operator!= ( const FEType f2) const
inline

Tests inequality

Definition at line 281 of file fe_type.h.

282  {
283  return !(*this == f2);
284  }
bool libMesh::FEType::operator< ( const FEType f2) const
inline

An ordering to make FEType useful as a std::map key

Definition at line 289 of file fe_type.h.

References family, inf_map, order, radial_family, and radial_order.

290  {
291  if (order != f2.order)
292  return (order < f2.order);
293  if (family != f2.family)
294  return (family < f2.family);
295 
296 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
297  if (radial_order != f2.radial_order)
298  return (radial_order < f2.radial_order);
299  if (radial_family != f2.radial_family)
300  return (radial_family < f2.radial_family);
301  if (inf_map != f2.inf_map)
302  return (inf_map < f2.inf_map);
303 #endif // ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
304  return false;
305  }
FEFamily family
Definition: fe_type.h:206
OrderWrapper radial_order
Definition: fe_type.h:238
OrderWrapper order
Definition: fe_type.h:200
InfMapType inf_map
Definition: fe_type.h:259
FEFamily radial_family
Definition: fe_type.h:251
bool libMesh::FEType::operator== ( const FEType f2) const
inline

Tests equality

Definition at line 266 of file fe_type.h.

References family, inf_map, order, radial_family, and radial_order.

267  {
268  return (order == f2.order
269  && family == f2.family
270 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
271  && radial_order == f2.radial_order
272  && radial_family == f2.radial_family
273  && inf_map == f2.inf_map
274 #endif // ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
275  );
276  }
FEFamily family
Definition: fe_type.h:206
OrderWrapper radial_order
Definition: fe_type.h:238
OrderWrapper order
Definition: fe_type.h:200
InfMapType inf_map
Definition: fe_type.h:259
FEFamily radial_family
Definition: fe_type.h:251

Member Data Documentation

FEFamily libMesh::FEType::family

The type of finite element. Valid types are LAGRANGE, HIERARCHIC, etc...

The type of approximation in radial direction. Valid types are JACOBI_20_00, JACOBI_30_00, etc...

Definition at line 206 of file fe_type.h.

Referenced by libMesh::DofMap::_dof_indices(), libMesh::FEMSystem::assembly(), libMesh::FEMap::build(), libMesh::FETransformationBase< OutputShape >::build(), libMesh::FEAbstract::build(), libMesh::FEGenericBase< OutputType >::build(), libMesh::FEMContext::build_new_fe(), libMesh::WrappedFunction< Output >::component(), libMesh::FEInterface::compute_constraints(), libMesh::GMVIO::copy_nodal_solution(), default_quadrature_rule(), libMesh::DofMap::distribute_dofs(), libMesh::DofMap::distribute_local_dofs_node_major(), libMesh::DofMap::distribute_local_dofs_var_major(), libMesh::DofMap::dof_indices(), libMesh::JumpErrorEstimator::estimate_error(), libMesh::FEInterface::extra_hanging_dofs(), libMesh::FEMContext::FEMContext(), libMesh::FEInterface::field_type(), libMesh::FEAbstract::get_family(), libMesh::System::get_info(), libMesh::DifferentiableSystem::have_first_order_scalar_vars(), libMesh::DifferentiableSystem::have_second_order_scalar_vars(), libMesh::FEInterface::max_order(), libMesh::Variable::n_components(), libMesh::FEInterface::n_vec_dim(), libMesh::DofMap::old_dof_indices(), libMesh::WrappedFunction< Output >::operator()(), libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::operator()(), libMesh::BoundaryProjectSolution::operator()(), operator<(), operator==(), libMesh::System::read_header(), libMesh::DofMap::reinit(), and libMesh::System::write_header().

InfMapType libMesh::FEType::inf_map

The coordinate mapping type of the infinite element. When the infinite elements are defined over a surface with a separable coordinate system (sphere, spheroid, ellipsoid), the infinite elements may take advantage of this fact.

Definition at line 259 of file fe_type.h.

Referenced by libMesh::FEGenericBase< OutputType >::build_InfFE(), libMesh::System::get_info(), libMesh::FEInterface::ifem_inverse_map(), libMesh::FEInterface::ifem_map(), libMesh::FEInterface::ifem_nodal_soln(), libMesh::InfFE< Dim, T_radial, T_map >::InfFE(), operator<(), operator==(), libMesh::System::read_header(), and libMesh::System::write_header().

OrderWrapper libMesh::FEType::order

The approximation order of the element.

The approximation order in radial direction of the infinite element.

Definition at line 200 of file fe_type.h.

Referenced by libMesh::DofMap::_dof_indices(), libMesh::FEMContext::build_new_fe(), libMesh::FEGenericBase< OutputType >::coarsened_dof_values(), libMesh::FEInterface::compute_data(), libMesh::FEGenericBase< OutputType >::compute_proj_constraints(), libMesh::DofMap::constrain_p_dofs(), libMesh::GMVIO::copy_nodal_solution(), libMesh::DofMap::distribute_dofs(), libMesh::DofMap::distribute_local_dofs_node_major(), libMesh::DofMap::distribute_local_dofs_var_major(), libMesh::DofMap::dof_indices(), libMesh::FEInterface::dofs_on_edge(), libMesh::FEInterface::dofs_on_side(), libMesh::FEMContext::FEMContext(), libMesh::System::get_info(), libMesh::FEAbstract::get_order(), libMesh::FESubdivision::init_shape_functions(), libMesh::Variable::n_components(), libMesh::FEInterface::n_dofs(), libMesh::FEInterface::n_dofs_at_node(), libMesh::FEInterface::n_dofs_per_elem(), libMesh::FEInterface::n_shape_functions(), libMesh::FEInterface::nodal_soln(), libMesh::DofMap::old_dof_indices(), libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::operator()(), libMesh::WeightedPatchRecoveryErrorEstimator::EstimateError::operator()(), libMesh::PatchRecoveryErrorEstimator::EstimateError::operator()(), operator<(), operator==(), libMesh::System::read_header(), libMesh::DofMap::reinit(), libMesh::DofMap::SCALAR_dof_indices(), libMesh::HPCoarsenTest::select_refinement(), libMesh::FEAbstract::set_fe_order(), and libMesh::FEInterface::shape().

FEFamily libMesh::FEType::radial_family

For InfFE, family contains the radial shape family, while base_family contains the approximation type in circumferential direction. Valid types are LAGRANGE, HIERARCHIC, etc...

Definition at line 251 of file fe_type.h.

Referenced by libMesh::FEGenericBase< OutputType >::build_InfFE(), libMesh::System::get_info(), libMesh::FEInterface::ifem_compute_data(), libMesh::FEInterface::ifem_nodal_soln(), libMesh::FEInterface::ifem_shape(), libMesh::InfFE< Dim, T_radial, T_map >::InfFE(), operator<(), operator==(), libMesh::System::read_header(), and libMesh::System::write_header().


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