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

Site Contact