SundanceRefIntegral.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_REFINTEGRAL_H
00043 #define SUNDANCE_REFINTEGRAL_H
00044 
00045 #include "SundanceDefs.hpp"
00046 #include "SundanceElementIntegral.hpp"
00047 
00048 
00049 namespace Sundance
00050 {
00051 
00052 using namespace Teuchos;
00053 
00054 /** 
00055  * RefIntegral represents the integrals of a product 
00056  * of basis functions (or their derivatives) over a reference cell. 
00057  * This object can be created once and for all up front, and then 
00058  * reused to produce the integrals of constant-coefficient weak forms
00059  * over simplicial cells by means of a linear transformation.
00060  * 
00061  * An instance of this object can represent either a one form,
00062  * \f[
00063  * W_{i\gamma} = \int_{T_R} D_\gamma \psi_i 
00064  * \f]  
00065  * or a two form,
00066  * \f[
00067  * W_{(ij)(\gamma\delta)} = \int_{T_R} D_\gamma \psi_i D_\delta \phi_j.
00068  * \f]  
00069  * In the notation for the two form, we have grouped together the index
00070  * pairs \f$(ij)\f$ and \f$(\gamma\delta)\f$ to indicate that we 
00071  * will think of this 4-tensor as a 2D array, with \f$(ij)\f$ combinations
00072  * running down rows and \f$(\gamma\delta)\f$ running across columns. 
00073  * We will consider the node number combination 
00074  * \f$(ij)\f$ to be a single index \f$k\f$, and 
00075  * \f$(\gamma\delta)\f$ to be a single index \f$\epsilon\f$. Thus, 
00076  * the storage and use of one-forms and two-forms is essentially identical.
00077  * 
00078  * This object's job in life is to be multiplied with a transformation
00079  * matrix \f$T_{\epsilon c}\f$to produce an array of local element matrices
00080  * \f$A_{kc}\f$,
00081  * \f[
00082  * A_{kc} = W_{k\epsilon} T_{\epsilon c}
00083  * \f]
00084  * The index \f$c\f$ is over cells (cells are processed in batches).
00085  * 
00086  * Physical storage is as a 1D vector stored in colum-major order. 
00087  */
00088 class RefIntegral : public ElementIntegral
00089 {
00090 public:
00091   /** Construct a reference zero-form */
00092   RefIntegral(int spatialDim,
00093     const CellType& maxCellType,
00094     int dim, 
00095     const CellType& cellType,
00096     const QuadratureFamily& quad_in,
00097     bool isInternalBdry,
00098     const ParametrizedCurve& globalCurve,
00099     const Mesh& mesh,
00100     int verb);
00101 
00102   /** Construct a reference one-form */
00103   RefIntegral(int spatialDim,
00104     const CellType& maxCellType,
00105     int dim, 
00106     const CellType& cellType,
00107     const BasisFamily& testBasis,
00108     int alpha,
00109     int testDerivOrder,
00110     const QuadratureFamily& quad_in,
00111     bool isInternalBdry,
00112     const ParametrizedCurve& globalCurve,
00113     const Mesh& mesh,
00114     int verb);
00115 
00116   /** Construct a reference two-form */
00117   RefIntegral(int spatialDim,
00118     const CellType& maxCellType,
00119     int dim,
00120     const CellType& cellType,
00121     const BasisFamily& testBasis,
00122     int alpha,
00123     int testDerivOrder,
00124     const BasisFamily& unkBasis,
00125     int beta,
00126     int unkDerivOrder,
00127     const QuadratureFamily& quad_in,
00128     bool isInternalBdry,
00129     const ParametrizedCurve& globalCurve,
00130     const Mesh& mesh,
00131     int verb);
00132 
00133   /** virtual dtor */
00134   virtual ~RefIntegral(){;}
00135       
00136   /** */
00137   void print(std::ostream& os) const ;
00138 
00139   /** */
00140   void transform(const CellJacobianBatch& JTrans,
00141     const CellJacobianBatch& JVol,
00142     const Array<int>& isLocalFlag,
00143     const Array<int>& facetNum,
00144     const RCP<Array<int> >& cellLIDs,
00145     const double& coeff,
00146     RCP<Array<double> >& A) const
00147     {
00148       if (order()==2) transformTwoForm(JTrans, JVol, facetNum, cellLIDs, coeff, A);
00149       else if (order()==1) transformOneForm(JTrans, JVol, facetNum, cellLIDs, coeff, A);
00150       else transformZeroForm(JVol, isLocalFlag, cellLIDs, coeff, A);
00151     }
00152 
00153   /** */
00154   void transformTwoForm(const CellJacobianBatch& JTrans,
00155     const CellJacobianBatch& JVol,
00156     const Array<int>& facetNum, 
00157     const RCP<Array<int> >& cellLIDs,
00158     const double& coeff,
00159     RCP<Array<double> >& A) const ;
00160 
00161   /** */
00162   void transformOneForm(const CellJacobianBatch& JTrans,
00163     const CellJacobianBatch& JVol,
00164     const Array<int>& facetNum, 
00165     const RCP<Array<int> >& cellLIDs,
00166     const double& coeff,
00167     RCP<Array<double> >& A) const ;
00168 
00169   /** */
00170   void transformZeroForm(const CellJacobianBatch& JVol,
00171     const Array<int>& isLocalFlag,
00172     const RCP<Array<int> >& cellLIDs,
00173     const double& coeff,
00174     RCP<Array<double> >& A) const ;
00175 
00176   /** */
00177   inline double& value(int facetCase, int testDerivDir, int testNode,
00178     int unkDerivDir, int unkNode)
00179     {return W_[facetCase][unkNode + nNodesUnk()*testNode 
00180         + nNodes()*(unkDerivDir 
00181           + nRefDerivUnk()*testDerivDir)];}
00182 
00183   /** */
00184   inline const double& value(int facetCase, 
00185     int testDerivDir, int testNode,
00186     int unkDerivDir, int unkNode) const 
00187     {
00188       return W_[facetCase][unkNode + nNodesUnk()*testNode 
00189         + nNodes()*(unkDerivDir 
00190           + nRefDerivUnk()*testDerivDir)];
00191     }
00192       
00193   /** */
00194   inline double& value(int facetCase, int testDerivDir, int testNode)
00195     {return W_[facetCase][nNodesTest()*testDerivDir + testNode];}
00196 
00197   /** */
00198   inline const double& value(int facetCase, 
00199     int testDerivDir, int testNode) const 
00200     {return W_[facetCase][nNodesTest()*testDerivDir + testNode];}
00201 
00202   static double& totalFlops() {static double rtn = 0; return rtn;}
00203 
00204 protected:
00205 
00206   static void addFlops(const double& flops) {totalFlops() += flops;}
00207       
00208 private:
00209 
00210   Array<Array<double> > W_;
00211 
00212   /** For ACI (ACI = Adaptive Cell Integration), store the reference integral values for one form
00213    * The indexes facet, quadPoints, nRefDerivTest , nNodesTest */
00214   Array<Array<Array<Array<double> > > > W_ACI_F1_;
00215 
00216   /** For ACI (ACI = Adaptive Cell Integration), store the reference integral values for two form
00217    * The indexes facet, quadPoints, nRefDerivTest , nNodesTest , nRefDerivUnk , nNodesUnk */
00218   Array<Array<Array<Array<Array<Array<double> > > > > > W_ACI_F2_;
00219 
00220   /** The quadrature family needed for special integration */
00221   QuadratureFamily quad_;
00222 
00223   /** The quadrature points*/
00224   Array<Point> quadPts_;
00225 
00226   /** The standard weights (in case of ACI these might change)*/
00227   Array<double> quadWeights_;
00228 };
00229 }
00230 
00231 #endif

Site Contact