SundanceGaussianQuadrature.hpp
Go to the documentation of this file.
00001 /* @HEADER@ */
00002 // ************************************************************************
00003 // 
00004 //                             Sundance
00005 //                 Copyright 2011 Sandia Corporation
00006 // 
00007 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
00008 // the U.S. Government retains certain rights in this software.
00009 //
00010 // Redistribution and use in source and binary forms, with or without
00011 // modification, are permitted provided that the following conditions are
00012 // met:
00013 //
00014 // 1. Redistributions of source code must retain the above copyright
00015 // notice, this list of conditions and the following disclaimer.
00016 //
00017 // 2. Redistributions in binary form must reproduce the above copyright
00018 // notice, this list of conditions and the following disclaimer in the
00019 // documentation and/or other materials provided with the distribution.
00020 //
00021 // 3. Neither the name of the Corporation nor the names of the
00022 // contributors may be used to endorse or promote products derived from
00023 // this software without specific prior written permission.
00024 //
00025 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
00026 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
00027 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
00028 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
00029 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
00030 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
00031 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
00032 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
00033 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
00034 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
00035 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
00036 //
00037 // Questions? Contact Kevin Long (kevin.long@ttu.edu)
00038 // 
00039 
00040 /* @HEADER@ */
00041 
00042 #ifndef SUNDANCE_GAUSSIANQUADRATURE_H
00043 #define SUNDANCE_GAUSSIANQUADRATURE_H
00044 
00045 #include "SundanceDefs.hpp"
00046 #include "Teuchos_RefCountPtr.hpp"
00047 #include "SundanceQuadratureFamilyBase.hpp"
00048 
00049 namespace Sundance
00050 {
00051 using namespace Teuchos;
00052 
00053 
00054 /** 
00055  * Family of optimal Gaussian integration rules, e.g., Gauss-Legendre on 
00056  * lines, Dunavant on triangles. 
00057  */
00058 class GaussianQuadrature : public QuadratureFamilyBase
00059 {
00060 public:
00061   /** */
00062   GaussianQuadrature(int order);
00063 
00064   /** */
00065   virtual ~GaussianQuadrature(){;}
00066 
00067 
00068   /** */
00069   virtual XMLObject toXML() const ;
00070 
00071   /** Describable interface */
00072   virtual std::string description() const 
00073     {return "GaussianQuadrature[order=" + Teuchos::toString(order()) 
00074         +  "]";}
00075 
00076   /* handleable boilerplate */
00077   GET_RCP(QuadratureFamilyStub);
00078 
00079   /** This methos is for the ACI integration */
00080   virtual void getAdaptedWeights(const CellType& cellType ,
00081                      int cellDim,
00082                                      int celLID ,
00083                                    int facetIndex ,
00084                                      const Mesh& mesh ,
00085                                      const ParametrizedCurve& globalCurve ,
00086                                      Array<Point>& quadPoints ,
00087                                      Array<double>& quadWeights ,
00088                                      bool &isCut) const ;
00089 
00090 protected:
00091   /** compute a rule for the reference line cell */
00092   virtual void getLineRule(Array<Point>& quadPoints,
00093     Array<double>& quadWeights) const ;
00094 
00095   /** compute a rule for the reference triangle cell */
00096   virtual void getTriangleRule(Array<Point>& quadPoints,
00097     Array<double>& quadWeights) const ;
00098 
00099   /** compute a rule for the reference quad cell */
00100   virtual void getQuadRule(Array<Point>& quadPoints,
00101     Array<double>& quadWeights) const ;
00102 
00103   /** compute a rule for the reference tet cell */
00104   virtual void getTetRule(Array<Point>& quadPoints,
00105     Array<double>& quadWeights) const ;
00106 
00107   /** compute a rule for the reference brick cell */
00108   virtual void getBrickRule(Array<Point>& quadPoints,
00109     Array<double>& quadWeights) const ;
00110 
00111 };
00112 }
00113 
00114 
00115 #endif

Site Contact