SundanceQuadratureFamilyBase.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_QUADRATUREFAMILYBASE_H
00043 #define SUNDANCE_QUADRATUREFAMILYBASE_H
00044 
00045 #include "SundanceDefs.hpp"
00046 #include "Teuchos_RefCountPtr.hpp"
00047 #include "Teuchos_Array.hpp"
00048 #include "SundanceQuadratureFamilyStub.hpp"
00049 #include "SundanceCellType.hpp"
00050 #include "SundancePoint.hpp"
00051 #include "SundanceParametrizedCurve.hpp"
00052 #include "SundanceMesh.hpp"
00053 
00054 namespace Sundance
00055 {
00056 using namespace Teuchos;
00057 
00058 /** 
00059  * QuadratureFamilyBase extends QuadratureFamilyStub to provide
00060  * an interface for getting quadrature points for a given cell type. 
00061  */
00062 class QuadratureFamilyBase : public QuadratureFamilyStub 
00063 {
00064 public:
00065   /** */
00066   QuadratureFamilyBase(int order) : QuadratureFamilyStub(order) {;}
00067 
00068   /** */
00069   virtual ~QuadratureFamilyBase(){;}
00070 
00071   /** Gets number of points associated with a particular cell type:
00072       WARNING: this is slow.  Call it once and store the result. 
00073       TODO: make it pure virtual and override with queries in
00074       the derived classes, making them supply the information.  */
00075   virtual int getNumPoints( const CellType &cellType ) const 
00076     {
00077       Array<Point> qp;
00078       Array<double> qw;
00079       this->getPoints(cellType,qp,qw);
00080       return qp.size();
00081     }
00082 
00083   /** Get the quadrature points and weights for the given cell type */
00084   virtual void getPoints(const CellType& cellType, 
00085     Array<Point>& quadPoints,
00086     Array<double>& quadWeights) const ;
00087       
00088   /** This method is used for the Adaptive Cell Integration, which returns
00089    * special quadrature weights for cells which are cut by the curve <br>
00090    * The quadPoints and the quadWeights are set to the default values.
00091    * (quadrature without the curve , no curve)
00092    * quadWeights will be changed depending on the curve, if the curve cuts the cell
00093    * isCut should be set by this method, if it is true then the quadWeights will be
00094    * used for quadrature of the cell */
00095   virtual void getAdaptedWeights(const CellType& cellType ,
00096   int cellDim,
00097   int celLID ,
00098     int facetIndex ,
00099     const Mesh& mesh ,
00100     const ParametrizedCurve& globalCurve ,
00101     Array<Point>& quadPoints ,
00102     Array<double>& quadWeights,
00103     bool& isCut) const ;
00104 
00105 protected:
00106 
00107   /** compute a rule for the reference line cell */
00108   virtual void getLineRule(Array<Point>& quadPoints,
00109     Array<double>& quadWeights) const ;
00110 
00111   /** compute a rule for the reference triangle cell */
00112   virtual void getTriangleRule(Array<Point>& quadPoints,
00113     Array<double>& quadWeights) const ;
00114 
00115   /** compute a rule for the reference quad cell */
00116   virtual void getQuadRule(Array<Point>& quadPoints,
00117     Array<double>& quadWeights) const ;
00118 
00119   /** compute a rule for the reference tet cell */
00120   virtual void getTetRule(Array<Point>& quadPoints,
00121     Array<double>& quadWeights) const ;
00122 
00123 
00124   /** compute a rule for the reference brick cell */
00125   virtual void getBrickRule(Array<Point>& quadPoints,
00126     Array<double>& quadWeights) const ;
00127 private:
00128 };
00129 
00130 }
00131 
00132 #endif
00133 

Site Contact