SundanceQuadratureIntegral.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_QUADRATUREINTEGRAL_H
00043 #define SUNDANCE_QUADRATUREINTEGRAL_H
00044 
00045 #include "SundanceDefs.hpp"
00046 #include "SundanceQuadratureIntegralBase.hpp"
00047 
00048 namespace Sundance
00049 {
00050 
00051 using namespace Teuchos;
00052 
00053 /** 
00054  *  
00055  *
00056  */
00057 class QuadratureIntegral 
00058   : public QuadratureIntegralBase
00059 {
00060 public:
00061   /** Construct a zero-form to be computed by quadrature */
00062   QuadratureIntegral(int spatialDim,
00063     const CellType& maxCellType,
00064     int dim, 
00065     const CellType& cellType,
00066     const QuadratureFamily& quad,
00067     bool isInternalBdry,
00068     const ParametrizedCurve& globalCurve,
00069     const Mesh& mesh,
00070     int verb);
00071 
00072   /** Construct a one form to be computed by quadrature */
00073   QuadratureIntegral(int spatialDim,
00074     const CellType& maxCellType,
00075     int dim, 
00076     const CellType& cellType,
00077     const BasisFamily& testBasis,
00078     int alpha,
00079     int testDerivOrder,
00080     const QuadratureFamily& quad,
00081     bool isInternalBdry,
00082     const ParametrizedCurve& globalCurve,
00083     const Mesh& mesh,
00084     int verb);
00085 
00086   /** Construct a two-form to be computed by quadrature */
00087   QuadratureIntegral(int spatialDim,
00088     const CellType& maxCellType,
00089     int dim,
00090     const CellType& cellType,
00091     const BasisFamily& testBasis,
00092     int alpha,
00093     int testDerivOrder,
00094     const BasisFamily& unkBasis,
00095     int beta,
00096     int unkDerivOrder,
00097     const QuadratureFamily& quad,
00098     bool isInternalBdry,
00099     const ParametrizedCurve& globalCurve,
00100     const Mesh& mesh,
00101     int verb);
00102 
00103   /** virtual dtor */
00104   virtual ~QuadratureIntegral(){;}
00105 
00106   /** */
00107   virtual void transformZeroForm(const CellJacobianBatch& JTrans,
00108          const CellJacobianBatch& JVol,
00109          const Array<int>& isLocalFlag,
00110          const Array<int>& facetIndex,
00111          const RCP<Array<int> >& cellLIDs,
00112          const double* const coeff,
00113          RCP<Array<double> >& A) const ;
00114       
00115   /** */
00116   virtual void transformTwoForm(const CellJacobianBatch& JTrans,
00117         const CellJacobianBatch& JVol,
00118         const Array<int>& facetIndex,
00119           const RCP<Array<int> >& cellLIDs,
00120         const double* const coeff,
00121         RCP<Array<double> >& A) const ;
00122       
00123   /** */
00124   void transformOneForm(const CellJacobianBatch& JTrans,
00125       const CellJacobianBatch& JVol,
00126       const Array<int>& facetIndex,
00127         const RCP<Array<int> >& cellLIDs,
00128       const double* const coeff,
00129       RCP<Array<double> >& A) const ;
00130 
00131 private:
00132 
00133   /** Do the integration by summing reference quantities over quadrature
00134    * points and then transforming the sum to physical quantities.  */
00135   void transformSummingFirst(int nCells,
00136     const Array<int>& facetIndex,
00137     const RCP<Array<int> >& cellLIDs,
00138     const double* const GPtr,
00139     const double* const coeff,
00140     RCP<Array<double> >& A) const ;
00141 
00142   /** Do the integration by transforming to physical coordinates 
00143    * at each quadrature point, and then summing */
00144   void transformSummingLast(int nCells,
00145     const Array<int>& facetIndex,
00146     const RCP<Array<int> >& cellLIDs,
00147     const double* const GPtr,
00148     const double* const coeff,
00149     RCP<Array<double> >& A) const ;
00150 
00151   /** Determine whether to do this batch of integrals using the
00152    * sum-first method or the sum-last method */
00153   bool useSumFirstMethod() const {return useSumFirstMethod_;}
00154       
00155   /** */
00156   inline double& wValue(int facetCase, 
00157     int q, int testDerivDir, int testNode,
00158     int unkDerivDir, int unkNode)
00159     {return W_[facetCase][unkNode
00160         + nNodesUnk()
00161         *(testNode + nNodesTest()
00162           *(unkDerivDir + nRefDerivUnk()
00163             *(testDerivDir + nRefDerivTest()*q)))];}
00164 
00165       
00166 
00167   /** */
00168   inline const double& wValue(int facetCase, 
00169     int q, 
00170     int testDerivDir, int testNode,
00171     int unkDerivDir, int unkNode) const 
00172     {
00173       return W_[facetCase][unkNode
00174         + nNodesUnk()
00175         *(testNode + nNodesTest()
00176           *(unkDerivDir + nRefDerivUnk()
00177             *(testDerivDir + nRefDerivTest()*q)))];
00178     }
00179       
00180   /** */
00181   inline double& wValue(int facetCase, 
00182     int q, int testDerivDir, int testNode)
00183     {return W_[facetCase][testNode + nNodesTest()*(testDerivDir + nRefDerivTest()*q)];}
00184 
00185 
00186   /** */
00187   inline const double& wValue(int facetCase, 
00188     int q, int testDerivDir, int testNode) const 
00189     {return W_[facetCase][testNode + nNodesTest()*(testDerivDir + nRefDerivTest()*q)];}
00190 
00191   /* */
00192   Array<Array<double> > W_;
00193 
00194   /* */
00195   bool useSumFirstMethod_;
00196       
00197   /** For ACI (ACI = Adaptive Cell Integration), store the reference integral values for one form
00198    * The indexes facet, quadPoints, nRefDerivTest , nNodesTest */
00199   Array<Array<Array<Array<double> > > > W_ACI_F1_;
00200 
00201   /** For ACI (ACI = Adaptive Cell Integration), store the reference integral values for two form
00202    * The indexes facet, quadPoints, nRefDerivTest , nNodesTest , nRefDerivUnk , nNodesUnk */
00203   Array<Array<Array<Array<Array<Array<double> > > > > > W_ACI_F2_;
00204 
00205   /** The quadrature family needed for special integration (ACI)*/
00206   QuadratureFamily quad_;
00207 
00208   /** The quadrature points*/
00209   Array < Array<Point> > quadPts_;
00210 
00211   /** The standard weights (in case of ACI these might change)*/
00212   Array < Array<double> > quadWeights_;
00213 };
00214 }
00215 
00216 
00217 #endif

Site Contact