libMesh::InfFE< Dim, T_radial, T_map >::Base Class Reference
#include <inf_fe.h>
Static Public Member Functions | |
| static Elem * | build_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.
- Date:
- 2003
Definition at line 188 of file inf_fe.h.
Constructor & Destructor Documentation
| libMesh::InfFE< Dim, T_radial, T_map >::Base::Base | ( | ) | [inline, private] |
Member Function Documentation
| Elem * libMesh::InfFE< Dim, T_radial, T_base >::Base::build_elem | ( | const Elem * | inf_elem | ) | [inline, static] |
Build the base element of an infinite element. Be careful, this method allocates memory! So be sure to delete the new element afterwards.
Definition at line 36 of file inf_fe_base_radial.C.
References libMesh::Elem::build_side().
Referenced by libMesh::InfFE< Dim, T_radial, T_map >::inverse_map(), libMesh::InfFE< Dim, T_radial, T_map >::map(), and libMesh::InfFE< Dim, T_radial, T_map >::update_base_elem().
| ElemType libMesh::InfFE< Dim, T_radial, T_base >::Base::get_elem_type | ( | const ElemType | type | ) | [inline, static] |
- Returns:
- the base element associated to
type. This is, for example,TRI3forINFPRISM6.
Definition at line 46 of file inf_fe_base_radial.C.
References libMeshEnums::EDGE2, libMeshEnums::EDGE3, libMesh::err, libMeshEnums::INFEDGE2, libMeshEnums::INFHEX16, libMeshEnums::INFHEX18, libMeshEnums::INFHEX8, libMeshEnums::INFPRISM12, libMeshEnums::INFPRISM6, libMeshEnums::INFQUAD4, libMeshEnums::INFQUAD6, libMeshEnums::INVALID_ELEM, libMeshEnums::QUAD4, libMeshEnums::QUAD8, libMeshEnums::QUAD9, libMeshEnums::TRI3, and libMeshEnums::TRI6.
Referenced by libMesh::InfFE< Dim, T_radial, T_map >::combine_base_radial(), 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(), and libMesh::InfFE< Dim, T_radial, T_map >::shape().
00047 { 00048 switch (type) 00049 { 00050 // 3D infinite elements: 00051 // with Dim=3 -> infinite elements on their own 00052 case INFHEX8: 00053 return QUAD4; 00054 00055 case INFHEX16: 00056 return QUAD8; 00057 00058 case INFHEX18: 00059 return QUAD9; 00060 00061 case INFPRISM6: 00062 return TRI3; 00063 00064 case INFPRISM12: 00065 return TRI6; 00066 00067 // 2D infinite elements: 00068 // with Dim=3 -> used as boundary condition, 00069 // with Dim=2 -> infinite elements on their own 00070 case INFQUAD4: 00071 return EDGE2; 00072 00073 case INFQUAD6: 00074 return EDGE3; 00075 00076 // 1D infinite elements: 00077 // with Dim=2 -> used as boundary condition, 00078 // with Dim=1 -> infinite elements on their own, 00079 // but no base element! 00080 case INFEDGE2: 00081 return INVALID_ELEM; 00082 00083 default: 00084 { 00085 libMesh::err << "ERROR: Unsupported element type!: " << type 00086 << std::endl; 00087 libmesh_error(); 00088 } 00089 } 00090 00091 00092 libmesh_error(); 00093 return INVALID_ELEM; 00094 }
| unsigned int libMesh::InfFE< Dim, T_radial, T_base >::Base::n_base_mapping_sf | ( | const ElemType | base_elem_type, | |
| const Order | base_mapping_order | |||
| ) | [inline, static] |
- Returns:
- the number of shape functions used in the mapping in the base element of type
base_elem_typemapped with orderbase_mapping_order
Definition at line 101 of file inf_fe_base_radial.C.
References libMesh::InfFE< Dim, T_radial, T_map >::n_shape_functions().
Referenced by libMesh::InfFE< Dim, T_radial, T_map >::init_face_shape_functions(), libMesh::InfFE< Dim, T_radial, T_map >::init_shape_functions(), and libMesh::InfFE< Dim, T_radial, T_map >::inverse_map().
00103 { 00104 if (Dim == 1) 00105 return 1; 00106 00107 else if (Dim == 2) 00108 return FE<1,LAGRANGE>::n_shape_functions (base_elem_type, 00109 base_mapping_order); 00110 else if (Dim == 3) 00111 return FE<2,LAGRANGE>::n_shape_functions (base_elem_type, 00112 base_mapping_order); 00113 else 00114 { 00115 // whoa, cool infinite element! 00116 libmesh_error(); 00117 return 0; 00118 } 00119 }
The documentation for this class was generated from the following files:
Site Created By: libMesh Developers
Last modified: February 05 2013 19:55:27 UTC
Hosted By: