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

Site Contact