SundanceCurveQuadratureIntegral.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_CURVEQUADRATUREINTEGRAL_H
00043 #define SUNDANCE_CURVEQUADRATUREINTEGRAL_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 CurveQuadratureIntegral : public ElementIntegral
00057 {
00058 public:
00059   /** Construct a zero-form to be computed by quadrature */
00060   CurveQuadratureIntegral(
00061     const CellType& maxCellType,
00062     const bool isConstantIntegral,
00063     const QuadratureFamily& quad,
00064     const ParametrizedCurve& globalCurve,
00065     const Mesh& mesh,
00066     int verb);
00067 
00068   /** Construct a one form to be computed by quadrature */
00069   CurveQuadratureIntegral(
00070     const CellType& maxCellType,
00071     const bool isConstantIntegral,
00072     const BasisFamily& testBasis,
00073     int alpha,
00074     int testDerivOrder,
00075     const QuadratureFamily& quad,
00076     const ParametrizedCurve& globalCurve,
00077     const Mesh& mesh,
00078     int verb);
00079 
00080   /** Construct a two-form to be computed by quadrature */
00081   CurveQuadratureIntegral(
00082     const CellType& maxCellType,
00083     const bool isConstantIntegral,
00084     const BasisFamily& testBasis,
00085     int alpha,
00086     int testDerivOrder,
00087     const BasisFamily& unkBasis,
00088     int beta,
00089     int unkDerivOrder,
00090     const QuadratureFamily& quad,
00091     const ParametrizedCurve& globalCurve,
00092     const Mesh& mesh,
00093     int verb);
00094 
00095   /** virtual dtor */
00096   virtual ~CurveQuadratureIntegral(){;}
00097       
00098      /** */
00099   virtual void transform(const CellJacobianBatch& JTrans,
00100     const CellJacobianBatch& JVol,
00101     const Array<int>& isLocalFlag,
00102     const Array<int>& facetNum,
00103     const RCP<Array<int> >& cellLIDs,
00104     const double constCoeff,
00105     const double* const coeff,
00106     RCP<Array<double> >& A) const 
00107     {
00108       if (order()==2) transformTwoForm(JTrans, JVol, facetNum, cellLIDs, constCoeff, coeff, A);
00109       else if (order()==1) transformOneForm(JTrans, JVol, facetNum, cellLIDs, constCoeff, coeff, A);
00110       else transformZeroForm(JTrans, JVol, isLocalFlag, facetNum, cellLIDs, constCoeff, coeff, A);
00111     }
00112 
00113   /** */
00114   virtual void transformZeroForm(const CellJacobianBatch& JTrans,
00115     const CellJacobianBatch& JVol,
00116     const Array<int>& isLocalFlag,
00117     const Array<int>& facetIndex,
00118     const RCP<Array<int> >& cellLIDs,
00119     const double constCoeff,
00120     const double* const coeff,
00121     RCP<Array<double> >& A) const ;
00122   
00123   /** */
00124   virtual void transformTwoForm(const CellJacobianBatch& JTrans,
00125     const CellJacobianBatch& JVol,
00126     const Array<int>& facetIndex,
00127     const RCP<Array<int> >& cellLIDs,
00128     const double constCoeff,
00129     const double* const coeff,
00130     RCP<Array<double> >& A) const ;
00131   
00132   /** */
00133   void transformOneForm(const CellJacobianBatch& JTrans,
00134     const CellJacobianBatch& JVol,
00135     const Array<int>& facetIndex,
00136     const RCP<Array<int> >& cellLIDs,
00137     const double constCoeff,
00138     const double* const coeff,
00139     RCP<Array<double> >& A) const ;
00140   
00141 private:
00142 
00143   /** updates the reference cell information (quadPoints,derivatives,normals)*/
00144   void updateRefCellInformation(int maxCellLID , const ParametrizedCurve& curve) const;
00145 
00146   /** updates the W_ for the given cell for one form Integral */
00147   void updateRefCellIntegralOneForm(int maxCellLID , int cellInBatch) const ;
00148 
00149   /** updates the W_ for the given cell for two form Integral */
00150   void updateRefCellIntegralTwoForm(int maxCellLID , int cellInBatch) const ;
00151 
00152   /** Do the integration by summing reference quantities over quadrature
00153    * points and then transforming the sum to physical quantities.  */
00154   void transformSummingFirst(int nCells,
00155   const CellJacobianBatch& JVol,
00156     const Array<int>& facetIndex,
00157     const RCP<Array<int> >& cellLIDs,
00158     const double constCoeff,
00159     const double* const GPtr,
00160     const double* const coeff,
00161     RCP<Array<double> >& A) const ;
00162 
00163   /** Do the integration by transforming to physical coordinates 
00164    * at each quadrature point, and then summing */
00165   void transformSummingLast(int nCells,
00166     const Array<int>& facetIndex,
00167     const RCP<Array<int> >& cellLIDs,
00168     const double* const GPtr,
00169     const double* const coeff,
00170     RCP<Array<double> >& A) const ;
00171 
00172   /** Determine whether to do this batch of integrals using the
00173    * sum-first method or the sum-last method */
00174   bool useSumFirstMethod() const {return useSumFirstMethod_;}
00175       
00176   /** */
00177   inline double& wValue(int q, int testDerivDir, int testNode,
00178     int unkDerivDir, int unkNode) const
00179     {return W_[unkNode
00180         + nNodesUnk()
00181         *(testNode + nNodesTest()
00182           *(unkDerivDir + nRefDerivUnk()
00183             *(testDerivDir + nRefDerivTest()*q)))];}
00184 
00185       
00186 
00187   /** */
00188   inline const double& wValue(int facetCase, 
00189     int q, 
00190     int testDerivDir, int testNode,
00191     int unkDerivDir, int unkNode) const 
00192     {
00193       return W_[unkNode
00194         + nNodesUnk()
00195         *(testNode + nNodesTest()
00196           *(unkDerivDir + nRefDerivUnk()
00197             *(testDerivDir + nRefDerivTest()*q)))];
00198     }
00199       
00200   /** */
00201   inline double& wValue(int q, int testDerivDir, int testNode) const
00202     {return W_[testNode + nNodesTest()*(testDerivDir + nRefDerivTest()*q)];}
00203 
00204   /** The quadrature family needed for curve integration */
00205   QuadratureFamily quad_;
00206 
00207   /** The quadrature points*/
00208   mutable Array<Point> quadPts_;
00209 
00210   /** The standard weights */
00211   Array<double> quadWeights_;
00212 
00213   /* this must be changeable, because each cell will have different reference cell values <br>
00214    * here we store the result from the reference basis evaluation*/
00215   mutable Array<double> W_;
00216 
00217   /* */
00218   bool useSumFirstMethod_;
00219 
00220   /** weather to use constant coefficients */
00221   bool useConstCoeff_;
00222 
00223   /** The derivative of the curve at the curve */
00224   mutable Array<Point> quadCurveDerivs_;
00225 
00226   /** The derivative of the curve at the curve */
00227   mutable Array<Point> quadCurveNormals_;
00228 
00229 };
00230 }
00231 
00232 
00233 #endif /* SUNDANCE_CURVEQUADRATUREINTEGRAL_H */

Site Contact