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

Site Contact