$extrastylesheet
fe_type.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 
00019 
00020 #ifndef LIBMESH_FE_TYPE_H
00021 #define LIBMESH_FE_TYPE_H
00022 
00023 // Local includes
00024 #include "libmesh/auto_ptr.h"
00025 #include "libmesh/libmesh_config.h"
00026 #include "libmesh/enum_order.h"
00027 #include "libmesh/enum_fe_family.h"
00028 #include "libmesh/enum_inf_map_type.h"
00029 
00030 // C++ includes
00031 
00032 namespace libMesh
00033 {
00034 
00035 // Forward declarations
00036 class QBase;
00037 
00038 
00043 class FEType
00044 {
00045 public:
00046 
00047 #ifndef LIBMESH_ENABLE_INFINITE_ELEMENTS
00048 
00053   FEType(const Order    o = FIRST,
00054          const FEFamily f = LAGRANGE) :
00055     order(o),
00056     family(f)
00057   {}
00058 
00059 
00060   //TODO:[BSK] Could these data types all be const?
00061   // [RHS] Order can't in the case of p refinement!
00065   Order order;
00066 
00071   FEFamily family;
00072 
00073 #else
00074 
00083   FEType(const Order      o  = FIRST,
00084          const FEFamily   f  = LAGRANGE,
00085          const Order      ro = THIRD,
00086          const FEFamily   rf = JACOBI_20_00,
00087          const InfMapType im = CARTESIAN) :
00088     order(o),
00089     radial_order(ro),
00090     family(f),
00091     radial_family(rf),
00092     inf_map(im)
00093   {}
00094 
00095 
00099   Order order;
00100 
00104   Order radial_order;
00105 
00110   FEFamily family;
00111 
00117   FEFamily radial_family;
00118 
00125   InfMapType inf_map;
00126 
00127 #endif // ifndef LIBMESH_ENABLE_INFINITE_ELEMENTS
00128 
00132   bool operator== (const FEType &f2) const
00133   {
00134     return (order == f2.order
00135             && family == f2.family
00136 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
00137             && radial_order == f2.radial_order
00138             && radial_family == f2.radial_family
00139             && inf_map == f2.inf_map
00140 #endif // ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
00141             );
00142   }
00143 
00147   bool operator!= (const FEType &f2) const
00148   {
00149     return !(*this == f2);
00150   }
00151 
00155   bool operator< (const FEType &f2) const
00156   {
00157     if (order != f2.order)
00158       return (order < f2.order);
00159     if (family != f2.family)
00160       return (family < f2.family);
00161 
00162 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
00163     if (radial_order != f2.radial_order)
00164       return (radial_order < f2.radial_order);
00165     if (radial_family != f2.radial_family)
00166       return (radial_family < f2.radial_family);
00167     if (inf_map != f2.inf_map)
00168       return (inf_map < f2.inf_map);
00169 #endif // ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
00170     return false;
00171   }
00172 
00179   Order default_quadrature_order () const;
00180 
00187   UniquePtr<QBase> default_quadrature_rule (const unsigned int dim,
00188                                             const int extraorder=0) const;
00189 
00190 
00191 private:
00192 
00193 };
00194 
00195 
00196 
00197 //-------------------------------------------------------------------
00198 // FEType inline methods
00199 inline
00200 Order FEType::default_quadrature_order () const
00201 {
00202   return static_cast<Order>(2*static_cast<unsigned int>(order) + 1);
00203 }
00204 
00205 } // namespace libMesh
00206 
00207 
00208 #endif // LIBMESH_FE_TYPE_H