libMesh
Static Public Member Functions | Private Member Functions | List of all members
libMesh::InfFEBase Class Reference

This nested class contains most of the static methods related to the base part of an infinite element. More...

#include <inf_fe.h>

Static Public Member Functions

static std::unique_ptr< const Elembuild_elem (const Elem *inf_elem)
 Build the base element of an infinite element. More...
 
static ElemType get_elem_type (const ElemType type)
 
static unsigned int n_base_mapping_sf (const Elem &base_elem, const Order base_mapping_order)
 

Private Member Functions

 InfFEBase ()
 Never use an object of this type. More...
 

Detailed Description

This nested class contains most of the static methods related to the base part of an infinite element.

Only static members are provided, for use within InfFE and InfFEMap.

Author
Daniel Dreyer
Date
2003

Definition at line 151 of file inf_fe.h.

Constructor & Destructor Documentation

◆ InfFEBase()

libMesh::InfFEBase::InfFEBase ( )
inlineprivate

Never use an object of this type.

Definition at line 158 of file inf_fe.h.

158 {}

Member Function Documentation

◆ build_elem()

std::unique_ptr< const Elem > libMesh::InfFEBase::build_elem ( const Elem inf_elem)
static

Build the base element of an infinite element.

Here "base" refers to the non-infinite/traditional geometric element face which is associated with this infinite element. A const base Elem is built from a const input InfElem.

Definition at line 38 of file inf_fe_base_radial.C.

References libMesh::Elem::build_side_ptr().

Referenced by InfFERadialTest::base_point(), libMesh::InfFE< Dim, T_radial, T_map >::compute_shape_indices(), libMesh::InfFEMap::inverse_map(), libMesh::InfFEMap::map(), and libMesh::InfFE< Dim, T_radial, T_map >::update_base_elem().

39 {
40  return inf_elem->build_side_ptr(0);
41 }

◆ get_elem_type()

ElemType libMesh::InfFEBase::get_elem_type ( const ElemType  type)
static
Returns
The base element associated to type. This is, for example, TRI3 for INFPRISM6.

Definition at line 45 of file inf_fe_base_radial.C.

References libMesh::EDGE2, libMesh::EDGE3, libMesh::Utility::enum_to_string(), libMesh::INFEDGE2, libMesh::INFHEX16, libMesh::INFHEX18, libMesh::INFHEX8, libMesh::INFPRISM12, libMesh::INFPRISM6, libMesh::INFQUAD4, libMesh::INFQUAD6, libMesh::INVALID_ELEM, libMesh::QUAD4, libMesh::QUAD8, libMesh::QUAD9, libMesh::TRI3, and libMesh::TRI6.

Referenced by libMesh::InfFE< Dim, T_radial, T_map >::compute_shape_functions(), libMesh::InfFE< Dim, T_radial, T_map >::compute_shape_indices(), 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(), libMesh::InfFE< Dim, T_radial, T_map >::shape(), and libMesh::InfFE< Dim, T_radial, T_map >::shape_deriv().

46 {
47  switch (type)
48  {
49  // 3D infinite elements:
50  // with Dim=3 -> infinite elements on their own
51  case INFHEX8:
52  return QUAD4;
53 
54  case INFHEX16:
55  return QUAD8;
56 
57  case INFHEX18:
58  return QUAD9;
59 
60  case INFPRISM6:
61  return TRI3;
62 
63  case INFPRISM12:
64  return TRI6;
65 
66  // 2D infinite elements:
67  // with Dim=3 -> used as boundary condition,
68  // with Dim=2 -> infinite elements on their own
69  case INFQUAD4:
70  return EDGE2;
71 
72  case INFQUAD6:
73  return EDGE3;
74 
75  // 1D infinite elements:
76  // with Dim=2 -> used as boundary condition,
77  // with Dim=1 -> infinite elements on their own,
78  // but no base element!
79  case INFEDGE2:
80  return INVALID_ELEM;
81 
82  default:
83  libmesh_error_msg("ERROR: Unsupported element type!: " << Utility::enum_to_string(type));
84  }
85 }
std::string enum_to_string(const T e)

◆ n_base_mapping_sf()

unsigned int libMesh::InfFEBase::n_base_mapping_sf ( const Elem base_elem,
const Order  base_mapping_order 
)
static
Returns
The number of shape functions used in the mapping in the base element of type base_elem_type mapped with order base_mapping_order

Definition at line 91 of file inf_fe_base_radial.C.

References libMesh::Elem::dim(), libMesh::FE< Dim, T >::n_shape_functions(), and libMesh::Elem::type().

93 {
94  switch (base_elem.dim())
95  {
96  case 0:
97  return 1;
98  case 1:
99  return FE<1,LAGRANGE>::n_shape_functions (base_elem.type(),
100  base_mapping_order);
101  case 2:
102  return FE<2,LAGRANGE>::n_shape_functions (base_elem.type(),
103  base_mapping_order);
104  default:
105  libmesh_error_msg("Unsupported base_elem dim = " << base_elem.dim());
106  }
107 }
virtual unsigned int n_shape_functions() const override
Definition: fe.C:75

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