libMesh
Public Member Functions | Static Public Member Functions | List of all members
libMesh::FETransformationBase< OutputShape > Class Template Referenceabstract

This class handles the computation of the shape functions in the physical domain. More...

#include <fe_base.h>

Inheritance diagram for libMesh::FETransformationBase< OutputShape >:
[legend]

Public Member Functions

 FETransformationBase ()=default
 
virtual ~FETransformationBase ()=default
 
virtual void init_map_phi (const FEGenericBase< OutputShape > &fe) const =0
 Pre-requests any necessary data from FEMap. More...
 
virtual void init_map_dphi (const FEGenericBase< OutputShape > &fe) const =0
 Pre-requests any necessary data from FEMap. More...
 
virtual void init_map_d2phi (const FEGenericBase< OutputShape > &fe) const =0
 Pre-requests any necessary data from FEMap. More...
 
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, bool add_p_level=true) const =0
 Evaluates shape functions in physical coordinates based on proper finite element transformation. More...
 
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
 Evaluates shape function gradients in physical coordinates based on proper finite element transformation. More...
 
virtual void map_d2phi (const unsigned int dim, 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
 Evaluates shape function Hessians in physical coordinates based on proper finite element transformation. More...
 
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
 Evaluates the shape function curl in physical coordinates based on proper finite element transformation. More...
 
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
 Evaluates the shape function divergence in physical coordinates based on proper finite element transformation. More...
 

Static Public Member Functions

static std::unique_ptr< FETransformationBase< OutputShape > > build (const FEType &type)
 Builds an FETransformation object based on the finite element type. More...
 

Detailed Description

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

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
Date
2012

Definition at line 54 of file fe_base.h.

Constructor & Destructor Documentation

◆ FETransformationBase()

template<typename OutputShape >
libMesh::FETransformationBase< OutputShape >::FETransformationBase ( )
default

◆ ~FETransformationBase()

template<typename OutputShape >
virtual libMesh::FETransformationBase< OutputShape >::~FETransformationBase ( )
virtualdefault

Member Function Documentation

◆ build()

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

Builds an FETransformation object based on the finite element type.

Definition at line 32 of file fe_transformation_base.C.

References libMesh::BERNSTEIN, libMesh::CLOUGH, libMesh::Utility::enum_to_string(), libMesh::FEType::family, libMesh::HERMITE, libMesh::HIERARCHIC, libMesh::HIERARCHIC_VEC, libMesh::JACOBI_20_00, libMesh::JACOBI_30_00, libMesh::L2_HIERARCHIC, libMesh::L2_HIERARCHIC_VEC, libMesh::L2_LAGRANGE, libMesh::L2_LAGRANGE_VEC, libMesh::L2_RAVIART_THOMAS, libMesh::LAGRANGE, libMesh::LAGRANGE_VEC, libMesh::MONOMIAL, libMesh::MONOMIAL_VEC, libMesh::NEDELEC_ONE, libMesh::RATIONAL_BERNSTEIN, libMesh::RAVIART_THOMAS, libMesh::SCALAR, libMesh::SIDE_HIERARCHIC, libMesh::SUBDIVISION, libMesh::SZABAB, and libMesh::XYZ.

33 {
34  switch (fe_type.family)
35  {
36  // H1 Conforming Elements
37  case LAGRANGE:
38  case HIERARCHIC:
39  case BERNSTEIN:
40  case SZABAB:
41  case CLOUGH: // PB: Really H2
42  case HERMITE: // PB: Really H2
43  case SUBDIVISION:
44  case LAGRANGE_VEC:
45  case HIERARCHIC_VEC:
46  case MONOMIAL: // PB: Shouldn't this be L2 conforming?
47  case MONOMIAL_VEC: // PB: Shouldn't this be L2 conforming?
48  case XYZ: // PB: Shouldn't this be L2 conforming?
49  case RATIONAL_BERNSTEIN:
50  case L2_HIERARCHIC:
51  case L2_HIERARCHIC_VEC:
52  case SIDE_HIERARCHIC:
53  case L2_LAGRANGE: // PB: Shouldn't this be L2 conforming?
54  case L2_LAGRANGE_VEC: // PB: Shouldn't this be L2 conforming?
55  case JACOBI_20_00: // PB: For infinite elements...
56  case JACOBI_30_00: // PB: For infinite elements...
57  return std::make_unique<H1FETransformation<OutputShape>>();
58 
59  // HCurl Conforming Elements
60  case NEDELEC_ONE:
61  return std::make_unique<HCurlFETransformation<OutputShape>>();
62 
63  // HDiv Conforming Elements
64  case RAVIART_THOMAS:
65  case L2_RAVIART_THOMAS:
66  return std::make_unique<HDivFETransformation<OutputShape>>();
67 
68  // L2 Conforming Elements
69 
70  // Other...
71  case SCALAR:
72  // Should never need this for SCALARs
73  return std::make_unique<H1FETransformation<OutputShape>>();
74 
75  default:
76  libmesh_error_msg("Unknown family = " << Utility::enum_to_string(fe_type.family));
77  }
78 }
std::string enum_to_string(const T e)

◆ init_map_d2phi()

template<typename OutputShape >
virtual void libMesh::FETransformationBase< OutputShape >::init_map_d2phi ( const FEGenericBase< OutputShape > &  fe) const
pure virtual

◆ init_map_dphi()

template<typename OutputShape >
virtual void libMesh::FETransformationBase< OutputShape >::init_map_dphi ( const FEGenericBase< OutputShape > &  fe) const
pure virtual

◆ init_map_phi()

template<typename OutputShape >
virtual void libMesh::FETransformationBase< OutputShape >::init_map_phi ( const FEGenericBase< OutputShape > &  fe) const
pure virtual

◆ map_curl()

template<typename OutputShape >
virtual void libMesh::FETransformationBase< OutputShape >::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< OutputShape >, libMesh::HDivFETransformation< OutputShape >, and libMesh::H1FETransformation< OutputShape >.

◆ map_d2phi()

template<typename OutputShape >
virtual void libMesh::FETransformationBase< OutputShape >::map_d2phi ( const unsigned int  dim,
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< OutputShape >, libMesh::HDivFETransformation< OutputShape >, and libMesh::H1FETransformation< OutputShape >.

◆ map_div()

template<typename OutputShape >
virtual void libMesh::FETransformationBase< OutputShape >::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::HDivFETransformation< OutputShape >, libMesh::HCurlFETransformation< OutputShape >, and libMesh::H1FETransformation< OutputShape >.

◆ map_dphi()

template<typename OutputShape >
virtual void libMesh::FETransformationBase< OutputShape >::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::H1FETransformation< OutputShape >, libMesh::HCurlFETransformation< OutputShape >, and libMesh::HDivFETransformation< OutputShape >.

◆ map_phi()

template<typename OutputShape >
virtual void libMesh::FETransformationBase< OutputShape >::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,
bool  add_p_level = true 
) const
pure virtual

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

Implemented in libMesh::H1FETransformation< OutputShape >, libMesh::HCurlFETransformation< OutputShape >, and libMesh::HDivFETransformation< OutputShape >.


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