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

Site Contact