libMesh::FETransformationBase< T > Class Template Referenceabstract

#include <fe_base.h>

Public Member Functions

 FETransformationBase ()
 
virtual ~FETransformationBase ()
 
virtual void map_phi (const unsigned int dim, const Elem *const elem, const std::vector< Point > &qp, const FEGenericBase< OutputShape > &fe, std::vector< std::vector< OutputShape > > &phi) const =0
 
virtual void map_dphi (const unsigned int dim, const Elem *const elem, const std::vector< Point > &qp, const FEGenericBase< OutputShape > &fe, std::vector< std::vector< typename FEGenericBase< OutputShape >::OutputGradient > > &dphi, std::vector< std::vector< OutputShape > > &dphidx, std::vector< std::vector< OutputShape > > &dphidy, std::vector< std::vector< OutputShape > > &dphidz) const =0
 
virtual void map_d2phi (const unsigned int dim, const Elem *const elem, const std::vector< Point > &qp, const FEGenericBase< OutputShape > &fe, std::vector< std::vector< typename FEGenericBase< OutputShape >::OutputTensor > > &d2phi, std::vector< std::vector< OutputShape > > &d2phidx2, std::vector< std::vector< OutputShape > > &d2phidxdy, std::vector< std::vector< OutputShape > > &d2phidxdz, std::vector< std::vector< OutputShape > > &d2phidy2, std::vector< std::vector< OutputShape > > &d2phidydz, std::vector< std::vector< OutputShape > > &d2phidz2) const =0
 
virtual void map_curl (const unsigned int dim, const Elem *const elem, const std::vector< Point > &qp, const FEGenericBase< OutputShape > &fe, std::vector< std::vector< OutputShape > > &curl_phi) const =0
 
virtual void map_div (const unsigned int dim, const Elem *const elem, const std::vector< Point > &qp, const FEGenericBase< OutputShape > &fe, std::vector< std::vector< typename FEGenericBase< OutputShape >::OutputDivergence > > &div_phi) const =0
 

Static Public Member Functions

static AutoPtr
< FETransformationBase
< OutputShape > > 
build (const FEType &type)
 

Detailed Description

template<typename T>
class libMesh::FETransformationBase< T >

This class handles the computation of the shape functions in the physical domain. Derived classes implement the particular mapping: H1, HCurl, HDiv, or L2. This class assumes the FEGenericBase object has been initialized in the reference domain (i.e. init_shape_functions has been called).

Author
Paul T. Bauman, 2012

Definition at line 54 of file fe_base.h.

Constructor & Destructor Documentation

template<typename T>
libMesh::FETransformationBase< T >::FETransformationBase ( )
inline

Definition at line 45 of file fe_transformation_base.h.

45 {}
template<typename T>
virtual libMesh::FETransformationBase< T >::~FETransformationBase ( )
inlinevirtual

Definition at line 46 of file fe_transformation_base.h.

46 {}

Member Function Documentation

template<typename OutputShape >
AutoPtr< FETransformationBase< OutputShape > > libMesh::FETransformationBase< OutputShape >::build ( const FEType type)
static

Builds an FETransformation object based on the finite element type

Definition at line 26 of file fe_transformation_base.C.

References libMeshEnums::BERNSTEIN, libMeshEnums::CLOUGH, libMesh::FEType::family, libMeshEnums::HERMITE, libMeshEnums::HIERARCHIC, libMeshEnums::JACOBI_20_00, libMeshEnums::JACOBI_30_00, libMeshEnums::L2_HIERARCHIC, libMeshEnums::L2_LAGRANGE, libMeshEnums::LAGRANGE, libMeshEnums::LAGRANGE_VEC, libMeshEnums::MONOMIAL, libMeshEnums::NEDELEC_ONE, libMeshEnums::SCALAR, libMeshEnums::SZABAB, and libMeshEnums::XYZ.

