$extrastylesheet
libMesh::HCurlFETransformation< OutputShape > Class Template Reference

#include <hcurl_fe_transformation.h>

Inheritance diagram for libMesh::HCurlFETransformation< OutputShape >:

List of all members.

Public Member Functions

 HCurlFETransformation ()
virtual ~HCurlFETransformation ()
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
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
virtual void map_d2phi (const unsigned int, const Elem *const , 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
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
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
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 UniquePtr
< FETransformationBase
< OutputShape > > 
build (const FEType &type)

Detailed Description

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

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).

Author:
Paul T. Bauman, 2012

Definition at line 34 of file hcurl_fe_transformation.h.


Constructor & Destructor Documentation

template<typename OutputShape >
libMesh::HCurlFETransformation< OutputShape >::HCurlFETransformation ( ) [inline]

Definition at line 38 of file hcurl_fe_transformation.h.

    : FETransformationBase<OutputShape>(){}
template<typename OutputShape >
virtual libMesh::HCurlFETransformation< OutputShape >::~HCurlFETransformation ( ) [inline, virtual]

Definition at line 41 of file hcurl_fe_transformation.h.

{}

Member Function Documentation

template<typename OutputShape >
UniquePtr< FETransformationBase< OutputShape > > libMesh::FETransformationBase< OutputShape >::build ( const FEType type) [static, inherited]

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> >();
}
template<typename OutputShape >
void libMesh::HCurlFETransformation< 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 [virtual]

Evaluates the shape function curl in physical coordinates based on HCurl conforming finite element transformation. In 2-D, the transformation is $ \nabla \times \phi = J^{-1} * \nabla \times \hat{\phi} $ where $ J = \det( dx/d\xi ) $ In 3-D, the transformation is $ \nabla \times \phi = J^{-1} dx/d\xi \nabla \times \hat{\phi} $

Implements libMesh::FETransformationBase< OutputShape >.

Definition at line 126 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(), and libMesh::Real.

{
  switch (dim)
    {
    case 0:
    case 1:
      libmesh_error_msg("These element transformations only make sense in 2D and 3D.");

    case 2:
      {
        const std::vector<std::vector<OutputShape> >& dphi_dxi = fe.get_dphidxi();
        const std::vector<std::vector<OutputShape> >& dphi_deta = fe.get_dphideta();

        const std::vector<Real>& J = fe.get_fe_map().get_jacobian();

        // FIXME: I don't think this is valid for 2D elements in 3D space
        /* In 2D: curl(phi) = J^{-1} * curl(\hat{phi}) */
        for (unsigned int i=0; i<curl_phi.size(); i++)
          for (unsigned int p=0; p<curl_phi[i].size(); p++)
            {
              curl_phi[i][p].slice(0) = curl_phi[i][p].slice(1) = 0.0;

              curl_phi[i][p].slice(2) = ( dphi_dxi[i][p].slice(1) - dphi_deta[i][p].slice(0) )/J[p];
            }

        break;
      }

    case 3:
      {
        const std::vector<std::vector<OutputShape> >& dphi_dxi = fe.get_dphidxi();
        const std::vector<std::vector<OutputShape> >& dphi_deta = fe.get_dphideta();
        const std::vector<std::vector<OutputShape> >& dphi_dzeta = fe.get_dphidzeta();

        const std::vector<RealGradient>& dxyz_dxi   = fe.get_fe_map().get_dxyzdxi();
        const std::vector<RealGradient>& dxyz_deta  = fe.get_fe_map().get_dxyzdeta();
        const std::vector<RealGradient>& dxyz_dzeta = fe.get_fe_map().get_dxyzdzeta();

        const std::vector<Real>& J = fe.get_fe_map().get_jacobian();

        for (unsigned int i=0; i<curl_phi.size(); i++)
          for (unsigned int p=0; p<curl_phi[i].size(); p++)
            {
              Real dx_dxi   = dxyz_dxi[p](0);
              Real dx_deta  = dxyz_deta[p](0);
              Real dx_dzeta = dxyz_dzeta[p](0);

              Real dy_dxi   = dxyz_dxi[p](1);
              Real dy_deta  = dxyz_deta[p](1);
              Real dy_dzeta = dxyz_dzeta[p](1);

              Real dz_dxi   = dxyz_dxi[p](2);
              Real dz_deta  = dxyz_deta[p](2);
              Real dz_dzeta = dxyz_dzeta[p](2);

              const Real inv_jac = 1.0/J[p];

              /* In 3D: curl(phi) = J^{-1} dx/dxi * curl(\hat{phi})

                 dx/dxi = [  dx/dxi  dx/deta  dx/dzeta
                 dy/dxi  dy/deta  dy/dzeta
                 dz/dxi  dz/deta  dz/dzeta ]

                 curl(u) = [ du_z/deta  - du_y/dzeta
                 du_x/dzeta - du_z/dxi
                 du_y/dxi   - du_x/deta ]
              */
              curl_phi[i][p].slice(0) = inv_jac*( dx_dxi*( dphi_deta[i][p].slice(2)  -
                                                           dphi_dzeta[i][p].slice(1)   ) +
                                                  dx_deta*( dphi_dzeta[i][p].slice(0) -
                                                            dphi_dxi[i][p].slice(2)     ) +
                                                  dx_dzeta*( dphi_dxi[i][p].slice(1) -
                                                             dphi_deta[i][p].slice(0)    ) );

              curl_phi[i][p].slice(1) = inv_jac*( dy_dxi*( dphi_deta[i][p].slice(2) -
                                                           dphi_dzeta[i][p].slice(1)  ) +
                                                  dy_deta*( dphi_dzeta[i][p].slice(0)-
                                                            dphi_dxi[i][p].slice(2)    ) +
                                                  dy_dzeta*( dphi_dxi[i][p].slice(1) -
                                                             dphi_deta[i][p].slice(0)   ) );

              curl_phi[i][p].slice(2) = inv_jac*( dz_dxi*( dphi_deta[i][p].slice(2) -
                                                           dphi_dzeta[i][p].slice(1)   ) +
                                                  dz_deta*( dphi_dzeta[i][p].slice(0) -
                                                            dphi_dxi[i][p].slice(2)     ) +
                                                  dz_dzeta*( dphi_dxi[i][p].slice(1) -
                                                             dphi_deta[i][p].slice(0)    ) );
            }

        break;
      }

    default:
      libmesh_error_msg("Invalid dim = " << dim);
    } // switch(dim)
}
template<>
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 240 of file hcurl_fe_transformation.C.

