Main Page   Namespace List   Class Hierarchy   Compound List   File List   Namespace Members   Compound Members   File Members   Related Pages   Examples  

GaussLegendre.h

Go to the documentation of this file.
00001 #ifndef GAUSSLEGENDRE_H
00002 #define GAUSSLEGENDRE_H
00003 
00004 #include "SundanceDefs.h"
00005 #include "QuadratureFamilyBase.h"
00006 #include "TSFHashtable.h"
00007 #include "XMLObject.h"
00008 
00009 
00010 namespace Sundance
00011 {
00012 
00013   using namespace TSF;
00014   using std::string;
00015 
00016   using std::ostream;
00017 
00018   /** \ingroup UserLevelFE
00019    * Abstract, geometry-independent Gauss-Legendre quadrature rule
00020    */
00021 
00022   class GaussLegendre : public QuadratureFamilyBase
00023     {
00024     public:
00025       /** \name User-level methods */
00026       //@{
00027       /** create a Gauss-Legendre quadrature rule of a given
00028        * polynomial order */
00029       GaussLegendre(int polyOrder) : QuadratureFamilyBase("GaussLegendre",
00030                                                           polyOrder) {;}
00031       //@}
00032 
00033       /** \name Developer-only methods */
00034       //@{
00035       /** create from an XML specification */
00036       GaussLegendre(const XMLObject& xml);
00037       /** TUVD */
00038       virtual ~GaussLegendre(){;}
00039 
00040       /** return an identifying string */
00041       virtual const string& typeName() const {return typeName_;}
00042 
00043       //@}
00044     protected:
00045 
00046       /** compute a rule for the reference line cell */
00047       virtual QuadraturePoints getLineRule() const ;
00048 
00049       /** compute a rule for the reference triangle cell */
00050       virtual QuadraturePoints getTriangleRule() const ;
00051 
00052       /** compute a rule for the reference quad cell */
00053       virtual QuadraturePoints getQuadRule() const ;
00054 
00055       /** compute a rule for the reference brick cell */
00056       virtual QuadraturePoints getBrickRule() const ;
00057 
00058       /** compute a rule for the reference tet cell */
00059       virtual QuadraturePoints getTetRule() const ;
00060 
00061       static string typeName_;
00062       static TSFHashtable<int, QuadraturePoints> lineRules_;
00063       static TSFHashtable<int, QuadraturePoints> triangleRules_;
00064       static TSFHashtable<int, QuadraturePoints> quadRules_;
00065       static TSFHashtable<int, QuadraturePoints> brickRules_;
00066       static TSFHashtable<int, QuadraturePoints> tetRules_;
00067     };
00068 
00069   extern "C"
00070   {
00071     /** \relates GaussLegendre creator function used in dynamic loading */
00072     QuadratureFamilyBase* createXMLGaussLegendre(const XMLObject& xml);
00073   }
00074 
00075 }
00076 #endif

Contact:
Kevin Long (krlong@ca.sandia.gov)


Documentation generated by