#include <fe_abstract.h>
Public Member Functions | |
virtual | ~FEAbstract () |
virtual void | reinit (const Elem *elem, const std::vector< Point > *const pts=libmesh_nullptr, const std::vector< Real > *const weights=libmesh_nullptr)=0 |
virtual void | reinit (const Elem *elem, const unsigned int side, const Real tolerance=TOLERANCE, const std::vector< Point > *const pts=libmesh_nullptr, const std::vector< Real > *const weights=libmesh_nullptr)=0 |
virtual void | edge_reinit (const Elem *elem, const unsigned int edge, const Real tolerance=TOLERANCE, const std::vector< Point > *pts=libmesh_nullptr, const std::vector< Real > *weights=libmesh_nullptr)=0 |
virtual void | side_map (const Elem *elem, const Elem *side, const unsigned int s, const std::vector< Point > &reference_side_points, std::vector< Point > &reference_points)=0 |
unsigned int | get_dim () const |
const std::vector< Point > & | get_xyz () const |
const std::vector< Real > & | get_JxW () const |
const std::vector< RealGradient > & | get_dxyzdxi () const |
const std::vector< RealGradient > & | get_dxyzdeta () const |
const std::vector< RealGradient > & | get_dxyzdzeta () const |
const std::vector< RealGradient > & | get_d2xyzdxi2 () const |
const std::vector< RealGradient > & | get_d2xyzdeta2 () const |
const std::vector< RealGradient > & | get_d2xyzdzeta2 () const |
const std::vector< RealGradient > & | get_d2xyzdxideta () const |
const std::vector< RealGradient > & | get_d2xyzdxidzeta () const |
const std::vector< RealGradient > & | get_d2xyzdetadzeta () const |
const std::vector< Real > & | get_dxidx () const |
const std::vector< Real > & | get_dxidy () const |
const std::vector< Real > & | get_dxidz () const |
const std::vector< Real > & | get_detadx () const |
const std::vector< Real > & | get_detady () const |
const std::vector< Real > & | get_detadz () const |
const std::vector< Real > & | get_dzetadx () const |
const std::vector< Real > & | get_dzetady () const |
const std::vector< Real > & | get_dzetadz () const |
const std::vector< std::vector< Point > > & | get_tangents () const |
const std::vector< Point > & | get_normals () const |
const std::vector< Real > & | get_curvatures () const |
virtual void | attach_quadrature_rule (QBase *q)=0 |
virtual unsigned int | n_shape_functions () const =0 |
virtual unsigned int | n_quadrature_points () const =0 |
ElemType | get_type () const |
unsigned int | get_p_level () const |
FEType | get_fe_type () const |
Order | get_order () const |
void | set_fe_order (int new_order) |
virtual FEContinuity | get_continuity () const =0 |
virtual bool | is_hierarchic () const =0 |
FEFamily | get_family () const |
const FEMap & | get_fe_map () const |
void | print_JxW (std::ostream &os) const |
virtual void | print_phi (std::ostream &os) const =0 |
virtual void | print_dphi (std::ostream &os) const =0 |
virtual void | print_d2phi (std::ostream &os) const =0 |
void | print_xyz (std::ostream &os) const |
void | print_info (std::ostream &os) const |
Static Public Member Functions | |
static std::unique_ptr< FEAbstract > | build (const unsigned int dim, const FEType &type) |
static bool | on_reference_element (const Point &p, const ElemType t, const Real eps=TOLERANCE) |
static void | get_refspace_nodes (const ElemType t, std::vector< Point > &nodes) |
static void | compute_node_constraints (NodeConstraints &constraints, const Elem *elem) |
static void | compute_periodic_node_constraints (NodeConstraints &constraints, const PeriodicBoundaries &boundaries, const MeshBase &mesh, const PointLocatorBase *point_locator, const Elem *elem) |
static std::string | get_info () |
static void | print_info (std::ostream &out=libMesh::out) |
static unsigned int | n_objects () |
static void | enable_print_counter_info () |
static void | disable_print_counter_info () |
Protected Types | |
typedef std::map< std::string, std::pair< unsigned int, unsigned int > > | Counts |
Protected Member Functions | |
FEAbstract (const unsigned int dim, const FEType &fet) | |
virtual void | compute_shape_functions (const Elem *, const std::vector< Point > &)=0 |
virtual bool | shapes_need_reinit () const =0 |
void | increment_constructor_count (const std::string &name) |
void | increment_destructor_count (const std::string &name) |
Protected Attributes | |
std::unique_ptr< FEMap > | _fe_map |
const unsigned int | dim |
bool | calculations_started |
bool | calculate_phi |
bool | calculate_dphi |
bool | calculate_d2phi |
bool | calculate_curl_phi |
bool | calculate_div_phi |
bool | calculate_dphiref |
FEType | fe_type |
ElemType | elem_type |
unsigned int | _p_level |
QBase * | qrule |
bool | shapes_on_quadrature |
Static Protected Attributes | |
static Counts | _counts |
static Threads::atomic< unsigned int > | _n_objects |
static Threads::spin_mutex | _mutex |
static bool | _enable_print_counter = true |
Friends | |
std::ostream & | operator<< (std::ostream &os, const FEAbstract &fe) |
This class forms the foundation from which generic finite elements may be derived. In the current implementation, the templated derived class FE
offers a wide variety of commonly used finite element concepts. Check there for details. Use the FEAbstract::build()
method to create an object of any of the derived classes.
All calls to static members of the FE
classes should be requested through the FEInterface
. This interface class approximates runtime polymorphism for the templated finite element classes. Even internal library classes, like DofMap
, request the number of DOFs through this interface class. This approach also enables the co-existence of various element-based schemes.
Definition at line 93 of file fe_abstract.h.
|
protectedinherited |
Data structure to log the information. The log is identified by the class name.
Definition at line 119 of file reference_counter.h.
|
inlineprotected |
Constructor. Optionally initializes required data structures. Protected so that this base class cannot be explicitly instantiated.
Definition at line 608 of file fe_abstract.h.
|
inlinevirtual |
|
pure virtual |
Provides the class with the quadrature rule. Implement this in derived classes.
Implemented in libMesh::FESubdivision, libMesh::InfFE< Dim, T_radial, T_map >, libMesh::FE< Dim, T >, libMesh::FE< 2, SUBDIVISION >, libMesh::FE< Dim, HIERARCHIC >, libMesh::FE< Dim, SCALAR >, libMesh::FE< Dim, L2_LAGRANGE >, libMesh::FE< Dim, NEDELEC_ONE >, libMesh::FE< Dim, HERMITE >, libMesh::FE< Dim, CLOUGH >, libMesh::FE< Dim, MONOMIAL >, libMesh::FE< Dim, XYZ >, libMesh::FE< Dim, LAGRANGE >, libMesh::FE< Dim, L2_HIERARCHIC >, and libMesh::FE< Dim, LAGRANGE_VEC >.
Referenced by get_curvatures().
|
static |
Builds a specific finite element type.
Definition at line 44 of file fe_abstract.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::SCALAR, libMesh::SUBDIVISION, libMesh::SZABAB, and libMesh::XYZ.
Referenced by libMesh::DGFEMContext::DGFEMContext(), libMesh::FEMContext::init_internal_data(), and libMesh::DofMap::use_coupled_neighbor_dofs().
|
static |
Computes the nodal constraint contributions (for non-conforming adapted meshes), using Lagrange geometry
Definition at line 793 of file fe_abstract.C.
References std::abs(), libMesh::Elem::build_side_ptr(), libMesh::Elem::default_order(), libMesh::Elem::dim(), fe_type, libMesh::FEInterface::inverse_map(), libMesh::LAGRANGE, libMesh::Elem::level(), libmesh_nullptr, libMesh::FEInterface::n_dofs(), libMesh::Elem::neighbor_ptr(), libMesh::Elem::parent(), libMesh::Real, libMesh::remote_elem, libMesh::FEInterface::shape(), libMesh::Elem::side_index_range(), libMesh::Threads::spin_mtx, and libMesh::Elem::subactive().
|
static |
Computes the node position constraint equation contributions (for meshes with periodic boundary conditions)
Definition at line 936 of file fe_abstract.C.
References libMesh::Elem::active(), libMesh::PeriodicBoundaries::boundary(), libMesh::BoundaryInfo::boundary_ids(), libMesh::Elem::build_side_ptr(), libMesh::Elem::default_order(), libMesh::Elem::dim(), fe_type, libMesh::MeshBase::get_boundary_info(), libMesh::PeriodicBoundaryBase::get_corresponding_pos(), libMesh::invalid_uint, libMesh::FEInterface::inverse_map(), libMesh::LAGRANGE, libMesh::Elem::level(), libMesh::FEInterface::n_dofs(), libMesh::PeriodicBoundaries::neighbor(), libMesh::Elem::neighbor_ptr(), libMesh::PeriodicBoundaryBase::pairedboundary, libMesh::Real, libMesh::FEInterface::shape(), libMesh::Elem::side_index_range(), libMesh::BoundaryInfo::side_with_boundary_id(), and libMesh::Threads::spin_mtx.
|
protectedpure virtual |
After having updated the jacobian and the transformation from local to global coordinates in FEMap::compute_map()
, the first derivatives of the shape functions are transformed to global coordinates, giving dphi
, dphidx
, dphidy
, and dphidz
. This method should rarely be re-defined in derived classes, but still should be usable for children. Therefore, keep it protected. This needs to be implemented in the derived class since this function depends on whether the shape functions are vector-valued or not.
Implemented in libMesh::FEXYZ< Dim >, libMesh::InfFE< Dim, T_radial, T_map >, libMesh::FEGenericBase< OutputType >, and libMesh::FEGenericBase< FEOutputType< T >::type >.
Referenced by get_fe_map().
|
staticinherited |
Definition at line 107 of file reference_counter.C.
References libMesh::ReferenceCounter::_enable_print_counter.
Referenced by libMesh::LibMeshInit::LibMeshInit(), and libMesh::ReferenceCounter::n_objects().
|
pure virtual |
Reinitializes all the physical element-dependent data based on the edge
of the element elem
. The tolerance
parameter is passed to the involved call to inverse_map()
. By default the element data are computed at the quadrature points specified by the quadrature rule qrule
, but any set of points on the reference edge element may be specified in the optional argument pts
.
Implemented in libMesh::InfFE< Dim, T_radial, T_map >, libMesh::FE< Dim, T >, libMesh::FE< 2, SUBDIVISION >, libMesh::FE< Dim, HIERARCHIC >, libMesh::FE< Dim, SCALAR >, libMesh::FE< Dim, L2_LAGRANGE >, libMesh::FE< Dim, NEDELEC_ONE >, libMesh::FE< Dim, HERMITE >, libMesh::FE< Dim, CLOUGH >, libMesh::FE< Dim, MONOMIAL >, libMesh::FE< Dim, XYZ >, libMesh::FE< Dim, LAGRANGE >, libMesh::FE< Dim, L2_HIERARCHIC >, and libMesh::FE< Dim, LAGRANGE_VEC >.
|
staticinherited |
Methods to enable/disable the reference counter output from print_info()
Definition at line 101 of file reference_counter.C.
References libMesh::ReferenceCounter::_enable_print_counter.
Referenced by libMesh::ReferenceCounter::n_objects().
|
pure virtual |
Implemented in libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::InfFE< Dim, T_radial, T_map >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< 2, SUBDIVISION >, libMesh::FE< Dim, HIERARCHIC >, libMesh::FE< Dim, SCALAR >, libMesh::FE< Dim, L2_LAGRANGE >, libMesh::FE< Dim, NEDELEC_ONE >, libMesh::FE< Dim, HERMITE >, libMesh::FE< Dim, CLOUGH >, libMesh::FE< Dim, MONOMIAL >, libMesh::FE< Dim, XYZ >, libMesh::FE< Dim, LAGRANGE >, libMesh::FE< Dim, L2_HIERARCHIC >, libMesh::FE< Dim, LAGRANGE_VEC >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, and libMesh::FE< Dim, T >.
Referenced by libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::operator()(), and set_fe_order().
|
inline |
Definition at line 383 of file fe_abstract.h.
References _fe_map, attach_quadrature_rule(), n_quadrature_points(), and n_shape_functions().
|
inline |
Definition at line 270 of file fe_abstract.h.
References _fe_map.
|
inline |
Definition at line 300 of file fe_abstract.h.
References _fe_map.
|
inline |
Definition at line 264 of file fe_abstract.h.
References _fe_map.
|
inline |
Definition at line 286 of file fe_abstract.h.
References _fe_map.
|
inline |
Definition at line 294 of file fe_abstract.h.
References _fe_map.
|
inline |
Definition at line 278 of file fe_abstract.h.
References _fe_map.
|
inline |
Definition at line 330 of file fe_abstract.h.
References _fe_map.
|
inline |
Definition at line 337 of file fe_abstract.h.
References _fe_map.
|
inline |
Definition at line 344 of file fe_abstract.h.
References _fe_map.
|
inline |
|
inline |
Definition at line 309 of file fe_abstract.h.
References _fe_map.
|
inline |
Definition at line 316 of file fe_abstract.h.
References _fe_map.
|
inline |
Definition at line 323 of file fe_abstract.h.
References _fe_map.
|
inline |
Definition at line 251 of file fe_abstract.h.
References _fe_map.
|
inline |
Definition at line 244 of file fe_abstract.h.
References _fe_map.
|
inline |
Definition at line 258 of file fe_abstract.h.
References _fe_map.
|
inline |
Definition at line 351 of file fe_abstract.h.
References _fe_map.
|
inline |
Definition at line 358 of file fe_abstract.h.
References _fe_map.
|
inline |
Definition at line 365 of file fe_abstract.h.
References _fe_map.
|
inline |
Definition at line 447 of file fe_abstract.h.
References libMesh::FEType::family, and fe_type.
Referenced by libMesh::FE< Dim, T >::FE().
|
inline |
Definition at line 452 of file fe_abstract.h.
References _fe_map, compute_shape_functions(), operator<<, print_d2phi(), print_dphi(), print_info(), print_JxW(), print_phi(), and print_xyz().
Referenced by libMesh::HCurlFETransformation< OutputShape >::init_map_d2phi(), libMesh::H1FETransformation< OutputShape >::init_map_d2phi(), libMesh::HCurlFETransformation< OutputShape >::init_map_dphi(), libMesh::H1FETransformation< OutputShape >::init_map_dphi(), libMesh::HCurlFETransformation< OutputShape >::init_map_phi(), libMesh::HCurlFETransformation< OutputShape >::map_curl(), libMesh::H1FETransformation< OutputShape >::map_curl(), libMesh::H1FETransformation< OutputShape >::map_d2phi(), libMesh::H1FETransformation< OutputShape >::map_div(), libMesh::H1FETransformation< OutputShape >::map_dphi(), and libMesh::HCurlFETransformation< OutputShape >::map_phi().
|
inline |
Definition at line 421 of file fe_abstract.h.
References fe_type.
Referenced by libMesh::FEMContext::build_new_fe(), libMesh::HCurlFETransformation< OutputShape >::map_phi(), libMesh::H1FETransformation< OutputShape >::map_phi(), libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::operator()(), and libMesh::JumpErrorEstimator::reinit_sides().
|
staticinherited |
Gets a string containing the reference information.
Definition at line 47 of file reference_counter.C.
References libMesh::ReferenceCounter::_counts, and libMesh::Quality::name().
Referenced by libMesh::ReferenceCounter::print_info().
|
inline |
Definition at line 237 of file fe_abstract.h.
References _fe_map.
Referenced by libMesh::ExactSolution::_compute_error(), libMesh::DiscontinuityMeasure::boundary_side_integration(), libMesh::KellyErrorEstimator::boundary_side_integration(), libMesh::System::calculate_norm(), libMesh::ExactErrorEstimator::find_squared_element_error(), libMesh::FEMSystem::init_context(), libMesh::LaplacianErrorEstimator::internal_side_integration(), libMesh::DiscontinuityMeasure::internal_side_integration(), libMesh::KellyErrorEstimator::internal_side_integration(), and libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::operator()().
|
inline |
Definition at line 377 of file fe_abstract.h.
References _fe_map.
Referenced by libMesh::KellyErrorEstimator::boundary_side_integration(), libMesh::ParsedFEMFunction< Output >::eval_args(), libMesh::ParsedFEMFunction< Output >::init_context(), libMesh::KellyErrorEstimator::init_context(), and libMesh::KellyErrorEstimator::internal_side_integration().
|
inline |
Definition at line 426 of file fe_abstract.h.
References _p_level, fe_type, and libMesh::FEType::order.
|
inline |
Definition at line 416 of file fe_abstract.h.
References _p_level.
|
static |
nodes
based on the element type. Definition at line 256 of file fe_abstract.C.
References libMesh::EDGE2, libMesh::EDGE3, libMesh::HEX20, libMesh::HEX27, libMesh::HEX8, libMesh::PRISM15, libMesh::PRISM18, libMesh::PRISM6, libMesh::PYRAMID13, libMesh::PYRAMID14, libMesh::PYRAMID5, libMesh::QUAD4, libMesh::QUAD8, libMesh::QUAD9, libMesh::QUADSHELL4, libMesh::QUADSHELL8, libMesh::TET10, libMesh::TET4, libMesh::TRI3, libMesh::TRI6, and libMesh::TRISHELL3.
|
inline |
Definition at line 371 of file fe_abstract.h.
References _fe_map.
|
inline |
Definition at line 410 of file fe_abstract.h.
References elem_type.
Referenced by libMesh::InfFE< Dim, T_radial, T_map >::reinit().
|
inline |
xyz
spatial locations of the quadrature points on the element. Definition at line 230 of file fe_abstract.h.
References _fe_map.
Referenced by libMesh::ExactSolution::_compute_error(), libMesh::DiscontinuityMeasure::boundary_side_integration(), libMesh::KellyErrorEstimator::boundary_side_integration(), libMesh::JumpErrorEstimator::estimate_error(), libMesh::ParsedFEMFunction< Output >::eval_args(), libMesh::ExactErrorEstimator::find_squared_element_error(), libMesh::ParsedFEMFunction< Output >::init_context(), libMesh::DGFEMContext::neighbor_side_fe_reinit(), libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::operator()(), and libMesh::InfFE< Dim, T_radial, T_map >::reinit().
|
inlineprotectedinherited |
Increments the construction counter. Should be called in the constructor of any derived class that will be reference counted.
Definition at line 185 of file reference_counter.h.
References libMesh::ReferenceCounter::_counts, libMesh::Quality::name(), and libMesh::Threads::spin_mtx.
Referenced by libMesh::ReferenceCounter::n_objects(), and libMesh::ReferenceCountedObject< RBParametrized >::ReferenceCountedObject().
|
inlineprotectedinherited |
Increments the destruction counter. Should be called in the destructor of any derived class that will be reference counted.
Definition at line 198 of file reference_counter.h.
References libMesh::ReferenceCounter::_counts, libMesh::Quality::name(), and libMesh::Threads::spin_mtx.
Referenced by libMesh::ReferenceCounter::n_objects(), and libMesh::ReferenceCountedObject< RBParametrized >::~ReferenceCountedObject().
|
pure virtual |
true
if the finite element's higher order shape functions are hierarchic Implemented in libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::InfFE< Dim, T_radial, T_map >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< 2, SUBDIVISION >, libMesh::FE< Dim, HIERARCHIC >, libMesh::FE< Dim, SCALAR >, libMesh::FE< Dim, L2_LAGRANGE >, libMesh::FE< Dim, NEDELEC_ONE >, libMesh::FE< Dim, HERMITE >, libMesh::FE< Dim, CLOUGH >, libMesh::FE< Dim, MONOMIAL >, libMesh::FE< Dim, XYZ >, libMesh::FE< Dim, LAGRANGE >, libMesh::FE< Dim, L2_HIERARCHIC >, libMesh::FE< Dim, LAGRANGE_VEC >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, and libMesh::FE< Dim, T >.
Referenced by set_fe_order().
|
inlinestaticinherited |
Prints the number of outstanding (created, but not yet destroyed) objects.
Definition at line 85 of file reference_counter.h.
References libMesh::ReferenceCounter::_n_objects, libMesh::ReferenceCounter::disable_print_counter_info(), libMesh::ReferenceCounter::enable_print_counter_info(), libMesh::ReferenceCounter::increment_constructor_count(), libMesh::ReferenceCounter::increment_destructor_count(), and libMesh::Quality::name().
Referenced by libMesh::System::hide_output(), and libMesh::LibMeshInit::LibMeshInit().
|
pure virtual |
Implemented in libMesh::InfFE< Dim, T_radial, T_map >, libMesh::FE< Dim, T >, libMesh::FE< 2, SUBDIVISION >, libMesh::FE< Dim, HIERARCHIC >, libMesh::FE< Dim, SCALAR >, libMesh::FE< Dim, L2_LAGRANGE >, libMesh::FE< Dim, NEDELEC_ONE >, libMesh::FE< Dim, HERMITE >, libMesh::FE< Dim, CLOUGH >, libMesh::FE< Dim, MONOMIAL >, libMesh::FE< Dim, XYZ >, libMesh::FE< Dim, LAGRANGE >, libMesh::FE< Dim, L2_HIERARCHIC >, and libMesh::FE< Dim, LAGRANGE_VEC >.
Referenced by libMesh::DiscontinuityMeasure::boundary_side_integration(), libMesh::KellyErrorEstimator::boundary_side_integration(), get_curvatures(), libMesh::LaplacianErrorEstimator::internal_side_integration(), libMesh::DiscontinuityMeasure::internal_side_integration(), and libMesh::KellyErrorEstimator::internal_side_integration().
|
pure virtual |
Implemented in libMesh::InfFE< Dim, T_radial, T_map >, libMesh::FE< Dim, T >, libMesh::FE< 2, SUBDIVISION >, libMesh::FE< Dim, HIERARCHIC >, libMesh::FE< Dim, SCALAR >, libMesh::FE< Dim, L2_LAGRANGE >, libMesh::FE< Dim, NEDELEC_ONE >, libMesh::FE< Dim, HERMITE >, libMesh::FE< Dim, CLOUGH >, libMesh::FE< Dim, MONOMIAL >, libMesh::FE< Dim, XYZ >, libMesh::FE< Dim, LAGRANGE >, libMesh::FE< Dim, L2_HIERARCHIC >, and libMesh::FE< Dim, LAGRANGE_VEC >.
Referenced by get_curvatures().
|
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 here the parameter eps
can be specified to indicate a tolerance. For example, Definition at line 554 of file fe_abstract.C.
References libMesh::EDGE2, libMesh::EDGE3, libMesh::EDGE4, libMesh::HEX20, libMesh::HEX27, libMesh::HEX8, libMesh::INFHEX16, libMesh::INFHEX18, libMesh::INFHEX8, libMesh::INFPRISM12, libMesh::INFPRISM6, libMesh::NODEELEM, libMesh::PRISM15, libMesh::PRISM18, libMesh::PRISM6, libMesh::PYRAMID13, libMesh::PYRAMID14, libMesh::PYRAMID5, libMesh::QUAD4, libMesh::QUAD8, libMesh::QUAD9, libMesh::QUADSHELL4, libMesh::QUADSHELL8, libMesh::Real, libMesh::TET10, libMesh::TET4, libMesh::TRI3, libMesh::TRI6, and libMesh::TRISHELL3.
Referenced by libMesh::FEInterface::ifem_on_reference_element(), libMesh::FE< Dim, T >::inverse_map(), and libMesh::FEInterface::on_reference_element().
|
pure virtual |
Prints the value of each shape function's second derivatives at each quadrature point. Implement in derived class since this depends on whether the element is vector-valued or not.
Implemented in libMesh::FEGenericBase< OutputType >, and libMesh::FEGenericBase< FEOutputType< T >::type >.
Referenced by get_fe_map().
|
pure virtual |
Prints the value of each shape function's derivative at each quadrature point. Implement in derived class since this depends on whether the element is vector-valued or not.
Implemented in libMesh::FEGenericBase< OutputType >, and libMesh::FEGenericBase< FEOutputType< T >::type >.
Referenced by get_fe_map(), and print_info().
|
staticinherited |
Prints the reference information, by default to libMesh::out
.
Definition at line 88 of file reference_counter.C.
References libMesh::ReferenceCounter::_enable_print_counter, and libMesh::ReferenceCounter::get_info().
Referenced by libMesh::LibMeshInit::LibMeshInit().
void libMesh::FEAbstract::print_info | ( | std::ostream & | os | ) | const |
Prints all the relevant information about the current element.
Definition at line 766 of file fe_abstract.C.
References print_dphi(), print_JxW(), print_phi(), and print_xyz().
Referenced by get_fe_map(), and libMesh::operator<<().
void libMesh::FEAbstract::print_JxW | ( | std::ostream & | os | ) | const |
Prints the Jacobian times the weight for each quadrature point.
Definition at line 753 of file fe_abstract.C.
References _fe_map.
Referenced by get_fe_map(), and print_info().
|
pure virtual |
Prints the value of each shape function at each quadrature point. Implement in derived class since this depends on whether the element is vector-valued or not.
Implemented in libMesh::FEGenericBase< OutputType >, and libMesh::FEGenericBase< FEOutputType< T >::type >.
Referenced by get_fe_map(), and print_info().
void libMesh::FEAbstract::print_xyz | ( | std::ostream & | os | ) | const |
Prints the spatial location of each quadrature point (on the physical element).
Definition at line 760 of file fe_abstract.C.
References _fe_map.
Referenced by get_fe_map(), and print_info().
|
pure virtual |
This is at the core of this class. Use this for each new element in the mesh. Reinitializes the requested physical element-dependent data based on the current element elem
. By default the element data are computed at the quadrature points specified by the quadrature rule qrule
, but any set of points on the reference element may be specified in the optional argument pts
.
get_phi()
or get_d2phi()
have been called, so all such accessors should be called before the first reinit()
. Implemented in libMesh::FEXYZ< Dim >, libMesh::FESubdivision, libMesh::InfFE< Dim, T_radial, T_map >, libMesh::FE< Dim, T >, libMesh::FE< 2, SUBDIVISION >, libMesh::FE< Dim, HIERARCHIC >, libMesh::FE< Dim, SCALAR >, libMesh::FE< Dim, L2_LAGRANGE >, libMesh::FE< Dim, NEDELEC_ONE >, libMesh::FE< Dim, HERMITE >, libMesh::FE< Dim, CLOUGH >, libMesh::FE< Dim, MONOMIAL >, libMesh::FE< Dim, XYZ >, libMesh::FE< Dim, LAGRANGE >, libMesh::FE< Dim, L2_HIERARCHIC >, and libMesh::FE< Dim, LAGRANGE_VEC >.
Referenced by libMesh::ExactSolution::_compute_error(), libMesh::FEMContext::build_new_fe(), libMesh::System::calculate_norm(), libMesh::ExactErrorEstimator::find_squared_element_error(), and libMesh::JumpErrorEstimator::reinit_sides().
|
pure virtual |
Reinitializes all the physical element-dependent data based on the side
of the element elem
. The tolerance
parameter is passed to the involved call to inverse_map()
. By default the element data are computed at the quadrature points specified by the quadrature rule qrule
, but any set of points on the reference side element may be specified in the optional argument pts
.
Implemented in libMesh::FEXYZ< Dim >, libMesh::FESubdivision, libMesh::InfFE< Dim, T_radial, T_map >, libMesh::FE< Dim, T >, libMesh::FE< 2, SUBDIVISION >, libMesh::FE< Dim, HIERARCHIC >, libMesh::FE< Dim, SCALAR >, libMesh::FE< Dim, L2_LAGRANGE >, libMesh::FE< Dim, NEDELEC_ONE >, libMesh::FE< Dim, HERMITE >, libMesh::FE< Dim, CLOUGH >, libMesh::FE< Dim, MONOMIAL >, libMesh::FE< Dim, XYZ >, libMesh::FE< Dim, LAGRANGE >, libMesh::FE< Dim, L2_HIERARCHIC >, libMesh::FE< Dim, LAGRANGE_VEC >, libMesh::FEXYZ< Dim >, and libMesh::FEXYZ< Dim >.
|
inline |
Sets the base FE order of the finite element.
Definition at line 431 of file fe_abstract.h.
References fe_type, get_continuity(), is_hierarchic(), and libMesh::FEType::order.
|
protectedpure virtual |
true
when the shape functions (for this FEFamily
) depend on the particular element, and therefore needs to be re-initialized for each new element. false
otherwise. Implemented in libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::InfFE< Dim, T_radial, T_map >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< 2, SUBDIVISION >, libMesh::FE< Dim, HIERARCHIC >, libMesh::FE< Dim, SCALAR >, libMesh::FE< Dim, L2_LAGRANGE >, libMesh::FE< Dim, NEDELEC_ONE >, libMesh::FE< Dim, HERMITE >, libMesh::FE< Dim, CLOUGH >, libMesh::FE< Dim, MONOMIAL >, libMesh::FE< Dim, XYZ >, libMesh::FE< Dim, LAGRANGE >, libMesh::FE< Dim, L2_HIERARCHIC >, libMesh::FE< Dim, LAGRANGE_VEC >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, and libMesh::FE< Dim, T >.
|
pure virtual |
Computes the reference space quadrature points on the side of an element based on the side quadrature points.
Implemented in libMesh::FE< Dim, T >, libMesh::InfFE< Dim, T_radial, T_map >, libMesh::FE< Dim, T >, libMesh::FE< 2, SUBDIVISION >, libMesh::FE< Dim, HIERARCHIC >, libMesh::FE< Dim, SCALAR >, libMesh::FE< Dim, L2_LAGRANGE >, libMesh::FE< Dim, NEDELEC_ONE >, libMesh::FE< Dim, HERMITE >, libMesh::FE< Dim, CLOUGH >, libMesh::FE< Dim, MONOMIAL >, libMesh::FE< Dim, XYZ >, libMesh::FE< Dim, LAGRANGE >, libMesh::FE< Dim, L2_HIERARCHIC >, and libMesh::FE< Dim, LAGRANGE_VEC >.
|
friend |
Same as above, but allows you to print to a stream.
Definition at line 782 of file fe_abstract.C.
Referenced by get_fe_map().
|
staticprotectedinherited |
Actually holds the data.
Definition at line 124 of file reference_counter.h.
Referenced by libMesh::ReferenceCounter::get_info(), libMesh::ReferenceCounter::increment_constructor_count(), and libMesh::ReferenceCounter::increment_destructor_count().
|
staticprotectedinherited |
Flag to control whether reference count information is printed when print_info is called.
Definition at line 143 of file reference_counter.h.
Referenced by libMesh::ReferenceCounter::disable_print_counter_info(), libMesh::ReferenceCounter::enable_print_counter_info(), and libMesh::ReferenceCounter::print_info().
|
protected |
Definition at line 517 of file fe_abstract.h.
Referenced by libMesh::FESubdivision::attach_quadrature_rule(), libMesh::InfFE< Dim, T_radial, T_map >::combine_base_radial(), libMesh::InfFE< Dim, T_radial, T_map >::compute_shape_functions(), get_curvatures(), get_d2xyzdeta2(), get_d2xyzdetadzeta(), get_d2xyzdxi2(), get_d2xyzdxideta(), get_d2xyzdxidzeta(), get_d2xyzdzeta2(), get_detadx(), get_detady(), get_detadz(), get_dxidx(), get_dxidy(), get_dxidz(), get_dxyzdeta(), get_dxyzdxi(), get_dxyzdzeta(), get_dzetadx(), get_dzetady(), get_dzetadz(), get_fe_map(), get_JxW(), get_normals(), get_tangents(), get_xyz(), libMesh::InfFE< Dim, T_radial, T_map >::init_shape_functions(), libMesh::FESubdivision::init_shape_functions(), print_JxW(), print_xyz(), and libMesh::InfFE< Dim, T_radial, T_map >::reinit().
|
staticprotectedinherited |
Mutual exclusion object to enable thread-safe reference counting.
Definition at line 137 of file reference_counter.h.
|
staticprotectedinherited |
The number of objects. Print the reference count information when the number returns to 0.
Definition at line 132 of file reference_counter.h.
Referenced by libMesh::ReferenceCounter::n_objects(), libMesh::ReferenceCounter::ReferenceCounter(), and libMesh::ReferenceCounter::~ReferenceCounter().
|
protected |
The p refinement level the current data structures are set up for.
Definition at line 579 of file fe_abstract.h.
Referenced by get_order(), and get_p_level().
|
mutableprotected |
Should we calculate shape function curls?
Definition at line 549 of file fe_abstract.h.
Referenced by libMesh::FEGenericBase< FEOutputType< T >::type >::get_curl_phi().
|
mutableprotected |
Should we calculate shape function hessians?
Definition at line 544 of file fe_abstract.h.
Referenced by libMesh::FESubdivision::attach_quadrature_rule(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phi(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phideta2(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidetadzeta(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidx2(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidxdy(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidxdz(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidxi2(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidxideta(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidxidzeta(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidy2(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidydz(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidz2(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidzeta2(), and libMesh::InfFE< Dim, T_radial, T_map >::reinit().
|
mutableprotected |
Should we calculate shape function divergences?
Definition at line 554 of file fe_abstract.h.
Referenced by libMesh::FEGenericBase< FEOutputType< T >::type >::get_div_phi().
|
mutableprotected |
Should we calculate shape function gradients?
Definition at line 539 of file fe_abstract.h.
Referenced by libMesh::FEGenericBase< FEOutputType< T >::type >::get_dphi(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_dphidx(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_dphidy(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_dphidz(), and libMesh::InfFE< Dim, T_radial, T_map >::reinit().
|
mutableprotected |
Should we calculate reference shape function gradients?
Definition at line 559 of file fe_abstract.h.
Referenced by libMesh::FEGenericBase< FEOutputType< T >::type >::get_curl_phi(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phi(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phideta2(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidetadzeta(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidx2(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidxdy(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidxdz(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidxi2(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidxideta(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidxidzeta(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidy2(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidydz(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidz2(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidzeta2(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_div_phi(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_dphi(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_dphideta(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_dphidx(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_dphidxi(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_dphidy(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_dphidz(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_dphidzeta(), and libMesh::InfFE< Dim, T_radial, T_map >::reinit().
|
mutableprotected |
Should we calculate shape functions?
Definition at line 534 of file fe_abstract.h.
Referenced by libMesh::FEGenericBase< FEOutputType< T >::type >::get_phi(), and libMesh::InfFE< Dim, T_radial, T_map >::reinit().
|
mutableprotected |
Have calculations with this object already been started? Then all get_* functions should already have been called.
Definition at line 529 of file fe_abstract.h.
Referenced by libMesh::FEGenericBase< FEOutputType< T >::type >::get_curl_phi(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phi(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phideta2(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidetadzeta(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidx2(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidxdy(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidxdz(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidxi2(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidxideta(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidxidzeta(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidy2(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidydz(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidz2(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidzeta2(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_div_phi(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_dphi(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_dphideta(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_dphidx(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_dphidxi(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_dphidy(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_dphidz(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_dphidzeta(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_phi(), and libMesh::FESubdivision::init_shape_functions().
|
protected |
The dimensionality of the object
Definition at line 523 of file fe_abstract.h.
Referenced by libMesh::FESubdivision::attach_quadrature_rule(), libMesh::InfFE< Dim, T_radial, T_map >::compute_shape_functions(), get_dim(), and libMesh::InfFE< Dim, T_radial, T_map >::reinit().
|
protected |
The element type the current data structures are set up for.
Definition at line 573 of file fe_abstract.h.
Referenced by libMesh::FESubdivision::attach_quadrature_rule(), get_type(), and libMesh::InfFE< Dim, T_radial, T_map >::reinit().
|
protected |
The finite element type for this object.
Definition at line 567 of file fe_abstract.h.
Referenced by libMesh::InfFE< Dim, T_radial, T_map >::attach_quadrature_rule(), libMesh::InfFE< Dim, T_radial, T_map >::combine_base_radial(), compute_node_constraints(), compute_periodic_node_constraints(), get_family(), get_fe_type(), get_order(), libMesh::InfFE< Dim, T_radial, T_map >::InfFE(), libMesh::InfFE< Dim, T_radial, T_map >::init_radial_shape_functions(), libMesh::InfFE< Dim, T_radial, T_map >::init_shape_functions(), libMesh::FESubdivision::init_shape_functions(), libMesh::InfFE< Dim, T_radial, T_map >::reinit(), and set_fe_order().
|
protected |
A pointer to the quadrature rule employed
Definition at line 584 of file fe_abstract.h.
Referenced by libMesh::InfFE< Dim, T_radial, T_map >::attach_quadrature_rule(), and libMesh::FESubdivision::attach_quadrature_rule().
|
protected |
A flag indicating if current data structures correspond to quadrature rule points
Definition at line 590 of file fe_abstract.h.
Referenced by libMesh::FESubdivision::attach_quadrature_rule().