Poster of Linux kernelThe best gift for a Linux geek
FEType

FEType

Section: C Library Functions (3) Updated: Thu Apr 7 2011
Local index Up
 

NAME

FEType -  

SYNOPSIS


#include <fe_type.h>  

Public Member Functions


FEType (const Order o=FIRST, const FEFamily f=LAGRANGE)

FEType (const Order o=FIRST, const FEFamily f=LAGRANGE, const Order ro=THIRD, const FEFamily rf=JACOBI_20_00, const InfMapType im=CARTESIAN)

bool operator== (const FEType &f2) const

bool operator< (const FEType &f2) const

Order default_quadrature_order () const

AutoPtr< QBase > default_quadrature_rule (const unsigned int dim, const int extraorder=0) const
 

Public Attributes


Order order

FEFamily family

Order radial_order

FEFamily radial_family

InfMapType inf_map
 

Detailed Description

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

Definition at line 42 of file fe_type.h.  

Constructor & Destructor Documentation

 

FEType::FEType (const Ordero = FIRST, const FEFamilyf = LAGRANGE) [inline]Constructor. Optionally takes the approximation Order and the finite element family FEFamily

Definition at line 52 of file fe_type.h.

                                      :
    order(o),
    family(f)
  {}
 

FEType::FEType (const Ordero = FIRST, const FEFamilyf = LAGRANGE, const Orderro = THIRD, const FEFamilyrf = JACOBI_20_00, const InfMapTypeim = 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 82 of file fe_type.h.

                                          :
    order(o),
    radial_order(ro),
    family(f),
    radial_family(rf),
    inf_map(im)
  {}
 

Member Function Documentation

 

Order 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 191 of file fe_type.h.

References order.

Referenced by FEBase::compute_periodic_constraints(), FEBase::compute_proj_constraints(), default_quadrature_rule(), JumpErrorEstimator::estimate_error(), and Elem::volume().

{
  return static_cast<Order>(2*static_cast<unsigned int>(order) + 1);
}
 

AutoPtr< QBase > FEType::default_quadrature_rule (const unsigned intdim, const intextraorder = 0) constReturns:

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 30 of file fe_type.C.

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

Referenced by ExactSolution::_compute_error(), UniformRefinementEstimator::_estimate_error(), System::calculate_norm(), FEBase::coarsened_dof_values(), ExactErrorEstimator::estimate_error(), FEMContext::FEMContext(), System::ProjectVector::operator()(), PatchRecoveryErrorEstimator::EstimateError::operator()(), and HPCoarsenTest::select_refinement().

{

  // Clough elements have at least piecewise cubic functions
  if (family == CLOUGH)
    {
      // this seems ridiculous but for some reason gcc 3.3.5 wants
      // this when using complex numbers (spetersen 04/20/06)
      const unsigned int seven = 7;

      return AutoPtr<QBase>
        (new QClough(dim,
                     static_cast<Order>
                     (std::max(static_cast<unsigned int>
                               (this->default_quadrature_order()), seven + extraorder))));
    }
  
  return AutoPtr<QBase>
    (new QGauss(dim, static_cast<Order>(this->default_quadrature_order()
                                        + extraorder)));
}
 

bool FEType::operator< (const FEType &f2) const [inline]ordering to make FEType useful as a std::map key

Definition at line 146 of file fe_type.h.

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

  {
    if (order != f2.order)
      return (order < f2.order);
    if (family != f2.family)
      return (family < f2.family);

#ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
    if (radial_order != f2.radial_order)
      return (radial_order < f2.radial_order);
    if (radial_family != f2.radial_family)
      return (radial_family < f2.radial_family);
    if (inf_map != f2.inf_map)
      return (inf_map < f2.inf_map);
#endif // ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
    return false;
  }
 

bool FEType::operator== (const FEType &f2) const [inline]equality

Definition at line 131 of file fe_type.h.

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

  {
    return (order == f2.order
            && family == f2.family
#ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
            && radial_order == f2.radial_order
            && radial_family == f2.radial_family
            && inf_map == f2.inf_map
#endif // ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
            );
  }
 

Member Data Documentation

 

FEFamily FEType::familyThe 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 70 of file fe_type.h.

Referenced by FEBase::build(), FEInterface::compute_constraints(), GMVIO::copy_nodal_solution(), default_quadrature_rule(), DofMap::distribute_local_dofs_var_major(), DofMap::dof_indices(), FEInterface::dofs_on_edge(), FEInterface::dofs_on_side(), FEInterface::extra_hanging_dofs(), FE< Dim, T >::FE(), FEBase::get_family(), System::get_info(), FEInterface::inverse_map(), FEInterface::max_order(), FEInterface::n_dofs(), FEInterface::n_dofs_at_node(), FEInterface::n_dofs_per_elem(), FEInterface::n_shape_functions(), FEInterface::nodal_soln(), DofMap::old_dof_indices(), System::ProjectVector::operator()(), operator<(), operator==(), System::project_vector(), System::read_header(), DofMap::reinit(), DofMap::SCALAR_dof_indices(), FEInterface::shape(), and System::write_header().  

InfMapType FEType::inf_mapThe 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 124 of file fe_type.h.

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

Order FEType::orderThe approximation order of the element.

The approximation order in radial direction of the infinite element.

Definition at line 64 of file fe_type.h.

Referenced by FEBase::coarsened_dof_values(), FEInterface::compute_data(), DofMap::constrain_p_dofs(), GMVIO::copy_nodal_solution(), default_quadrature_order(), DofMap::distribute_local_dofs_node_major(), DofMap::distribute_local_dofs_var_major(), DofMap::dof_indices(), FEInterface::dofs_on_edge(), FEInterface::dofs_on_side(), FEMContext::FEMContext(), System::get_info(), FEBase::get_order(), FEInterface::n_dofs(), FEInterface::n_dofs_at_node(), FEInterface::n_dofs_per_elem(), FEInterface::n_shape_functions(), FEInterface::nodal_soln(), DofMap::old_dof_indices(), System::ProjectVector::operator()(), PatchRecoveryErrorEstimator::EstimateError::operator()(), operator<(), operator==(), System::read_header(), DofMap::reinit(), DofMap::SCALAR_dof_indices(), HPCoarsenTest::select_refinement(), FEInterface::shape(), and System::write_header().  

FEFamily FEType::radial_familyFor 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 116 of file fe_type.h.

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

Order FEType::radial_orderThe approximation order in the base of the infinite element.

Definition at line 103 of file fe_type.h.

Referenced by InfFE< Dim, T_radial, T_map >::compute_data(), InfFE< Dim, T_radial, T_map >::compute_shape_indices(), System::get_info(), InfFE< Dim, T_radial, T_map >::n_dofs(), InfFE< Dim, T_radial, T_map >::n_dofs_at_node(), InfFE< Dim, T_radial, T_map >::n_dofs_per_elem(), operator<(), operator==(), System::read_header(), InfFE< Dim, T_radial, T_map >::reinit(), InfFE< Dim, T_radial, T_map >::shape(), and System::write_header().

 

Author

Generated automatically by Doxygen for libMesh from the source code.


 

Index

NAME
SYNOPSIS
Public Member Functions
Public Attributes
Detailed Description
Constructor & Destructor Documentation
FEType::FEType (const Ordero = FIRST, const FEFamilyf = LAGRANGE) [inline]Constructor. Optionally takes the approximation Order and the finite element family FEFamily
FEType::FEType (const Ordero = FIRST, const FEFamilyf = LAGRANGE, const Orderro = THIRD, const FEFamilyrf = JACOBI_20_00, const InfMapTypeim = 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.
Member Function Documentation
Order FEType::default_quadrature_order () const [inline]Returns:
AutoPtr< QBase > FEType::default_quadrature_rule (const unsigned intdim, const intextraorder = 0) constReturns:
bool FEType::operator< (const FEType &f2) const [inline]ordering to make FEType useful as a std::map key
bool FEType::operator== (const FEType &f2) const [inline]equality
Member Data Documentation
FEFamily FEType::familyThe type of finite element. Valid types are LAGRANGE, HIERARCHIC, etc...
InfMapType FEType::inf_mapThe 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.
Order FEType::orderThe approximation order of the element.
FEFamily FEType::radial_familyFor InfFE, family contains the radial shape family, while base_family contains the approximation type in circumferential direction. Valid types are LAGRANGE, HIERARCHIC, etc...
Order FEType::radial_orderThe approximation order in the base of the infinite element.
Author

This document was created by man2html, using the manual pages.
Time: 21:46:04 GMT, April 16, 2011