Poster of Linux kernelThe best gift for a Linux geek
InfFE::Base

InfFE::Base

Section: C Library Functions (3) Updated: Thu Apr 7 2011
Local index Up
 

NAME

InfFE::Base -  

SYNOPSIS


#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 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

Version:

Revision:

3391

Definition at line 188 of file inf_fe.h.  

Constructor & Destructor Documentation

 

template<unsigned int Dim, FEFamily T_radial, InfMapType T_map> InfFE< Dim, T_radial, T_map >::Base::Base () [inline, private]Never use an object of this type.

Definition at line 195 of file inf_fe.h.

{}
 

Member Function Documentation

 

template<unsigned int Dim, FEFamily T_radial, InfMapType T_base> Elem * 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 34 of file inf_fe_base_radial.C.

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

Referenced by InfFE< Dim, T_radial, T_map >::inverse_map(), and InfFE< Dim, T_radial, T_map >::map().

{ 
  AutoPtr<Elem> ape(inf_elem->build_side(0)); 
  return ape.release(); 
}
 

template<unsigned int Dim, FEFamily T_radial, InfMapType T_base> ElemType InfFE< Dim, T_radial, T_base >::Base::get_elem_type (const ElemTypetype) [static]Returns:

the base element associated to type. This is, for example, TRI3 for INFPRISM6.

Definition at line 44 of file inf_fe_base_radial.C.

References libMeshEnums::EDGE2, libMeshEnums::EDGE3, 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 InfFE< Dim, T_radial, T_map >::compute_shape_indices(), InfFE< Dim, T_radial, T_map >::n_dofs(), InfFE< Dim, T_radial, T_map >::n_dofs_at_node(), InfFE< Dim, T_radial, T_map >::n_dofs_per_elem(), and InfFE< Dim, T_radial, T_map >::shape().

{
  switch (type)
    {
      // 3D infinite elements:
      // with Dim=3 -> infinite elements on their own
      case INFHEX8:
          return QUAD4;

      case INFHEX16:
          return QUAD8;

      case INFHEX18:
          return QUAD9;
                 
      case INFPRISM6:
          return TRI3;

      case INFPRISM12:
          return TRI6;

      // 2D infinite elements:
      // with Dim=3 -> used as boundary condition,
      // with Dim=2 -> infinite elements on their own
      case INFQUAD4:
          return EDGE2;

      case INFQUAD6:
          return EDGE3;

      // 1D infinite elements:
      // with Dim=2 -> used as boundary condition,
      // with Dim=1 -> infinite elements on their own,
      //               but no base element!
      case INFEDGE2:
          return INVALID_ELEM;

      default:
        {
          std::cerr << 'ERROR: Unsupported element type!: ' << type
                    << std::endl;
          libmesh_error();
        }
    }


  libmesh_error();
  return INVALID_ELEM;
}
 

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

References InfFE< Dim, T_radial, T_map >::n_shape_functions().

Referenced by InfFE< Dim, T_radial, T_map >::init_face_shape_functions(), and InfFE< Dim, T_radial, T_map >::inverse_map().

{
  if (Dim == 1)
    return 1;
  
  else if (Dim == 2)
    return FE<1,LAGRANGE>::n_shape_functions (base_elem_type,
                                              base_mapping_order);
  else if (Dim == 3)
    return FE<2,LAGRANGE>::n_shape_functions (base_elem_type,
                                              base_mapping_order);
  else
    {
      // whoa, cool infinite element!
      libmesh_error();
      return 0;
    }
}

 

Author

Generated automatically by Doxygen for libMesh from the source code.


 

Index

NAME
SYNOPSIS
Static Public Member Functions
Private Member Functions
Detailed Description
template<unsigned int Dim, FEFamily T_radial, InfMapType T_map> class InfFE< Dim, T_radial, T_map >::Base
Constructor & Destructor Documentation
template<unsigned int Dim, FEFamily T_radial, InfMapType T_map> InfFE< Dim, T_radial, T_map >::Base::Base () [inline, private]Never use an object of this type.
Member Function Documentation
template<unsigned int Dim, FEFamily T_radial, InfMapType T_base> Elem * 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.
template<unsigned int Dim, FEFamily T_radial, InfMapType T_base> ElemType InfFE< Dim, T_radial, T_base >::Base::get_elem_type (const ElemTypetype) [static]Returns:
template<unsigned int Dim, FEFamily T_radial, InfMapType T_base> unsigned int InfFE< Dim, T_radial, T_base >::Base::n_base_mapping_sf (const ElemTypebase_elem_type, const Orderbase_mapping_order) [static]Returns:
Author

This document was created by man2html, using the manual pages.
Time: 21:48:29 GMT, April 16, 2011