libMesh::FEType Class Reference
#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 |
| 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 43 of file fe_type.h.
Constructor & Destructor Documentation
| libMesh::FEType::FEType | ( | const Order | o = FIRST, |
|
| const FEFamily | f = LAGRANGE | |||
| ) | [inline] |
| libMesh::FEType::FEType | ( | const Order | o = FIRST, |
|
| const FEFamily | f = LAGRANGE, |
|||
| const Order | 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 83 of file fe_type.h.
00087 : 00088 order(o), 00089 radial_order(ro), 00090 family(f), 00091 radial_family(rf), 00092 inf_map(im) 00093 {}
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 degreeorderand is based on integrating the mass matrix for such an element exactly.
Definition at line 200 of file fe_type.h.
References order.
Referenced by libMesh::FEGenericBase< OutputType >::compute_periodic_constraints(), libMesh::FEGenericBase< OutputType >::compute_proj_constraints(), default_quadrature_rule(), libMesh::JumpErrorEstimator::estimate_error(), and libMesh::Elem::volume().
| AutoPtr< 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 libMeshEnums::CLOUGH, default_quadrature_order(), family, and std::max().
Referenced by libMesh::ExactSolution::_compute_error(), libMesh::System::calculate_norm(), libMesh::FEGenericBase< OutputType >::coarsened_dof_values(), libMesh::ExactErrorEstimator::estimate_error(), libMesh::FEMContext::FEMContext(), libMesh::WeightedPatchRecoveryErrorEstimator::EstimateError::operator()(), libMesh::BoundaryProjectSolution::operator()(), libMesh::ProjectSolution::operator()(), libMesh::PatchRecoveryErrorEstimator::EstimateError::operator()(), and libMesh::HPCoarsenTest::select_refinement().
00033 { 00034 00035 // Clough elements have at least piecewise cubic functions 00036 if (family == CLOUGH) 00037 { 00038 // this seems ridiculous but for some reason gcc 3.3.5 wants 00039 // this when using complex numbers (spetersen 04/20/06) 00040 const unsigned int seven = 7; 00041 00042 return AutoPtr<QBase> 00043 (new QClough(dim, 00044 static_cast<Order> 00045 (std::max(static_cast<unsigned int> 00046 (this->default_quadrature_order()), seven + extraorder)))); 00047 } 00048 00049 return AutoPtr<QBase> 00050 (new QGauss(dim, static_cast<Order>(this->default_quadrature_order() 00051 + extraorder))); 00052 }
| bool libMesh::FEType::operator!= | ( | const FEType & | f2 | ) | const [inline] |
| bool libMesh::FEType::operator< | ( | const FEType & | f2 | ) | const [inline] |
An ordering to make FEType useful as a std::map key
Definition at line 155 of file fe_type.h.
References family, inf_map, order, radial_family, and radial_order.
00156 { 00157 if (order != f2.order) 00158 return (order < f2.order); 00159 if (family != f2.family) 00160 return (family < f2.family); 00161 00162 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 00163 if (radial_order != f2.radial_order) 00164 return (radial_order < f2.radial_order); 00165 if (radial_family != f2.radial_family) 00166 return (radial_family < f2.radial_family); 00167 if (inf_map != f2.inf_map) 00168 return (inf_map < f2.inf_map); 00169 #endif // ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 00170 return false; 00171 }
| bool libMesh::FEType::operator== | ( | const FEType & | f2 | ) | const [inline] |
Tests equality
Definition at line 132 of file fe_type.h.
References family, inf_map, order, radial_family, and radial_order.
00133 { 00134 return (order == f2.order 00135 && family == f2.family 00136 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 00137 && radial_order == f2.radial_order 00138 && radial_family == f2.radial_family 00139 && inf_map == f2.inf_map 00140 #endif // ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 00141 ); 00142 }
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 71 of file fe_type.h.
Referenced by libMesh::FETransformationBase< OutputShape >::build(), libMesh::FEMap::build(), libMesh::FEGenericBase< OutputType >::build(), libMesh::FEAbstract::build(), libMesh::WrappedFunction< Output >::component(), libMesh::FEInterface::compute_constraints(), libMesh::GMVIO::copy_nodal_solution(), default_quadrature_rule(), libMesh::DofMap::distribute_local_dofs_node_major(), libMesh::DofMap::distribute_local_dofs_var_major(), libMesh::DofMap::dof_indices(), libMesh::FEInterface::extra_hanging_dofs(), libMesh::FEInterface::field_type(), libMesh::FEAbstract::get_family(), libMesh::System::get_info(), libMesh::FEInterface::max_order(), libMesh::Variable::n_components(), libMesh::FEInterface::n_vec_dim(), libMesh::DofMap::old_dof_indices(), libMesh::WrappedFunction< Output >::operator()(), libMesh::BoundaryProjectSolution::operator()(), libMesh::ProjectFEMSolution::operator()(), libMesh::ProjectSolution::operator()(), operator<(), operator==(), libMesh::System::project_vector(), libMesh::System::read_header(), libMesh::System::read_parallel_data(), libMesh::System::read_serialized_vectors(), libMesh::DofMap::reinit(), libMesh::System::write_header(), libMesh::System::write_parallel_data(), libMesh::System::write_serialized_vector(), and libMesh::System::write_serialized_vectors().
| 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 125 of file fe_type.h.
Referenced by libMesh::FEGenericBase< OutputType >::build_InfFE(), libMesh::System::get_info(), libMesh::FEInterface::ifem_inverse_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().
| Order libMesh::FEType::order |
The approximation order of the element.
The approximation order in radial direction of the infinite element.
Definition at line 65 of file fe_type.h.
Referenced by libMesh::FEGenericBase< OutputType >::coarsened_dof_values(), libMesh::FEInterface::compute_data(), libMesh::FEXYZ< Dim >::compute_face_values(), libMesh::FEXYZ< Dim >::compute_shape_functions(), libMesh::DofMap::constrain_p_dofs(), libMesh::GMVIO::copy_nodal_solution(), default_quadrature_order(), 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::FE< Dim, T >::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::FE< Dim, T >::n_shape_functions(), libMesh::FEInterface::nodal_soln(), libMesh::DofMap::old_dof_indices(), libMesh::WeightedPatchRecoveryErrorEstimator::EstimateError::operator()(), libMesh::PatchRecoveryErrorEstimator::EstimateError::operator()(), operator<(), operator==(), libMesh::System::read_header(), libMesh::System::read_SCALAR_dofs(), libMesh::DofMap::reinit(), libMesh::HPCoarsenTest::select_refinement(), libMesh::FEInterface::shape(), and libMesh::System::write_header().
| 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 117 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 approximation order in the base of the infinite element.
Definition at line 104 of file fe_type.h.
Referenced by libMesh::InfFE< Dim, T_radial, T_map >::attach_quadrature_rule(), libMesh::InfFE< Dim, T_radial, T_map >::combine_base_radial(), libMesh::InfFE< Dim, T_radial, T_map >::compute_data(), libMesh::InfFE< Dim, T_radial, T_map >::compute_shape_indices(), libMesh::System::get_info(), libMesh::InfFE< Dim, T_radial, T_map >::n_dofs(), libMesh::InfFE< Dim, T_radial, T_map >::n_dofs_at_node(), libMesh::InfFE< Dim, T_radial, T_map >::n_dofs_per_elem(), operator<(), operator==(), libMesh::System::read_header(), libMesh::InfFE< Dim, T_radial, T_map >::reinit(), libMesh::InfFE< Dim, T_radial, T_map >::shape(), and libMesh::System::write_header().
The documentation for this class was generated from the following files:
Site Created By: libMesh Developers
Last modified: February 05 2013 19:55:25 UTC
Hosted By: