SundanceGaussianQuadrature.hpp
Go to the documentation of this file.
00001 /* @HEADER@ */
00002 // ************************************************************************
00003 // 
00004 //                              Sundance
00005 //                 Copyright (2005) Sandia Corporation
00006 // 
00007 // Copyright (year first published) Sandia Corporation.  Under the terms 
00008 // of Contract DE-AC04-94AL85000 with Sandia Corporation, the U.S. Government 
00009 // retains certain rights in this software.
00010 // 
00011 // This library is free software; you can redistribute it and/or modify
00012 // it under the terms of the GNU Lesser General Public License as
00013 // published by the Free Software Foundation; either version 2.1 of the
00014 // License, or (at your option) any later version.
00015 //  
00016 // This library is distributed in the hope that it will be useful, but
00017 // WITHOUT ANY WARRANTY; without even the implied warranty of
00018 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00019 // Lesser General Public License for more details.
00020 //                                                                                 
00021 // You should have received a copy of the GNU Lesser General Public
00022 // License along with this library; if not, write to the Free Software
00023 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307
00024 // USA                                                                                
00025 // Questions? Contact Kevin Long (krlong@sandia.gov), 
00026 // Sandia National Laboratories, Livermore, California, USA
00027 // 
00028 // ************************************************************************
00029 /* @HEADER@ */
00030 
00031 #ifndef SUNDANCE_GAUSSIANQUADRATURE_H
00032 #define SUNDANCE_GAUSSIANQUADRATURE_H
00033 
00034 #include "SundanceDefs.hpp"
00035 #include "Teuchos_RefCountPtr.hpp"
00036 #include "SundanceQuadratureFamilyBase.hpp"
00037 
00038 namespace Sundance
00039 {
00040 using namespace Teuchos;
00041 
00042 
00043 /** 
00044  * Family of optimal Gaussian integration rules, e.g., Gauss-Legendre on 
00045  * lines, Dunavant on triangles. 
00046  */
00047 class GaussianQuadrature : public QuadratureFamilyBase
00048 {
00049 public:
00050   /** */
00051   GaussianQuadrature(int order);
00052 
00053   /** */
00054   virtual ~GaussianQuadrature(){;}
00055 
00056 
00057   /** */
00058   virtual XMLObject toXML() const ;
00059 
00060   /** Describable interface */
00061   virtual std::string description() const 
00062     {return "GaussianQuadrature[order=" + Teuchos::toString(order()) 
00063         +  "]";}
00064 
00065   /* handleable boilerplate */
00066   GET_RCP(QuadratureFamilyStub);
00067 
00068   /** This methos is for the ACI integration */
00069   virtual void getAdaptedWeights(const CellType& cellType ,
00070                      int cellDim,
00071                                      int celLID ,
00072                                    int facetIndex ,
00073                                      const Mesh& mesh ,
00074                                      const ParametrizedCurve& globalCurve ,
00075                                      Array<Point>& quadPoints ,
00076                                      Array<double>& quadWeights ,
00077                                      bool &isCut) const ;
00078 
00079 protected:
00080   /** compute a rule for the reference line cell */
00081   virtual void getLineRule(Array<Point>& quadPoints,
00082     Array<double>& quadWeights) const ;
00083 
00084   /** compute a rule for the reference triangle cell */
00085   virtual void getTriangleRule(Array<Point>& quadPoints,
00086     Array<double>& quadWeights) const ;
00087 
00088   /** compute a rule for the reference quad cell */
00089   virtual void getQuadRule(Array<Point>& quadPoints,
00090     Array<double>& quadWeights) const ;
00091 
00092   /** compute a rule for the reference tet cell */
00093   virtual void getTetRule(Array<Point>& quadPoints,
00094     Array<double>& quadWeights) const ;
00095 
00096   /** compute a rule for the reference brick cell */
00097   virtual void getBrickRule(Array<Point>& quadPoints,
00098     Array<double>& quadWeights) const ;
00099 
00100 };
00101 }
00102 
00103 
00104 #endif

Site Contact