fe.h
Go to the documentation of this file.00001 // The libMesh Finite Element Library. 00002 // Copyright (C) 2002-2012 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 00019 00020 #ifndef LIBMESH_FE_H 00021 #define LIBMESH_FE_H 00022 00023 // Local includes 00024 #include "libmesh/fe_base.h" 00025 #include "libmesh/libmesh.h" 00026 00027 // C++ includes 00028 #include <cstddef> 00029 00030 namespace libMesh 00031 { 00032 00033 // forward declarations 00034 class DofConstraints; 00035 class DofMap; 00036 00037 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 00038 00039 template <unsigned int friend_Dim, FEFamily friend_T_radial, InfMapType friend_T_map> 00040 class InfFE; 00041 00042 #endif 00043 00044 00048 template <FEFamily T> 00049 struct FEOutputType 00050 { 00051 typedef Real type; 00052 }; 00053 00054 00058 template<> 00059 struct FEOutputType<LAGRANGE_VEC> 00060 { 00061 typedef RealGradient type; 00062 }; 00063 00064 template<> 00065 struct FEOutputType<NEDELEC_ONE> 00066 { 00067 typedef RealGradient type; 00068 }; 00069 00070 00088 //------------------------------------------------------------- 00089 // FE class definition 00090 template <unsigned int Dim, FEFamily T> 00091 class FE : public FEGenericBase<typename FEOutputType<T>::type> 00092 { 00093 public: 00094 00098 explicit 00099 FE(const FEType& fet); 00100 00101 typedef typename 00102 FEGenericBase<typename FEOutputType<T>::type>::OutputShape 00103 OutputShape; 00104 00113 static OutputShape shape(const ElemType t, 00114 const Order o, 00115 const unsigned int i, 00116 const Point& p); 00117 00126 static OutputShape shape(const Elem* elem, 00127 const Order o, 00128 const unsigned int i, 00129 const Point& p); 00130 00138 static OutputShape shape_deriv(const ElemType t, 00139 const Order o, 00140 const unsigned int i, 00141 const unsigned int j, 00142 const Point& p); 00143 00150 static OutputShape shape_deriv(const Elem* elem, 00151 const Order o, 00152 const unsigned int i, 00153 const unsigned int j, 00154 const Point& p); 00155 00174 static OutputShape shape_second_deriv(const ElemType t, 00175 const Order o, 00176 const unsigned int i, 00177 const unsigned int j, 00178 const Point& p); 00179 00198 static OutputShape shape_second_deriv(const Elem* elem, 00199 const Order o, 00200 const unsigned int i, 00201 const unsigned int j, 00202 const Point& p); 00203 00210 static void nodal_soln(const Elem* elem, const Order o, 00211 const std::vector<Number>& elem_soln, 00212 std::vector<Number>& nodal_soln); 00213 00218 virtual unsigned int n_shape_functions () const; 00219 00226 static unsigned int n_shape_functions (const ElemType t, 00227 const Order o) 00228 { return FE<Dim,T>::n_dofs (t,o); } 00229 00236 static unsigned int n_dofs(const ElemType t, 00237 const Order o); 00238 00245 static unsigned int n_dofs_at_node(const ElemType t, 00246 const Order o, 00247 const unsigned int n); 00248 00255 static unsigned int n_dofs_per_elem(const ElemType t, 00256 const Order o); 00257 00261 virtual FEContinuity get_continuity() const; 00262 00267 virtual bool is_hierarchic() const; 00268 00275 static void dofs_on_side(const Elem* const elem, 00276 const Order o, 00277 unsigned int s, 00278 std::vector<unsigned int>& di); 00285 static void dofs_on_edge(const Elem* const elem, 00286 const Order o, 00287 unsigned int e, 00288 std::vector<unsigned int>& di); 00289 00299 static Point inverse_map (const Elem* elem, 00300 const Point& p, 00301 const Real tolerance = TOLERANCE, 00302 const bool secure = true); 00303 00314 static void inverse_map (const Elem* elem, 00315 const std::vector<Point>& physical_points, 00316 std::vector<Point>& reference_points, 00317 const Real tolerance = TOLERANCE, 00318 const bool secure = true); 00319 00330 virtual void reinit (const Elem* elem, 00331 const std::vector<Point>* const pts = NULL, 00332 const std::vector<Real>* const weights = NULL); 00333 00343 virtual void reinit (const Elem* elem, 00344 const unsigned int side, 00345 const Real tolerance = TOLERANCE, 00346 const std::vector<Point>* const pts = NULL, 00347 const std::vector<Real>* const weights = NULL); 00348 00358 virtual void edge_reinit (const Elem* elem, 00359 const unsigned int edge, 00360 const Real tolerance = TOLERANCE, 00361 const std::vector<Point>* const pts = NULL, 00362 const std::vector<Real>* const weights = NULL); 00363 00368 virtual void side_map (const Elem* elem, 00369 const Elem* side, 00370 const unsigned int s, 00371 const std::vector<Point>& reference_side_points, 00372 std::vector<Point>& reference_points); 00373 00379 virtual void attach_quadrature_rule (QBase* q); 00380 00386 virtual unsigned int n_quadrature_points () const; 00387 00388 #ifdef LIBMESH_ENABLE_AMR 00389 00395 static void compute_constraints (DofConstraints &constraints, 00396 DofMap &dof_map, 00397 const unsigned int variable_number, 00398 const Elem* elem); 00399 #endif // #ifdef LIBMESH_ENABLE_AMR 00400 00407 virtual bool shapes_need_reinit() const; 00408 00413 static Point map (const Elem* elem, 00414 const Point& reference_point); 00415 00420 static Point map_xi (const Elem* elem, 00421 const Point& reference_point); 00422 00427 static Point map_eta (const Elem* elem, 00428 const Point& reference_point); 00429 00434 static Point map_zeta (const Elem* elem, 00435 const Point& reference_point); 00436 00437 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 00438 00442 template <unsigned int friend_Dim, FEFamily friend_T_radial, InfMapType friend_T_map> 00443 friend class InfFE; 00444 #endif 00445 00446 protected: 00447 00455 virtual void init_shape_functions(const std::vector<Point>& qp, 00456 const Elem* e); 00457 00458 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 00459 00464 virtual void init_base_shape_functions(const std::vector<Point>& qp, 00465 const Elem* e); 00466 00467 #endif 00468 00473 std::vector<Point> cached_nodes; 00474 00478 ElemType last_side; 00479 00480 unsigned int last_edge; 00481 }; 00482 00483 00484 00493 //------------------------------------------------------------- 00494 // FEHierarchic class definition 00495 template <unsigned int Dim> 00496 class FEClough : public FE<Dim,CLOUGH> 00497 { 00498 public: 00499 00504 explicit 00505 FEClough(const FEType& fet); 00506 }; 00507 00508 00509 00518 //------------------------------------------------------------- 00519 // FEHierarchic class definition 00520 template <unsigned int Dim> 00521 class FEHermite : public FE<Dim,HERMITE> 00522 { 00523 public: 00524 00529 explicit 00530 FEHermite(const FEType& fet); 00531 00535 static Real hermite_raw_shape_second_deriv(const unsigned int basis_num, 00536 const Real xi); 00537 static Real hermite_raw_shape_deriv(const unsigned int basis_num, 00538 const Real xi); 00539 static Real hermite_raw_shape(const unsigned int basis_num, 00540 const Real xi); 00541 }; 00542 00543 00544 00553 //------------------------------------------------------------- 00554 // FEHierarchic class definition 00555 template <unsigned int Dim> 00556 class FEHierarchic : public FE<Dim,HIERARCHIC> 00557 { 00558 public: 00559 00564 explicit 00565 FEHierarchic(const FEType& fet); 00566 }; 00567 00568 00569 00578 //------------------------------------------------------------- 00579 // FEL2Hierarchic class definition 00580 template <unsigned int Dim> 00581 class FEL2Hierarchic : public FE<Dim,L2_HIERARCHIC> 00582 { 00583 public: 00584 00589 explicit 00590 FEL2Hierarchic(const FEType& fet); 00591 }; 00592 00593 00594 00603 //------------------------------------------------------------- 00604 // FELagrange class definition 00605 template <unsigned int Dim> 00606 class FELagrange : public FE<Dim,LAGRANGE> 00607 { 00608 public: 00609 00614 explicit 00615 FELagrange(const FEType& fet); 00616 }; 00617 00618 00622 //------------------------------------------------------------- 00623 // FEL2Lagrange class definition 00624 template <unsigned int Dim> 00625 class FEL2Lagrange : public FE<Dim,L2_LAGRANGE> 00626 { 00627 public: 00628 00633 explicit 00634 FEL2Lagrange(const FEType& fet); 00635 }; 00636 00637 00646 //------------------------------------------------------------- 00647 // FEMonomial class definition 00648 template <unsigned int Dim> 00649 class FEMonomial : public FE<Dim,MONOMIAL> 00650 { 00651 public: 00652 00657 explicit 00658 FEMonomial(const FEType& fet); 00659 }; 00660 00661 00662 //------------------------------------------------------------- 00663 // FEScalar class definition 00664 template <unsigned int Dim> 00665 class FEScalar : public FE<Dim,SCALAR> 00666 { 00667 public: 00668 00675 explicit 00676 FEScalar(const FEType& fet); 00677 }; 00678 00679 00689 //------------------------------------------------------------- 00690 // FEXYZ class definition 00691 template <unsigned int Dim> 00692 class FEXYZ : public FE<Dim,XYZ> 00693 { 00694 public: 00695 00700 explicit 00701 FEXYZ(const FEType& fet); 00702 00707 virtual void reinit (const Elem* elem, 00708 const std::vector<Point>* const pts = NULL, 00709 const std::vector<Real>* const weights = NULL) 00710 { FE<Dim,XYZ>::reinit (elem, pts, weights); } 00711 00716 virtual void reinit (const Elem* elem, 00717 const unsigned int side, 00718 const Real tolerance = TOLERANCE, 00719 const std::vector<Point>* const pts = NULL, 00720 const std::vector<Real>* const weights = NULL); 00721 00722 00723 protected: 00724 00732 virtual void init_shape_functions(const std::vector<Point>& qp, 00733 const Elem* e); 00734 00745 virtual void compute_shape_functions(const Elem* elem, const std::vector<Point>& qp); 00746 00750 void compute_face_values (const Elem* elem, 00751 const Elem* side, 00752 const std::vector<Real>& weights); 00753 }; 00754 00755 //------------------------------------------------------------- 00756 // FELagrangeVec class definition 00757 template <unsigned int Dim> 00758 class FELagrangeVec : public FE<Dim,LAGRANGE_VEC> 00759 { 00760 public: 00761 00766 explicit 00767 FELagrangeVec(const FEType& fet); 00768 00769 }; 00770 00771 //------------------------------------------------------------- 00772 // FENedelecOne class definition 00773 template <unsigned int Dim> 00774 class FENedelecOne : public FE<Dim,NEDELEC_ONE> 00775 { 00776 public: 00777 00782 explicit 00783 FENedelecOne(const FEType& fet); 00784 00785 }; 00786 00787 00788 00789 00793 namespace FiniteElements 00794 { 00799 typedef FEClough<2> FEClough2D; 00800 00805 typedef FE<1,HIERARCHIC> FEHierarchic1D; 00806 00811 typedef FE<2,HIERARCHIC> FEHierarchic2D; 00812 00817 typedef FE<3,HIERARCHIC> FEHierarchic3D; 00818 00819 00824 typedef FE<1,L2_HIERARCHIC> FEL2Hierarchic1D; 00825 00830 typedef FE<2,L2_HIERARCHIC> FEL2Hierarchic2D; 00831 00836 typedef FE<3,L2_HIERARCHIC> FEL2Hierarchic3D; 00837 00838 00843 typedef FE<1,LAGRANGE> FELagrange1D; 00844 00849 typedef FE<2,LAGRANGE> FELagrange2D; 00850 00855 typedef FE<3,LAGRANGE> FELagrange3D; 00856 00857 00862 typedef FE<1,L2_LAGRANGE> FEL2Lagrange1D; 00863 00868 typedef FE<2,L2_LAGRANGE> FEL2Lagrange2D; 00869 00874 typedef FE<3,L2_LAGRANGE> FEL2Lagrange3D; 00875 00876 00881 typedef FE<1,MONOMIAL> FEMonomial1D; 00882 00887 typedef FE<2,MONOMIAL> FEMonomial2D; 00888 00893 typedef FE<3,MONOMIAL> FEMonomial3D; 00894 00895 } 00896 00897 00898 00899 00900 // ------------------------------------------------------------ 00901 // FE class inline members 00902 template <unsigned int Dim, FEFamily T> 00903 inline 00904 FE<Dim,T>::FE (const FEType& fet) : 00905 FEGenericBase<typename FEOutputType<T>::type> (Dim,fet), 00906 last_side(INVALID_ELEM), 00907 last_edge(libMesh::invalid_uint) 00908 { 00909 // Sanity check. Make sure the 00910 // Family specified in the template instantiation 00911 // matches the one in the FEType object 00912 libmesh_assert_equal_to (T, this->get_family()); 00913 } 00914 00915 00916 00917 // ------------------------------------------------------------ 00918 // FEClough class inline members 00919 template <unsigned int Dim> 00920 inline 00921 FEClough<Dim>::FEClough (const FEType& fet) : 00922 FE<Dim,CLOUGH> (fet) 00923 { 00924 } 00925 00926 00927 00928 // ------------------------------------------------------------ 00929 // FEHermite class inline members 00930 template <unsigned int Dim> 00931 inline 00932 FEHermite<Dim>::FEHermite (const FEType& fet) : 00933 FE<Dim,HERMITE> (fet) 00934 { 00935 } 00936 00937 00938 00939 // ------------------------------------------------------------ 00940 // FEHierarchic class inline members 00941 template <unsigned int Dim> 00942 inline 00943 FEHierarchic<Dim>::FEHierarchic (const FEType& fet) : 00944 FE<Dim,HIERARCHIC> (fet) 00945 { 00946 } 00947 00948 00949 00950 // ------------------------------------------------------------ 00951 // FEL2Hierarchic class inline members 00952 template <unsigned int Dim> 00953 inline 00954 FEL2Hierarchic<Dim>::FEL2Hierarchic (const FEType& fet) : 00955 FE<Dim,L2_HIERARCHIC> (fet) 00956 { 00957 } 00958 00959 00960 00961 // ------------------------------------------------------------ 00962 // FELagrange class inline members 00963 template <unsigned int Dim> 00964 inline 00965 FELagrange<Dim>::FELagrange (const FEType& fet) : 00966 FE<Dim,LAGRANGE> (fet) 00967 { 00968 } 00969 00970 // ------------------------------------------------------------ 00971 // FELagrangeVec class inline members 00972 template <unsigned int Dim> 00973 inline 00974 FELagrangeVec<Dim>::FELagrangeVec (const FEType& fet) : 00975 FE<Dim,LAGRANGE_VEC> (fet) 00976 { 00977 } 00978 00979 // ------------------------------------------------------------ 00980 // FEL2Lagrange class inline members 00981 template <unsigned int Dim> 00982 inline 00983 FEL2Lagrange<Dim>::FEL2Lagrange (const FEType& fet) : 00984 FE<Dim,L2_LAGRANGE> (fet) 00985 { 00986 } 00987 00988 00989 00990 // ------------------------------------------------------------ 00991 // FEMonomial class inline members 00992 template <unsigned int Dim> 00993 inline 00994 FEMonomial<Dim>::FEMonomial (const FEType& fet) : 00995 FE<Dim,MONOMIAL> (fet) 00996 { 00997 } 00998 00999 01000 01001 01002 // ------------------------------------------------------------ 01003 // FEXYZ class inline members 01004 template <unsigned int Dim> 01005 inline 01006 FEXYZ<Dim>::FEXYZ (const FEType& fet) : 01007 FE<Dim,XYZ> (fet) 01008 { 01009 } 01010 01011 // ------------------------------------------------------------ 01012 // FEScalar class inline members 01013 template <unsigned int Dim> 01014 inline 01015 FEScalar<Dim>::FEScalar (const FEType& fet) : 01016 FE<Dim,SCALAR> (fet) 01017 { 01018 } 01019 01020 // ------------------------------------------------------------ 01021 // FENedelecOne class inline members 01022 template <unsigned int Dim> 01023 inline 01024 FENedelecOne<Dim>::FENedelecOne (const FEType& fet) : 01025 FE<Dim,NEDELEC_ONE> (fet) 01026 { 01027 } 01028 01029 } // namespace libMesh 01030 01031 #endif // LIBMESH_FE_H
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:46 UTC
Hosted By: