SundanceQuadratureFamily.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_QUADRATUREFAMILY_H
00032 #define SUNDANCE_QUADRATUREFAMILY_H
00033 
00034 #include "SundanceDefs.hpp"
00035 #include "SundanceQuadratureFamilyBase.hpp"
00036 #include "PlayaHandle.hpp"
00037 
00038 namespace Sundance
00039 {
00040 
00041 /** 
00042  * QuadratureFamily is a geometry-independent specification of
00043  * a method by which quadrature is to be carried out. For example,
00044  * a GaussianQuadrature family will generate Gaussian
00045  * quadrature points on any cell type.
00046  */
00047 class QuadratureFamily : public Playa::Handle<QuadratureFamilyStub>
00048 {
00049 public:
00050   /* */
00051   HANDLE_CTORS(QuadratureFamily, QuadratureFamilyStub);
00052   /** */
00053   XMLObject toXML() const ;
00054 
00055   /** */
00056   int order() const ;
00057 
00058   /** Returns the number of points in a rule of the given cell type 
00059       WARNING: this is slow.  Call it once and store the result. 
00060       TODO: make it pure virtual and override with queries in
00061       the derived classes, making them supply the information.  */
00062   int getNumPoints( const CellType& cellType ) const;
00063 
00064 
00065   /** Get the quadrature points and weights for the given cell type */
00066   void getPoints(const CellType& cellType, 
00067     Array<Point>& quadPoints,
00068     Array<double>& quadWeights) const ;
00069 
00070   /** Get quadrature points and weights for integration on a facet of a cell */
00071   void getFacetPoints(const CellType& cellType, 
00072     int facetDim,
00073     int facetIndex,
00074     Array<Point>& quadPoints,
00075     Array<double>& quadWeights) const ;
00076 
00077   /** Get the quadrature points and weights for the given cell type ,
00078    * which might be cut by a curve in the case of, see class for more info */
00079   void getAdaptedWeights(const CellType& cellType ,
00080                int cellDim,
00081                int celLID,
00082                int facetIndex ,
00083                  const Mesh& mesh ,
00084                  const ParametrizedCurve& globalCurve ,
00085                  Array<Point>& quadPoints ,
00086                  Array<double>& quadWeights ,
00087                  bool& isCut) const ;
00088 
00089 private:
00090   /** Get quad points for a facet of a line */
00091   void getLineFacetQuad(int facetDim,
00092     int facetIndex,
00093     Array<Point>& quadPoints,
00094     Array<double>& quadWeights) const ;
00095   /** Get quad points for a facet of a triangle */
00096   void getTriangleFacetQuad(int facetDim,
00097     int facetIndex,
00098     Array<Point>& quadPoints,
00099     Array<double>& quadWeights) const ;
00100 
00101   /** Get quad points for a facet of a quadlateral */
00102   void getQuadFacetQuad(int facetDim,
00103     int facetIndex,
00104     Array<Point>& quadPoints,
00105     Array<double>& quadWeights) const ;
00106 
00107   /** Get quad points for a facet of a tet */
00108   void getTetFacetQuad(int facetDim,
00109     int facetIndex,
00110     Array<Point>& quadPoints,
00111     Array<double>& quadWeights) const ;
00112 
00113   /** Get quad points for a facet of a Brick cell */
00114   void getBrickFacetQuad(int facetDim,
00115     int facetIndex,
00116     Array<Point>& quadPoints,
00117     Array<double>& quadWeights) const ;
00118 };
00119 
00120 
00121 /** \relates QuadratureFamily */
00122 void printQuad(std::ostream& os, 
00123   const Array<Point>& pts, const Array<double>& wgts);
00124 
00125 
00126 }
00127 
00128 #endif

Site Contact