$extrastylesheet
#include <fe_transformation_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 UniquePtr < FETransformationBase < OutputShape > > | build (const FEType &type) |
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).
Definition at line 41 of file fe_transformation_base.h.
| libMesh::FETransformationBase< OutputShape >::FETransformationBase | ( | ) | [inline] |
Definition at line 45 of file fe_transformation_base.h.
{}
| virtual libMesh::FETransformationBase< OutputShape >::~FETransformationBase | ( | ) | [inline, virtual] |
Definition at line 46 of file fe_transformation_base.h.
{}
| UniquePtr< 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 libMesh::BERNSTEIN, libMesh::CLOUGH, libMesh::FEType::family, libMesh::HERMITE, libMesh::HIERARCHIC, libMesh::JACOBI_20_00, libMesh::JACOBI_30_00, libMesh::L2_HIERARCHIC, libMesh::L2_LAGRANGE, libMesh::LAGRANGE, libMesh::LAGRANGE_VEC, libMesh::MONOMIAL, libMesh::NEDELEC_ONE, libMesh::SCALAR, libMesh::SUBDIVISION, libMesh::SZABAB, and libMesh::XYZ.
{
switch (fe_type.family)
{
// H1 Conforming Elements
case LAGRANGE:
case HIERARCHIC:
case BERNSTEIN:
case SZABAB:
case CLOUGH: // PB: Really H2
case HERMITE: // PB: Really H2
case SUBDIVISION:
case LAGRANGE_VEC:
case MONOMIAL: // PB: Shouldn't this be L2 conforming?
case XYZ: // PB: Shouldn't this be L2 conforming?
case L2_HIERARCHIC: // PB: Shouldn't this be L2 conforming?
case L2_LAGRANGE: // PB: Shouldn't this be L2 conforming?
case JACOBI_20_00: // PB: For infinite elements...
case JACOBI_30_00: // PB: For infinite elements...
return UniquePtr<FETransformationBase<OutputShape> >(new H1FETransformation<OutputShape>);
// HCurl Conforming Elements
case NEDELEC_ONE:
return UniquePtr<FETransformationBase<OutputShape> >(new HCurlFETransformation<OutputShape>);
// HDiv Conforming Elements
// L2 Conforming Elements
// Other...
case SCALAR:
// Should never need this for SCALARs
return UniquePtr<FETransformationBase<OutputShape> >(new H1FETransformation<OutputShape>);
default:
libmesh_error_msg("Unknown family = " << fe_type.family);
}
libmesh_error_msg("We'll never get here!");
return UniquePtr<FETransformationBase<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 >, and libMesh::H1FETransformation< OutputShape >.
| virtual void libMesh::FETransformationBase< OutputShape >::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< OutputShape >, and libMesh::H1FETransformation< 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::HCurlFETransformation< OutputShape >, and libMesh::H1FETransformation< 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::HCurlFETransformation< OutputShape >, and libMesh::H1FETransformation< 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 | ||
| ) | const [pure virtual] |
Evaluates shape functions in physical coordinates based on proper finite element transformation.
Implemented in libMesh::HCurlFETransformation< OutputShape >, and libMesh::H1FETransformation< OutputShape >.