{
  libmesh_error_msg("HCurl transformations only make sense for vector-valued elements.");
}
template<typename OutputShape >
virtual void libMesh::HCurlFETransformation< OutputShape >::map_d2phi ( const unsigned int  ,
const Elem const,
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 [inline, virtual]

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

Implements libMesh::FETransformationBase< OutputShape >.

Definition at line 75 of file hcurl_fe_transformation.h.

References libMesh::err.

  { libmesh_do_once( libMesh::err << "WARNING: Shape function Hessians for HCurl elements are not currently "
                     << "being computed!" << std::endl; ); return; }
template<typename OutputShape >
virtual void libMesh::HCurlFETransformation< OutputShape >::map_div ( const unsigned int  ,
const Elem const,
const std::vector< Point > &  ,
const FEGenericBase< OutputShape > &  ,
std::vector< std::vector< typename FEGenericBase< OutputShape >::OutputDivergence > > &   
) const [inline, virtual]

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

Implements libMesh::FETransformationBase< OutputShape >.

Definition at line 107 of file hcurl_fe_transformation.h.

References libMesh::err.

  { libmesh_do_once( libMesh::err << "WARNING: Shape function divergences for HCurl elements are not currently "
                     << "being computed!" << std::endl; ); return; }
template<typename OutputShape >
virtual void libMesh::HCurlFETransformation< OutputShape >::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 [inline, virtual]

Evaluates shape function gradients in physical coordinates for HCurl conforming elements.

Implements libMesh::FETransformationBase< OutputShape >.

Definition at line 59 of file hcurl_fe_transformation.h.

References libMesh::err.

  { libmesh_do_once( libMesh::err << "WARNING: Shape function gradients for HCurl elements are not currently "
                     << "being computed!" << std::endl; ); return; }
template<typename OutputShape >
void libMesh::HCurlFETransformation< 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 [virtual]

Evaluates shape functions in physical coordinates for HCurl conforming elements. In this case $ \phi = (dx/d\xi)^{-T} \hat{\phi} $, where $ (dx/d\xi)^{-T} $ is the inverse-transpose of the Jacobian matrix of the element map. Note here $ x, \xi $ are vectors.

Implements libMesh::FETransformationBase< OutputShape >.

Definition at line 25 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(), and libMesh::FEAbstract::get_fe_type().

{
  switch (dim)
    {
    case 0:
    case 1:
      libmesh_error_msg("These element transformations only make sense in 2D and 3D.");

    case 2:
      {
        const std::vector<Real>& dxidx_map = fe.get_fe_map().get_dxidx();
        const std::vector<Real>& dxidy_map = fe.get_fe_map().get_dxidy();

        const std::vector<Real>& detadx_map = fe.get_fe_map().get_detadx();
        const std::vector<Real>& detady_map = fe.get_fe_map().get_detady();

        // FIXME: Need to update for 2D elements in 3D space
        /* phi = (dx/dxi)^-T * \hat{phi}
           In 2D:
           (dx/dxi)^{-1} = [  dxi/dx   dxi/dy
           deta/dx  deta/dy ]

           so: dxi/dx^{-T} * \hat{phi} = [ dxi/dx  deta/dx   [ \hat{phi}_xi
           dxi/dy  deta/dy ]   \hat{phi}_eta ]

           or in indicial notation:  phi_j = xi_{i,j}*\hat{phi}_i */

        for (unsigned int i=0; i<phi.size(); i++)
          for (unsigned int p=0; p<phi[i].size(); p++)
            {
              // Need to temporarily cache reference shape functions
              // TODO: PB: Might be worth trying to build phi_ref separately to see
              //           if we can get vectorization
              OutputShape phi_ref;
              FEInterface::shape<OutputShape>(2, fe.get_fe_type(), elem, i, qp[p], phi_ref);

              phi[i][p](0) = dxidx_map[p]*phi_ref.slice(0) + detadx_map[p]*phi_ref.slice(1);

              phi[i][p](1) = dxidy_map[p]*phi_ref.slice(0) + detady_map[p]*phi_ref.slice(1);
            }

        break;
      }

    case 3:
      {
        const std::vector<Real>& dxidx_map = fe.get_fe_map().get_dxidx();
        const std::vector<Real>& dxidy_map = fe.get_fe_map().get_dxidy();
        const std::vector<Real>& dxidz_map = fe.get_fe_map().get_dxidz();

        const std::vector<Real>& detadx_map = fe.get_fe_map().get_detadx();
        const std::vector<Real>& detady_map = fe.get_fe_map().get_detady();
        const std::vector<Real>& detadz_map = fe.get_fe_map().get_detadz();

        const std::vector<Real>& dzetadx_map = fe.get_fe_map().get_dzetadx();
        const std::vector<Real>& dzetady_map = fe.get_fe_map().get_dzetady();
        const std::vector<Real>& dzetadz_map = fe.get_fe_map().get_dzetadz();

        /* phi = (dx/dxi)^-T * \hat{phi}
           In 3D:
           dx/dxi^-1 = [  dxi/dx    dxi/dy    dxi/dz
           deta/dx   deta/dy   deta/dz
           dzeta/dx  dzeta/dy  dzeta/dz]

           so: dxi/dx^-T * \hat{phi} = [ dxi/dx  deta/dx  dzeta/dx   [ \hat{phi}_xi
           dxi/dy  deta/dy  dzeta/dy     \hat{phi}_eta
           dxi/dz  deta/dz  dzeta/dz ]   \hat{phi}_zeta ]

           or in indicial notation:  phi_j = xi_{i,j}*\hat{phi}_i */

        for (unsigned int i=0; i<phi.size(); i++)
          for (unsigned int p=0; p<phi[i].size(); p++)
            {
              // Need to temporarily cache reference shape functions
              // TODO: PB: Might be worth trying to build phi_ref separately to see
              //           if we can get vectorization
              OutputShape phi_ref;
              FEInterface::shape<OutputShape>(3, fe.get_fe_type(), elem, i, qp[p], phi_ref);

              phi[i][p].slice(0) = dxidx_map[p]*phi_ref.slice(0) + detadx_map[p]*phi_ref.slice(1)
                + dzetadx_map[p]*phi_ref.slice(2);

              phi[i][p].slice(1) = dxidy_map[p]*phi_ref.slice(0) + detady_map[p]*phi_ref.slice(1)
                + dzetady_map[p]*phi_ref.slice(2);

              phi[i][p].slice(2) = dxidz_map[p]*phi_ref.slice(0) + detadz_map[p]*phi_ref.slice(1)
                + dzetadz_map[p]*phi_ref.slice(2);
            }
        break;
      }

    default:
      libmesh_error_msg("Invalid dim = " << dim);
    } // switch(dim)
}
template<>
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 230 of file hcurl_fe_transformation.C.

{
  libmesh_error_msg("HCurl transformations only make sense for vector-valued elements.");
}

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