SundanceFeketeQuadrature.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_FEKETEQUADRATURE_H
00032 #define SUNDANCE_FEKETEQUADRATURE_H
00033 
00034 #include "SundanceDefs.hpp"
00035 #include "Teuchos_RefCountPtr.hpp"
00036 #include "SundanceQuadratureFamilyBase.hpp"
00037 
00038 namespace Sundance
00039 {
00040 
00041 using namespace Teuchos;
00042 
00043 /**
00044  * Family of Fekete approach integration rules, suitable for integration
00045  * over arbitrary regions
00046  * (For tensor elements the Gauss-Lobatto points are Fekete points)
00047  */
00048 class FeketeQuadrature: public QuadratureFamilyBase
00049 {
00050 public:
00051   /** */
00052   FeketeQuadrature(int order);
00053 
00054   /** */
00055   virtual ~FeketeQuadrature()
00056   {;}
00057 
00058   /** */
00059   virtual XMLObject toXML() const;
00060 
00061   /** Describable interface */
00062   virtual std::string description() const
00063   {
00064     return "FeketeQuadrature[order=" + Teuchos::toString(order()) + "]";
00065   }
00066 
00067   /* handleable boilerplate */
00068   GET_RCP(QuadratureFamilyStub)
00069   ;
00070 
00071 protected:
00072   /** compute a rule for the reference line cell */
00073   virtual void getLineRule(Array<Point>& quadPoints,
00074       Array<double>& quadWeights) const;
00075 
00076   /** compute a rule for the reference triangle cell */
00077   virtual void getTriangleRule(Array<Point>& quadPoints,
00078       Array<double>& quadWeights) const;
00079 
00080   /** compute a rule for the reference quad cell */
00081   virtual void getQuadRule(Array<Point>& quadPoints,
00082       Array<double>& quadWeights) const;
00083 
00084   /** compute a rule for the reference brick cell */
00085   virtual void getBrickRule(Array<Point>& quadPoints,
00086       Array<double>& quadWeights) const;
00087 
00088   /**
00089    * Compute adapted weights according to curve
00090    *
00091    * @param cellType
00092    * @param cellDim
00093    * @param cellLID
00094    * @param facetIndex
00095    * @param mesh
00096    * @param globalCurve
00097    * @param quadPoints
00098    * @param quadWeights
00099    * @param changedWeights
00100    */
00101   virtual void getAdaptedWeights(const CellType& cellType, int cellDim,
00102       int cellLID, int facetIndex, const Mesh& mesh,
00103       const ParametrizedCurve& globalCurve,
00104       Array<Point>& quadPoints, Array<double>& quadWeights,
00105       bool &weightsChanged) const;
00106 
00107   virtual void getAdaptedTriangleWeights(int cellLID, const Mesh& mesh,
00108       const ParametrizedCurve& globalCurve, Array<Point>& quadPoints,
00109       Array<double>& quadWeights, bool& weightsChanged) const;
00110 
00111   virtual void getAdaptedQuadWeights(int cellLID, const Mesh& mesh,
00112       const ParametrizedCurve& globalCurve, Array<Point>& quadPoints,
00113       Array<double>& quadWeights, bool& weightsChanged) const;
00114 
00115   virtual void integrateRegion(const CellType& cellType, int cellDim,
00116       const int innerOrder, const Point& x, const Point& xref,
00117       const Point& vec1, const Point& vec2, const Point& vec1ref,
00118       const Point& vec2ref, const ParametrizedCurve& curve,
00119       Array<double>& integrals) const;
00120 
00121   virtual void evaluateAllBasisFunctions(const CellType cellType, const Point& q,
00122       Array<double>& result) const;
00123 
00124 private:
00125 
00126   mutable bool _hasBasisCoeffs;
00127   mutable Array<double> _basisCoeffs;
00128 };
00129 }
00130 
00131 #endif

Site Contact