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

Site Contact