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

Site Contact