libMesh::InfFE< Dim, T_radial, T_map >::Base Class Reference

#include <inf_fe.h>

Static Public Member Functions

static Elembuild_elem (const Elem *inf_elem)
 
static ElemType get_elem_type (const ElemType type)
 
static unsigned int n_base_mapping_sf (const ElemType base_elem_type, const Order base_mapping_order)
 

Private Member Functions

 Base ()
 

Detailed Description

template<unsigned int Dim, FEFamily T_radial, InfMapType T_map>
class libMesh::InfFE< Dim, T_radial, T_map >::Base

This nested class contains most of the static methods related to the base part of an infinite element. Only static members are provided, and these should only be accessible from within InfFE.

Author
Daniel Dreyer
Date
2003

Definition at line 185 of file inf_fe.h.

Constructor & Destructor Documentation

◆ Base()

template<unsigned int Dim, FEFamily T_radial, InfMapType T_map>
libMesh::InfFE< Dim, T_radial, T_map >::Base::Base ( )
inlineprivate

Never use an object of this type.

Definition at line 192 of file inf_fe.h.

192 {}

Member Function Documentation

◆ build_elem()

template<unsigned int Dim, FEFamily T_radial, InfMapType T_base>
Elem * libMesh::InfFE< Dim, T_radial, T_base >::Base::build_elem ( const Elem inf_elem)
static

Build the base element of an infinite element. Be careful, this method allocates memory! So be sure to delete the new element afterward.

Definition at line 36 of file inf_fe_base_radial.C.

References libMesh::Elem::build_side_ptr().

37 {
38  std::unique_ptr<const Elem> ape(inf_elem->build_side_ptr(0));
39 
40  // The incoming inf_elem is const, but this function is required to
41  // return a non-const Elem * so that it can be used by
42  // update_base_elem(). Therefore a const_cast seems to be
43  // unavoidable here.
44  return const_cast<Elem *>(ape.release());
45 }

◆ get_elem_type()

template<unsigned int Dim, FEFamily T_radial, InfMapType T_base>
ElemType libMesh::InfFE< Dim, T_radial, T_base >::Base::get_elem_type ( const ElemType  type)
static
Returns
The base element associated to type. This is, for example, TRI3 for INFPRISM6.

Definition at line 51 of file inf_fe_base_radial.C.

References libMesh::EDGE2, libMesh::EDGE3, 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.

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

◆ n_base_mapping_sf()

template<unsigned int Dim, FEFamily T_radial, InfMapType T_base>
unsigned int libMesh::InfFE< Dim, T_radial, T_base >::Base::n_base_mapping_sf ( const ElemType  base_elem_type,
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 98 of file inf_fe_base_radial.C.

References libMesh::FE< Dim, T >::n_shape_functions().

100 {
101  if (Dim == 1)
102  return 1;
103 
104  else if (Dim == 2)
105  return FE<1,LAGRANGE>::n_shape_functions (base_elem_type,
106  base_mapping_order);
107  else if (Dim == 3)
108  return FE<2,LAGRANGE>::n_shape_functions (base_elem_type,
109  base_mapping_order);
110  else
111  {
112  libmesh_error_msg("Unsupported Dim = " << Dim);
113  return 0;
114  }
115 }
virtual unsigned int n_shape_functions() const override
Definition: fe.C:50

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