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

00039       : FETransformationBase<OutputShape>(){};

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

Definition at line 41 of file hcurl_fe_transformation.h.

00041 {};


Member Function Documentation

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

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.

00027   {
00028     switch( fe_type.family )
00029       {
00030         /* H1 Conforming Elements */
00031       case LAGRANGE:
00032       case HIERARCHIC:
00033       case BERNSTEIN:
00034       case SZABAB:
00035       case CLOUGH: /* PB: Really H2 */
00036       case HERMITE: /* PB: Really H2 */
00037       case LAGRANGE_VEC:
00038       case MONOMIAL: /* PB: Shouldn't this be L2 conforming? */
00039       case XYZ: /* PB: Shouldn't this be L2 conforming? */
00040       case L2_HIERARCHIC: /* PB: Shouldn't this be L2 conforming? */
00041       case L2_LAGRANGE: /* PB: Shouldn't this be L2 conforming? */
00042       case JACOBI_20_00: /* PB: For infinite elements... */
00043       case JACOBI_30_00: /* PB: For infinite elements... */
00044         {
00045           AutoPtr<FETransformationBase<OutputShape> > ap( new H1FETransformation<OutputShape> );
00046           return ap;
00047         }
00048         /* HCurl Conforming Elements */
00049       case NEDELEC_ONE:
00050         {
00051           AutoPtr<FETransformationBase<OutputShape> > ap( new HCurlFETransformation<OutputShape> );
00052           return ap;
00053         }
00054 
00055         /* HDiv Conforming Elements */
00056         /* L2 Conforming Elements */
00057 
00058         /* Other... */
00059       case SCALAR:
00060         {
00061           // Should never need this for SCALARs
00062           AutoPtr<FETransformationBase<OutputShape> > ap( new H1FETransformation<OutputShape> );
00063           return ap;
00064         }
00065         
00066       default:
00067         libmesh_error();
00068       }
00069 
00070     // Should never get here...
00071     libmesh_error();
00072     AutoPtr<FETransformationBase<OutputShape> > ap( NULL );
00073     return ap;
00074   }

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 [inline]

Definition at line 254 of file hcurl_fe_transformation.C.

References libMesh::err.