27  {
28  switch( fe_type.family )
29  {
30  /* H1 Conforming Elements */
31  case LAGRANGE:
32  case HIERARCHIC:
33  case BERNSTEIN:
34  case SZABAB:
35  case CLOUGH: /* PB: Really H2 */
36  case HERMITE: /* PB: Really H2 */
37  case LAGRANGE_VEC:
38  case MONOMIAL: /* PB: Shouldn't this be L2 conforming? */
39  case XYZ: /* PB: Shouldn't this be L2 conforming? */
40  case L2_HIERARCHIC: /* PB: Shouldn't this be L2 conforming? */
41  case L2_LAGRANGE: /* PB: Shouldn't this be L2 conforming? */
42  case JACOBI_20_00: /* PB: For infinite elements... */
43  case JACOBI_30_00: /* PB: For infinite elements... */
44  {
45  AutoPtr<FETransformationBase<OutputShape> > ap( new H1FETransformation<OutputShape> );
46  return ap;
47  }
48  /* HCurl Conforming Elements */
49  case NEDELEC_ONE:
50  {
51  AutoPtr<FETransformationBase<OutputShape> > ap( new HCurlFETransformation<OutputShape> );
52  return ap;
53  }
54 
55  /* HDiv Conforming Elements */
56  /* L2 Conforming Elements */
57 
58  /* Other... */
59  case SCALAR:
60  {
61  // Should never need this for SCALARs
62  AutoPtr<FETransformationBase<OutputShape> > ap( new H1FETransformation<OutputShape> );
63  return ap;
64  }
65 
66  default:
67  libmesh_error();
68  }
69 
70  // Should never get here...
71  libmesh_error();
72  AutoPtr<FETransformationBase<OutputShape> > ap( NULL );
73  return ap;
74  }
template<typename T>
virtual void libMesh::FETransformationBase< T >::map_curl ( const unsigned int  dim,
const Elem *const  elem,
const std::vector< Point > &  qp,
const FEGenericBase< OutputShape > &  fe,
std::vector< std::vector< OutputShape > > &  curl_phi 
) const
pure virtual

Evaluates the shape function curl in physical coordinates based on proper finite element transformation.

Implemented in libMesh::HCurlFETransformation< T >, and libMesh::H1FETransformation< T >.

template<typename T>
virtual void libMesh::FETransformationBase< T >::map_d2phi ( const unsigned int  dim,
const Elem *const  elem,
const std::vector< Point > &  qp,
const FEGenericBase< OutputShape > &  fe,
std::vector< std::vector< typename FEGenericBase< OutputShape >::OutputTensor > > &  d2phi,
std::vector< std::vector< OutputShape > > &  d2phidx2,
std::vector< std::vector< OutputShape > > &  d2phidxdy,
std::vector< std::vector< OutputShape > > &  d2phidxdz,
std::vector< std::vector< OutputShape > > &  d2phidy2,
std::vector< std::vector< OutputShape > > &  d2phidydz,
std::vector< std::vector< OutputShape > > &  d2phidz2 
) const
pure virtual

Evaluates shape function Hessians in physical coordinates based on proper finite element transformation.

Implemented in libMesh::HCurlFETransformation< T >, and libMesh::H1FETransformation< T >.

template<typename T>
virtual void libMesh::FETransformationBase< T >::map_div ( const unsigned int  dim,
const Elem *const  elem,
const std::vector< Point > &  qp,
const FEGenericBase< OutputShape > &  fe,
std::vector< std::vector< typename FEGenericBase< OutputShape >::OutputDivergence > > &  div_phi 
) const
pure virtual

Evaluates the shape function divergence in physical coordinates based on proper finite element transformation.

Implemented in libMesh::HCurlFETransformation< T >, and libMesh::H1FETransformation< T >.

template<typename T>
virtual void libMesh::FETransformationBase< T >::map_dphi ( const unsigned int  dim,
const Elem *const  elem,
const std::vector< Point > &  qp,
const FEGenericBase< OutputShape > &  fe,
std::vector< std::vector< typename FEGenericBase< OutputShape >::OutputGradient > > &  dphi,
std::vector< std::vector< OutputShape > > &  dphidx,
std::vector< std::vector< OutputShape > > &  dphidy,
std::vector< std::vector< OutputShape > > &  dphidz 
) const
pure virtual

Evaluates shape function gradients in physical coordinates based on proper finite element transformation.

Implemented in libMesh::HCurlFETransformation< T >, and libMesh::H1FETransformation< T >.

template<typename T>
virtual void libMesh::FETransformationBase< T >::map_phi ( const unsigned int  dim,
const Elem *const  elem,
const std::vector< Point > &  qp,
const FEGenericBase< OutputShape > &  fe,
std::vector< std::vector< OutputShape > > &  phi 
) const
pure virtual

Evaluates shape functions in physical coordinates based on proper finite element transformation.

Implemented in libMesh::HCurlFETransformation< T >, and libMesh::H1FETransformation< T >.


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

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

Hosted By:
SourceForge.net Logo