libMesh::InfFE< friend_Dim, friend_T_radial, friend_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 friend_Dim, FEFamily friend_T_radial, InfMapType friend_T_map>
class libMesh::InfFE< friend_Dim, friend_T_radial, friend_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 188 of file inf_fe.h.

Constructor & Destructor Documentation

template<unsigned int friend_Dim, FEFamily friend_T_radial, InfMapType friend_T_map>
libMesh::InfFE< friend_Dim, friend_T_radial, friend_T_map >::Base::Base ( )
inlineprivate

Never use an object of this type.

Definition at line 195 of file inf_fe.h.

195 {}

Member Function Documentation

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 afterwards.

Definition at line 36 of file inf_fe_base_radial.C.

References libMesh::Elem::build_side(), and libMesh::AutoPtr< Tp >::release().

37 {
38  AutoPtr<Elem> ape(inf_elem->build_side(0));
39  return ape.release();
40 }
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 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.

47 {
48  switch (type)
49  {
50  // 3D infinite elements:
51  // with Dim=3 -> infinite elements on their own
52  case INFHEX8:
53  return QUAD4;
54 
55  case INFHEX16:
56  return QUAD8;
57 
58  case INFHEX18:
59  return QUAD9;
60 
61  case INFPRISM6:
62  return TRI3;
63 
64  case INFPRISM12:
65  return TRI6;
66 
67  // 2D infinite elements:
68  // with Dim=3 -> used as boundary condition,
69  // with Dim=2 -> infinite elements on their own
70  case INFQUAD4:
71  return EDGE2;
72 
73  case INFQUAD6:
74  return EDGE3;
75 
76  // 1D infinite elements:
77  // with Dim=2 -> used as boundary condition,
78  // with Dim=1 -> infinite elements on their own,
79  // but no base element!
80  case INFEDGE2:
81  return INVALID_ELEM;
82 
83  default:
84  {
85  libMesh::err << "ERROR: Unsupported element type!: " << type
86  << std::endl;
87  libmesh_error();
88  }
89  }
90 
91 
92  libmesh_error();
93  return INVALID_ELEM;
94 }
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 101 of file inf_fe_base_radial.C.

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

103 {
104  if (Dim == 1)
105  return 1;
106 
107  else if (Dim == 2)
108  return FE<1,LAGRANGE>::n_shape_functions (base_elem_type,
109  base_mapping_order);
110  else if (Dim == 3)
111  return FE<2,LAGRANGE>::n_shape_functions (base_elem_type,
112  base_mapping_order);
113  else
114  {
115  // whoa, cool infinite element!
116  libmesh_error();
117  return 0;
118  }
119 }

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

Site Created By: libMesh Developers
Last modified: February 07 2014 16:57:22 UTC

Hosted By:
SourceForge.net Logo