#include <fe_transformation_base.h>
Public Member Functions | |
HCurlFETransformation () | |
virtual | ~HCurlFETransformation () |
virtual void | init_map_phi (const FEGenericBase< OutputShape > &fe) const override |
virtual void | init_map_dphi (const FEGenericBase< OutputShape > &fe) const override |
virtual void | init_map_d2phi (const FEGenericBase< OutputShape > &fe) const override |
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 override |
virtual void | map_dphi (const unsigned int, const Elem *const, const std::vector< Point > &, const FEGenericBase< OutputShape > &, std::vector< std::vector< typename FEGenericBase< OutputShape >::OutputGradient >> &, std::vector< std::vector< OutputShape >> &, std::vector< std::vector< OutputShape >> &, std::vector< std::vector< OutputShape >> &) const override |
virtual void | map_d2phi (const unsigned int, const std::vector< Point > &, const FEGenericBase< OutputShape > &, std::vector< std::vector< typename FEGenericBase< OutputShape >::OutputTensor >> &, std::vector< std::vector< OutputShape >> &, std::vector< std::vector< OutputShape >> &, std::vector< std::vector< OutputShape >> &, std::vector< std::vector< OutputShape >> &, std::vector< std::vector< OutputShape >> &, std::vector< std::vector< OutputShape >> &) const override |
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 override |
virtual void | map_div (const unsigned int, const Elem *const, const std::vector< Point > &, const FEGenericBase< OutputShape > &, std::vector< std::vector< typename FEGenericBase< OutputShape >::OutputDivergence >> &) const override |
template<> | |
void | init_map_phi (const FEGenericBase< Real > &) const |
template<> | |
void | init_map_dphi (const FEGenericBase< Real > &) const |
template<> | |
void | init_map_d2phi (const FEGenericBase< Real > &) const |
template<> | |
void | map_phi (const unsigned int, const Elem *const, const std::vector< Point > &, const FEGenericBase< Real > &, std::vector< std::vector< Real >> &) const |
template<> | |
void | map_curl (const unsigned int, const Elem *const, const std::vector< Point > &, const FEGenericBase< Real > &, std::vector< std::vector< Real >> &) const |
Static Public Member Functions | |
static std::unique_ptr< FETransformationBase< OutputShape > > | build (const FEType &type) |
This class handles the computation of the shape functions in the physical domain for HCurl conforming elements. This class assumes the FEGenericBase
object has been initialized in the reference domain (i.e. init_shape_functions
has been called).
Definition at line 29 of file fe_transformation_base.h.
|
inline |
Definition at line 40 of file hcurl_fe_transformation.h.
|
inlinevirtual |
Definition at line 43 of file hcurl_fe_transformation.h.
|
staticinherited |
Builds an FETransformation object based on the finite element type
Definition at line 28 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.
|
overridevirtual |
Pre-requests any necessary data from FEMap
Implements libMesh::FETransformationBase< OutputShape >.
Definition at line 42 of file hcurl_fe_transformation.C.
References libMesh::FEMap::get_d2xidxyz2(), libMesh::FEMap::get_dxidx(), and libMesh::FEAbstract::get_fe_map().
void libMesh::HCurlFETransformation< Real >::init_map_d2phi | ( | const FEGenericBase< Real > & | ) | const |
Definition at line 267 of file hcurl_fe_transformation.C.
|
overridevirtual |
Pre-requests any necessary data from FEMap
Implements libMesh::FETransformationBase< OutputShape >.
Definition at line 34 of file hcurl_fe_transformation.C.
References libMesh::FEMap::get_dxidx(), and libMesh::FEAbstract::get_fe_map().
void libMesh::HCurlFETransformation< Real >::init_map_dphi | ( | const FEGenericBase< Real > & | ) | const |
Definition at line 261 of file hcurl_fe_transformation.C.
|
overridevirtual |
Pre-requests any necessary data from FEMap
Implements libMesh::FETransformationBase< OutputShape >.
Definition at line 26 of file hcurl_fe_transformation.C.
References libMesh::FEMap::get_dxidx(), and libMesh::FEAbstract::get_fe_map().
void libMesh::HCurlFETransformation< Real >::init_map_phi | ( | const FEGenericBase< Real > & | ) | const |
Definition at line 255 of file hcurl_fe_transformation.C.
|
overridevirtual |
Evaluates the shape function curl in physical coordinates based on conforming finite element transformation. In 2-D, the transformation is where In 3-D, the transformation is
Implements libMesh::FETransformationBase< OutputShape >.
Definition at line 152 of file hcurl_fe_transformation.C.
References libMesh::FEGenericBase< OutputType >::get_dphideta(), libMesh::FEGenericBase< OutputType >::get_dphidxi(), libMesh::FEGenericBase< OutputType >::get_dphidzeta(), libMesh::FEMap::get_dxyzdeta(), libMesh::FEMap::get_dxyzdxi(), libMesh::FEMap::get_dxyzdzeta(), libMesh::FEAbstract::get_fe_map(), libMesh::FEMap::get_jacobian(), libMesh::index_range(), and libMesh::Real.
void libMesh::HCurlFETransformation< Real >::map_curl | ( | const unsigned int | , |
const Elem * | const, | ||
const std::vector< Point > & | , | ||
const FEGenericBase< Real > & | , | ||
std::vector< std::vector< Real >> & | |||
) | const |
Definition at line 283 of file hcurl_fe_transformation.C.
|
inlineoverridevirtual |
Evaluates shape function Hessians in physical coordinates based on conforming finite element transformation.
Implements libMesh::FETransformationBase< OutputShape >.
Definition at line 95 of file hcurl_fe_transformation.h.
|
inlineoverridevirtual |
Evaluates the shape function divergence in physical coordinates based on conforming finite element transformation.
Implements libMesh::FETransformationBase< OutputShape >.
Definition at line 127 of file hcurl_fe_transformation.h.
|
inlineoverridevirtual |
Evaluates shape function gradients in physical coordinates for conforming elements.
Implements libMesh::FETransformationBase< OutputShape >.
Definition at line 78 of file hcurl_fe_transformation.h.
|
overridevirtual |
Evaluates shape functions in physical coordinates for conforming elements. In this case , where is the inverse-transpose of the Jacobian matrix of the element map.
Implements libMesh::FETransformationBase< OutputShape >.
Definition at line 53 of file hcurl_fe_transformation.C.
References libMesh::FEMap::get_detadx(), libMesh::FEMap::get_detady(), libMesh::FEMap::get_detadz(), libMesh::FEMap::get_dxidx(), libMesh::FEMap::get_dxidy(), libMesh::FEMap::get_dxidz(), libMesh::FEMap::get_dzetadx(), libMesh::FEMap::get_dzetady(), libMesh::FEMap::get_dzetadz(), libMesh::FEAbstract::get_fe_map(), libMesh::FEAbstract::get_fe_type(), and libMesh::index_range().
void libMesh::HCurlFETransformation< Real >::map_phi | ( | const unsigned int | , |
const Elem * | const, | ||
const std::vector< Point > & | , | ||
const FEGenericBase< Real > & | , | ||
std::vector< std::vector< Real >> & | |||
) | const |
Definition at line 273 of file hcurl_fe_transformation.C.