00259 {
00260   libMesh::err << "HCurl transformations only make sense for vector-valued elements."
00261                << std::endl;
00262   libmesh_error();
00263 }

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 [inline, 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 132 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.

00137   {
00138     switch(dim)
00139       {
00140         // These element transformations only make sense in 2D and 3D
00141       case 0: 
00142       case 1:
00143         {
00144           libmesh_error();
00145         }
00146         
00147       case 2:
00148         {
00149           const std::vector<std::vector<OutputShape> >& dphi_dxi = fe.get_dphidxi();
00150           const std::vector<std::vector<OutputShape> >& dphi_deta = fe.get_dphideta();
00151 
00152           const std::vector<Real>& J = fe.get_fe_map().get_jacobian();
00153 
00154           // FIXME: I don't think this is valid for 2D elements in 3D space
00155           /* In 2D: curl(phi) = J^{-1} * curl(\hat{phi}) */
00156           for (unsigned int i=0; i<curl_phi.size(); i++)
00157             for (unsigned int p=0; p<curl_phi[i].size(); p++)
00158               {
00159                 curl_phi[i][p].slice(0) = curl_phi[i][p].slice(1) = 0.0;
00160 
00161                 curl_phi[i][p].slice(2) = ( dphi_dxi[i][p].slice(1) - dphi_deta[i][p].slice(0) )/J[p];
00162               }
00163 
00164           break;
00165         }
00166 
00167       case 3:
00168         {
00169           const std::vector<std::vector<OutputShape> >& dphi_dxi = fe.get_dphidxi();
00170           const std::vector<std::vector<OutputShape> >& dphi_deta = fe.get_dphideta();
00171           const std::vector<std::vector<OutputShape> >& dphi_dzeta = fe.get_dphidzeta();
00172           
00173           const std::vector<RealGradient>& dxyz_dxi   = fe.get_fe_map().get_dxyzdxi();
00174           const std::vector<RealGradient>& dxyz_deta  = fe.get_fe_map().get_dxyzdeta();
00175           const std::vector<RealGradient>& dxyz_dzeta = fe.get_fe_map().get_dxyzdzeta();
00176 
00177           const std::vector<Real>& J = fe.get_fe_map().get_jacobian();
00178 
00179           for (unsigned int i=0; i<curl_phi.size(); i++)
00180             for (unsigned int p=0; p<curl_phi[i].size(); p++)
00181               {
00182                 Real dx_dxi   = dxyz_dxi[p](0);
00183                 Real dx_deta  = dxyz_deta[p](0);
00184                 Real dx_dzeta = dxyz_dzeta[p](0);
00185 
00186                 Real dy_dxi   = dxyz_dxi[p](1);
00187                 Real dy_deta  = dxyz_deta[p](1);
00188                 Real dy_dzeta = dxyz_dzeta[p](1);
00189 
00190                 Real dz_dxi   = dxyz_dxi[p](2);
00191                 Real dz_deta  = dxyz_deta[p](2);
00192                 Real dz_dzeta = dxyz_dzeta[p](2);
00193 
00194                 const Real inv_jac = 1.0/J[p];
00195 
00196                 /* In 3D: curl(phi) = J^{-1} dx/dxi * curl(\hat{phi}) 
00197 
00198                    dx/dxi = [  dx/dxi  dx/deta  dx/dzeta
00199                                dy/dxi  dy/deta  dy/dzeta
00200                                dz/dxi  dz/deta  dz/dzeta ]
00201 
00202                    curl(u) = [ du_z/deta  - du_y/dzeta
00203                                du_x/dzeta - du_z/dxi
00204                                du_y/dxi   - du_x/deta ]
00205                  */
00206                 curl_phi[i][p].slice(0) = inv_jac*( dx_dxi*( dphi_deta[i][p].slice(2)  - 
00207                                                              dphi_dzeta[i][p].slice(1)   ) + 
00208                                                     dx_deta*( dphi_dzeta[i][p].slice(0) - 
00209                                                               dphi_dxi[i][p].slice(2)     ) +
00210                                                     dx_dzeta*( dphi_dxi[i][p].slice(1) - 
00211                                                                dphi_deta[i][p].slice(0)    ) );
00212 
00213                 curl_phi[i][p].slice(1) = inv_jac*( dy_dxi*( dphi_deta[i][p].slice(2) - 
00214                                                              dphi_dzeta[i][p].slice(1)  ) +
00215                                                     dy_deta*( dphi_dzeta[i][p].slice(0)- 
00216                                                               dphi_dxi[i][p].slice(2)    ) +
00217                                                     dy_dzeta*( dphi_dxi[i][p].slice(1) - 
00218                                                                dphi_deta[i][p].slice(0)   ) );
00219 
00220                 curl_phi[i][p].slice(2) = inv_jac*( dz_dxi*( dphi_deta[i][p].slice(2) - 
00221                                                              dphi_dzeta[i][p].slice(1)   ) +
00222                                                     dz_deta*( dphi_dzeta[i][p].slice(0) -
00223                                                               dphi_dxi[i][p].slice(2)     ) +
00224                                                     dz_dzeta*( dphi_dxi[i][p].slice(1) - 
00225                                                                dphi_deta[i][p].slice(0)    ) );
00226               }
00227           
00228           break;
00229         }
00230         
00231       default:
00232         libmesh_error();
00233         
00234       } // switch(dim)
00235 
00236     return;
00237   }

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.

00086     { libmesh_do_once( libMesh::err << "WARNING: Shape function Hessians for HCurl elements are not currently "
00087                        << "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.

00112     { libmesh_do_once( libMesh::err << "WARNING: Shape function divergences for HCurl elements are not currently "
00113                        << "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.

00067     { libmesh_do_once( libMesh::err << "WARNING: Shape function gradients for HCurl elements are not currently "
00068                        << "being computed!" << std::endl; ); return; }

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 [inline]

Definition at line 242 of file hcurl_fe_transformation.C.

References libMesh::err.

00247 {
00248   libMesh::err << "HCurl transformations only make sense for vector-valued elements."
00249                << std::endl;
00250   libmesh_error();
00251 }

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 [inline, 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().

00030   {
00031     switch(dim)
00032       {
00033         // These element transformations only make sense in 2D and 3D
00034       case 0: 
00035       case 1:
00036         {
00037           libmesh_error();
00038         }
00039         
00040       case 2:
00041         {
00042           const std::vector<Real>& dxidx_map = fe.get_fe_map().get_dxidx();
00043           const std::vector<Real>& dxidy_map = fe.get_fe_map().get_dxidy();
00044 
00045           const std::vector<Real>& detadx_map = fe.get_fe_map().get_detadx();
00046           const std::vector<Real>& detady_map = fe.get_fe_map().get_detady();
00047 
00048           // FIXME: Need to update for 2D elements in 3D space
00049           /* phi = (dx/dxi)^-T * \hat{phi} 
00050              In 2D:
00051                      (dx/dxi)^{-1} = [  dxi/dx   dxi/dy 
00052                                         deta/dx  deta/dy ]
00053 
00054                  so: dxi/dx^{-T} * \hat{phi} = [ dxi/dx  deta/dx   [ \hat{phi}_xi
00055                                                  dxi/dy  deta/dy ]   \hat{phi}_eta ]
00056 
00057               or in indicial notation:  phi_j = xi_{i,j}*\hat{phi}_i */
00058 
00059           for (unsigned int i=0; i<phi.size(); i++)
00060             for (unsigned int p=0; p<phi[i].size(); p++)
00061               {
00062                 // Need to temporarily cache reference shape functions
00063                 // TODO: PB: Might be worth trying to build phi_ref separately to see
00064                 //           if we can get vectorization
00065                 OutputShape phi_ref;
00066                 FEInterface::shape<OutputShape>(2, fe.get_fe_type(), elem, i, qp[p], phi_ref);
00067 
00068                 phi[i][p](0) = dxidx_map[p]*phi_ref.slice(0) + detadx_map[p]*phi_ref.slice(1);
00069 
00070                 phi[i][p](1) = dxidy_map[p]*phi_ref.slice(0) + detady_map[p]*phi_ref.slice(1);
00071               }
00072 
00073           break;
00074         }
00075 
00076       case 3:
00077         {
00078           const std::vector<Real>& dxidx_map = fe.get_fe_map().get_dxidx();
00079           const std::vector<Real>& dxidy_map = fe.get_fe_map().get_dxidy();
00080           const std::vector<Real>& dxidz_map = fe.get_fe_map().get_dxidz();
00081 
00082           const std::vector<Real>& detadx_map = fe.get_fe_map().get_detadx();
00083           const std::vector<Real>& detady_map = fe.get_fe_map().get_detady();
00084           const std::vector<Real>& detadz_map = fe.get_fe_map().get_detadz();
00085 
00086           const std::vector<Real>& dzetadx_map = fe.get_fe_map().get_dzetadx();
00087           const std::vector<Real>& dzetady_map = fe.get_fe_map().get_dzetady();
00088           const std::vector<Real>& dzetadz_map = fe.get_fe_map().get_dzetadz();
00089 
00090           /* phi = (dx/dxi)^-T * \hat{phi} 
00091              In 3D:
00092                   dx/dxi^-1 = [  dxi/dx    dxi/dy    dxi/dz
00093                                  deta/dx   deta/dy   deta/dz
00094                                  dzeta/dx  dzeta/dy  dzeta/dz]
00095 
00096                  so: dxi/dx^-T * \hat{phi} = [ dxi/dx  deta/dx  dzeta/dx   [ \hat{phi}_xi
00097                                                dxi/dy  deta/dy  dzeta/dy     \hat{phi}_eta 
00098                                                dxi/dz  deta/dz  dzeta/dz ]   \hat{phi}_zeta ]
00099 
00100               or in indicial notation:  phi_j = xi_{i,j}*\hat{phi}_i */
00101 
00102           for (unsigned int i=0; i<phi.size(); i++)
00103             for (unsigned int p=0; p<phi[i].size(); p++)
00104               {
00105                 // Need to temporarily cache reference shape functions
00106                 // TODO: PB: Might be worth trying to build phi_ref separately to see
00107                 //           if we can get vectorization
00108                 OutputShape phi_ref;
00109                 FEInterface::shape<OutputShape>(3, fe.get_fe_type(), elem, i, qp[p], phi_ref);
00110 
00111                 phi[i][p].slice(0) = dxidx_map[p]*phi_ref.slice(0) + detadx_map[p]*phi_ref.slice(1)
00112                   + dzetadx_map[p]*phi_ref.slice(2);
00113 
00114                 phi[i][p].slice(1) = dxidy_map[p]*phi_ref.slice(0) + detady_map[p]*phi_ref.slice(1)
00115                   + dzetady_map[p]*phi_ref.slice(2);
00116 
00117                 phi[i][p].slice(2) = dxidz_map[p]*phi_ref.slice(0) + detadz_map[p]*phi_ref.slice(1)
00118                   + dzetadz_map[p]*phi_ref.slice(2);
00119               }
00120           break;
00121         }
00122         
00123       default:
00124         libmesh_error();
00125 
00126       } // switch(dim)
00127 
00128     return;
00129   }


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

Site Created By: libMesh Developers
Last modified: February 05 2013 19:55:26 UTC

Hosted By:
SourceForge.net Logo