SundanceReducedIntegral.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_REDUCED_INTEGRAL_H
00032 #define SUNDANCE_REDUCED_INTEGRAL_H
00033 
00034 #include "SundanceDefs.hpp"
00035 #include "SundanceElementIntegral.hpp"
00036 
00037 namespace Sundance
00038 {
00039 
00040 using namespace Teuchos;
00041 
00042 /** 
00043  * 
00044  */
00045 class ReducedIntegral : public ElementIntegral
00046 {
00047 public:
00048   /** Construct a zero-form to be computed by reference integration
00049    * with coefficients that are piecewise constant */
00050   ReducedIntegral(int spatialDim,
00051     const CellType& maxCellType,
00052     int dim, 
00053     const CellType& cellType,
00054     const QuadratureFamily& quad,
00055     bool isInternalBdry,
00056     const ParametrizedCurve& globalCurve,
00057     const Mesh& mesh,
00058     int verb);
00059 
00060   /** Construct a one form to be computed by reference integration
00061    * with coefficients that are piecewise constant  */
00062   ReducedIntegral(int spatialDim,
00063     const CellType& maxCellType,
00064     int dim, 
00065     const CellType& cellType,
00066     const BasisFamily& testBasis,
00067     int alpha,
00068     int testDerivOrder,
00069     const QuadratureFamily& quad,
00070     bool isInternalBdry,
00071     const ParametrizedCurve& globalCurve,
00072     const Mesh& mesh,
00073     int verb);
00074 
00075   /** Construct a two-form to be computed by reference integration
00076    * with coefficients that are piecewise constant */
00077   ReducedIntegral(int spatialDim,
00078     const CellType& maxCellType,
00079     int dim,
00080     const CellType& cellType,
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     bool isInternalBdry,
00089     const ParametrizedCurve& globalCurve,
00090     const Mesh& mesh,
00091     int verb);
00092 
00093   /** virtual dtor */
00094   virtual ~ReducedIntegral(){;}
00095 
00096   /** */
00097   void transform(const CellJacobianBatch& JTrans,
00098     const CellJacobianBatch& JVol,
00099     const Array<int>& isLocalFlag,
00100     const Array<int>& facetNum,
00101     const RCP<Array<int> >& cellLIDs,
00102     const double* const coeffs,
00103     RCP<Array<double> >& A) const
00104     {
00105       if (order()==2) transformTwoForm(JTrans, JVol, facetNum, cellLIDs, coeffs, A);
00106       else if (order()==1) transformOneForm(JTrans, JVol, facetNum, cellLIDs, coeffs, A);
00107       else transformZeroForm(JTrans, JVol, isLocalFlag, facetNum,
00108         cellLIDs, coeffs, A);
00109     }
00110 
00111   /** */
00112   virtual void transformZeroForm(const CellJacobianBatch& JTrans,
00113     const CellJacobianBatch& JVol,
00114     const Array<int>& isLocalFlag,
00115     const Array<int>& facetIndex,
00116     const RCP<Array<int> >& cellLIDs,
00117     const double* const coeffs,
00118     RCP<Array<double> >& A) const ;
00119       
00120   /** */
00121   virtual void transformTwoForm(const CellJacobianBatch& JTrans,
00122     const CellJacobianBatch& JVol,
00123     const Array<int>& facetIndex,
00124     const RCP<Array<int> >& cellLIDs,
00125     const double* const coeffs,
00126     RCP<Array<double> >& A) const ;
00127       
00128   /** */
00129   void transformOneForm(const CellJacobianBatch& JTrans,
00130     const CellJacobianBatch& JVol,
00131     const Array<int>& facetIndex,
00132     const RCP<Array<int> >& cellLIDs,
00133     const double* const coeffs,
00134     RCP<Array<double> >& A) const ;
00135 
00136 private:
00137 
00138   
00139 
00140   /** */
00141   inline double& value(int facetCase, int testDerivDir, int testNode,
00142     int unkDerivDir, int unkNode)
00143     {return W_[facetCase][unkNode + nNodesUnk()*testNode 
00144         + nNodes()*(unkDerivDir 
00145           + nRefDerivUnk()*testDerivDir)];}
00146 
00147   /** */
00148   inline const double& value(int facetCase, 
00149     int testDerivDir, int testNode,
00150     int unkDerivDir, int unkNode) const 
00151     {
00152       return W_[facetCase][unkNode + nNodesUnk()*testNode 
00153         + nNodes()*(unkDerivDir 
00154           + nRefDerivUnk()*testDerivDir)];
00155     }
00156       
00157   /** */
00158   inline double& value(int facetCase, int testDerivDir, int testNode)
00159     {return W_[facetCase][nNodesTest()*testDerivDir + testNode];}
00160 
00161   /** */
00162   inline const double& value(int facetCase, 
00163     int testDerivDir, int testNode) const 
00164     {return W_[facetCase][nNodesTest()*testDerivDir + testNode];}
00165 
00166   static double& totalFlops() {static double rtn = 0; return rtn;}
00167 
00168 protected:
00169 
00170   static void addFlops(const double& flops) {totalFlops() += flops;}
00171       
00172 private:
00173 
00174   Array<Array<double> > W_;
00175 
00176 };
00177 }
00178 
00179 
00180 #endif

Site Contact