$extrastylesheet
hcurl_fe_transformation.h
Go to the documentation of this file.
00001 // The libMesh Finite Element Library.
00002 // Copyright (C) 2002-2014 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
00003 
00004 // This library is free software; you can redistribute it and/or
00005 // modify it under the terms of the GNU Lesser General Public
00006 // License as published by the Free Software Foundation; either
00007 // version 2.1 of the License, or (at your option) any later version.
00008 
00009 // This library is distributed in the hope that it will be useful,
00010 // but WITHOUT ANY WARRANTY; without even the implied warranty of
00011 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00012 // Lesser General Public License for more details.
00013 
00014 // You should have received a copy of the GNU Lesser General Public
00015 // License along with this library; if not, write to the Free Software
00016 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00017 
00018 #ifndef LIBMESH_HCURL_FE_TRANSFORMATION_H
00019 #define LIBMESH_HCURL_FE_TRANSFORMATION_H
00020 
00021 #include "libmesh/fe_transformation_base.h"
00022 
00023 namespace libMesh
00024 {
00025 
00033 template< typename OutputShape >
00034 class HCurlFETransformation : public FETransformationBase<OutputShape>
00035 {
00036 public:
00037 
00038   HCurlFETransformation()
00039     : FETransformationBase<OutputShape>(){}
00040 
00041   virtual ~HCurlFETransformation(){}
00042 
00049   virtual void map_phi( const unsigned int dim,
00050                         const Elem* const elem,
00051                         const std::vector<Point>& qp,
00052                         const FEGenericBase<OutputShape>& fe,
00053                         std::vector<std::vector<OutputShape> >& phi ) const;
00054 
00059   virtual void map_dphi( const unsigned int /*dim*/,
00060                          const Elem* const /*elem*/,
00061                          const std::vector<Point>& /*qp*/,
00062                          const FEGenericBase<OutputShape>& /*fe*/,
00063                          std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputGradient> >& /*dphi*/,
00064                          std::vector<std::vector<OutputShape> >& /*dphidx*/,
00065                          std::vector<std::vector<OutputShape> >& /*dphidy*/,
00066                          std::vector<std::vector<OutputShape> >& /*dphidz*/) const
00067   { libmesh_do_once( libMesh::err << "WARNING: Shape function gradients for HCurl elements are not currently "
00068                      << "being computed!" << std::endl; ); return; }
00069 
00070 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
00071 
00075   virtual void map_d2phi( const unsigned int /*dim*/,
00076                           const Elem* const /*elem*/,
00077                           const std::vector<Point>& /*qp*/,
00078                           const FEGenericBase<OutputShape>& /*fe*/,
00079                           std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputTensor> >& /*d2phi*/,
00080                           std::vector<std::vector<OutputShape> >& /*d2phidx2*/,
00081                           std::vector<std::vector<OutputShape> >& /*d2phidxdy*/,
00082                           std::vector<std::vector<OutputShape> >& /*d2phidxdz*/,
00083                           std::vector<std::vector<OutputShape> >& /*d2phidy2*/,
00084                           std::vector<std::vector<OutputShape> >& /*d2phidydz*/,
00085                           std::vector<std::vector<OutputShape> >& /*d2phidz2*/  ) const
00086   { libmesh_do_once( libMesh::err << "WARNING: Shape function Hessians for HCurl elements are not currently "
00087                      << "being computed!" << std::endl; ); return; }
00088 #endif //LIBMESH_ENABLE_SECOND_DERIVATIVES
00089 
00097   virtual void map_curl( const unsigned int dim,
00098                          const Elem* const elem,
00099                          const std::vector<Point>& qp,
00100                          const FEGenericBase<OutputShape>& fe,
00101                          std::vector<std::vector<OutputShape> >& curl_phi ) const;
00102 
00107   virtual void map_div( const unsigned int /*dim*/,
00108                         const Elem* const /*elem*/,
00109                         const std::vector<Point>& /*qp*/,
00110                         const FEGenericBase<OutputShape>& /*fe*/,
00111                         std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputDivergence> >& /*div_phi*/ ) const
00112   { libmesh_do_once( libMesh::err << "WARNING: Shape function divergences for HCurl elements are not currently "
00113                      << "being computed!" << std::endl; ); return; }
00114 
00115 }; // class HCurlFETransformation
00116 
00117 }
00118 
00119 #endif // LIBMESH_HCURL_FE_TRANSFORMATION_H