Intrepid
http://trilinos.sandia.gov/packages/docs/r10.10/packages/intrepid/test/Discretization/Integration/test_03.cpp
Go to the documentation of this file.
00001 // @HEADER
00002 // ************************************************************************
00003 //
00004 //                           Intrepid Package
00005 //                 Copyright (2007) Sandia Corporation
00006 //
00007 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
00008 // license for use of this work by or on behalf of the U.S. Government.
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 Pavel Bochev  (pbboche@sandia.gov)
00038 //                    Denis Ridzal  (dridzal@sandia.gov), or
00039 //                    Kara Peterson (kjpeter@sandia.gov)
00040 //
00041 // ************************************************************************
00042 // @HEADER
00043 
00044 
00051 #include "Intrepid_DefaultCubatureFactory.hpp"
00052 #include "Intrepid_Utils.hpp"
00053 #include "Teuchos_oblackholestream.hpp"
00054 #include "Teuchos_RCP.hpp"
00055 #include "Teuchos_GlobalMPISession.hpp"
00056 
00057 using namespace Intrepid;
00058 
00059 
00060 /*
00061   Monomial evaluation.
00062     in 1D, for point p(x)    : x^xDeg
00063     in 2D, for point p(x,y)  : x^xDeg * y^yDeg
00064     in 3D, for point p(x,y,z): x^xDeg * y^yDeg * z^zDeg
00065 */
00066 double computeMonomial(FieldContainer<double> & p, int xDeg, int yDeg=0, int zDeg=0) {
00067   double val = 1.0;
00068   int polydeg[3];
00069   polydeg[0] = xDeg; polydeg[1] = yDeg; polydeg[2] = zDeg;
00070   for (int i=0; i<p.dimension(0); i++) {
00071     val *= std::pow(p(i),polydeg[i]);
00072   }
00073   return val;
00074 }
00075 
00076 
00077 /*
00078   Computes integrals of monomials over a given reference cell.
00079 */
00080 double computeIntegral(shards::CellTopology & cellTopology, int cubDegree, int xDeg, int yDeg) {
00081  
00082   DefaultCubatureFactory<double>  cubFactory;                                         // create factory
00083   Teuchos::RCP<Cubature<double> > myCub = cubFactory.create(cellTopology, cubDegree); // create default cubature
00084 
00085   double val       = 0.0;
00086   int cubDim       = myCub->getDimension();
00087   int numCubPoints = myCub->getNumPoints();
00088 
00089   FieldContainer<double> point(cubDim);
00090   FieldContainer<double> cubPoints(numCubPoints, cubDim);
00091   FieldContainer<double> cubWeights(numCubPoints);
00092 
00093   myCub->getCubature(cubPoints, cubWeights);
00094 
00095   for (int i=0; i<numCubPoints; i++) {
00096     for (int j=0; j<cubDim; j++) {
00097       point(j) = cubPoints(i,j);
00098     }
00099     val += computeMonomial(point, xDeg, yDeg)*cubWeights(i);
00100   }
00101 
00102   return val;
00103 }
00104 
00105 
00106 int main(int argc, char *argv[]) {
00107 
00108   Teuchos::GlobalMPISession mpiSession(&argc, &argv);
00109 
00110   // This little trick lets us print to std::cout only if
00111   // a (dummy) command-line argument is provided.
00112   int iprint     = argc - 1;
00113   Teuchos::RCP<std::ostream> outStream;
00114   Teuchos::oblackholestream bhs; // outputs nothing
00115   if (iprint > 0)
00116     outStream = Teuchos::rcp(&std::cout, false);
00117   else
00118     outStream = Teuchos::rcp(&bhs, false);
00119 
00120   // Save the format state of the original std::cout.
00121   Teuchos::oblackholestream oldFormatState;
00122   oldFormatState.copyfmt(std::cout);
00123  
00124   *outStream \
00125   << "===============================================================================\n" \
00126   << "|                                                                             |\n" \
00127   << "|     Unit Test (CubatureDirect,CubatureTensor,DefaultCubatureFactory)        |\n" \
00128   << "|                                                                             |\n" \
00129   << "|     1) Computing integrals of monomials on reference cells in 2D            |\n" \
00130   << "|                                                                             |\n" \
00131   << "|  Questions? Contact  Pavel Bochev (pbboche@sandia.gov) or                   |\n" \
00132   << "|                      Denis Ridzal (dridzal@sandia.gov).                     |\n" \
00133   << "|                                                                             |\n" \
00134   << "|  Intrepid's website: http://trilinos.sandia.gov/packages/intrepid           |\n" \
00135   << "|  Trilinos website:   http://trilinos.sandia.gov                             |\n" \
00136   << "|                                                                             |\n" \
00137   << "===============================================================================\n"\
00138   << "| TEST 1: integrals of monomials in 2D                                        |\n"\
00139   << "===============================================================================\n";
00140 
00141   // internal variables:
00142   int                                      errorFlag = 0;
00143   int                                      polyCt = 0;
00144   int                                      offset = 0;
00145   Teuchos::Array< Teuchos::Array<double> > testInt;
00146   Teuchos::Array< Teuchos::Array<double> > analyticInt;
00147   Teuchos::Array<double>                   tmparray(1);
00148   double                                   reltol = 1.0e+03 * INTREPID_TOL;
00149   int                                      maxDeg[2];
00150   int                                      numPoly[2];
00151   maxDeg[0]  = INTREPID_CUBATURE_TRI_DEFAULT_MAX;
00152   maxDeg[1]  = INTREPID_CUBATURE_LINE_GAUSS_MAX;
00153   numPoly[0] = (INTREPID_CUBATURE_TRI_DEFAULT_MAX+1)*(INTREPID_CUBATURE_TRI_DEFAULT_MAX+2)/2;
00154   numPoly[1] = (INTREPID_CUBATURE_LINE_GAUSS_MAX+1)*(INTREPID_CUBATURE_LINE_GAUSS_MAX+2)/2;
00155 
00156   // get names of files with analytic values
00157   std::string basedir = "./data";
00158   std::stringstream namestream[2];
00159   std::string filename[2];
00160   namestream[0] << basedir << "/TRI_integrals" << ".dat";
00161   namestream[0] >> filename[0];
00162   namestream[1] << basedir << "/QUAD_integrals" << ".dat";
00163   namestream[1] >> filename[1];
00164 
00165   shards::CellTopology cellType[] = {shards::getCellTopologyData< shards::Triangle<> >(),
00166                                      shards::getCellTopologyData< shards::Quadrilateral<> >()};
00167 
00168   // compute and compare integrals
00169   try {
00170     for (int cellCt=0; cellCt < 2; cellCt++) {
00171       testInt.assign(numPoly[cellCt], tmparray);
00172       analyticInt.assign(numPoly[cellCt], tmparray);
00173       *outStream << "\nIntegrals of monomials on a reference " << cellType[cellCt].getBaseCellTopologyData()->name << ":\n";
00174       std::ifstream filecompare(&filename[cellCt][0]);
00175       // compute integrals
00176       for (int cubDeg=0; cubDeg <= maxDeg[cellCt]; cubDeg++) {
00177         polyCt = 0;
00178         testInt[cubDeg].resize((cubDeg+1)*(cubDeg+2)/2);
00179         for (int xDeg=0; xDeg <= cubDeg; xDeg++) {
00180           for (int yDeg=0; yDeg <= cubDeg-xDeg; yDeg++) {
00181             testInt[cubDeg][polyCt] = computeIntegral(cellType[cellCt], cubDeg, xDeg, yDeg);
00182             polyCt++; 
00183           }
00184         }
00185       }
00186       // get analytic values
00187       if (filecompare.is_open()) {
00188         getAnalytic(analyticInt, filecompare);
00189         // close file
00190         filecompare.close();
00191       }
00192       // perform comparison
00193       for (int cubDeg=0; cubDeg <= maxDeg[cellCt]; cubDeg++) {
00194         polyCt = 0;
00195         offset = 0;
00196         for (int xDeg=0; xDeg <= cubDeg; xDeg++) {
00197           for (int yDeg=0; yDeg <= cubDeg-xDeg; yDeg++) {
00198             double abstol = ( analyticInt[polyCt+offset][0] == 0.0 ? reltol : std::fabs(reltol*analyticInt[polyCt+offset][0]) );
00199             double absdiff = std::fabs(analyticInt[polyCt+offset][0] - testInt[cubDeg][polyCt]);
00200             *outStream << "Cubature order " << std::setw(2) << std::left << cubDeg << " integrating "
00201                        << "x^" << std::setw(2) << std::left << xDeg << " * y^" << std::setw(2) << yDeg << ":" << "   "
00202                        << std::scientific << std::setprecision(16) << testInt[cubDeg][polyCt] << "   " << analyticInt[polyCt+offset][0] << "   "
00203                        << std::setprecision(4) << absdiff << "   " << "<?" << "   " << abstol << "\n";
00204             if (absdiff > abstol) {
00205               errorFlag++;
00206               *outStream << std::right << std::setw(111) << "^^^^---FAILURE!\n";
00207             }
00208             polyCt++;
00209           }
00210           offset = offset + maxDeg[cellCt] - cubDeg;
00211         }
00212         *outStream << "\n";
00213       }
00214       *outStream << "\n";
00215     }  // end for cellCt
00216   }
00217   catch (std::logic_error err) {
00218     *outStream << err.what() << "\n";
00219     errorFlag = -1;
00220   };
00221 
00222 
00223   if (errorFlag != 0)
00224     std::cout << "End Result: TEST FAILED\n";
00225   else
00226     std::cout << "End Result: TEST PASSED\n";
00227 
00228   // reset format state of std::cout
00229   std::cout.copyfmt(oldFormatState);
00230 
00231   return errorFlag;
00232 }