Intrepid
/usr/src/RPM/BUILD/trilinos10-10.6.4/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 // This library is free software; you can redistribute it and/or modify
00011 // it under the terms of the GNU Lesser General Public License as
00012 // published by the Free Software Foundation; either version 2.1 of the
00013 // License, or (at your option) any later version.
00014 //
00015 // This library is distributed in the hope that it will be useful, but
00016 // WITHOUT ANY WARRANTY; without even the implied warranty of
00017 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00018 // Lesser General Public License for more details.
00019 //
00020 // You should have received a copy of the GNU Lesser General Public
00021 // License along with this library; if not, write to the Free Software
00022 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307
00023 // USA
00024 // Questions? Contact Pavel Bochev (pbboche@sandia.gov) or
00025 //                    Denis Ridzal (dridzal@sandia.gov).
00026 //
00027 // ************************************************************************
00028 // @HEADER
00029 
00030 
00037 #include "Intrepid_DefaultCubatureFactory.hpp"
00038 #include "Intrepid_Utils.hpp"
00039 #include "Teuchos_oblackholestream.hpp"
00040 #include "Teuchos_RCP.hpp"
00041 #include "Teuchos_GlobalMPISession.hpp"
00042 
00043 using namespace Intrepid;
00044 
00045 
00046 /*
00047   Monomial evaluation.
00048     in 1D, for point p(x)    : x^xDeg
00049     in 2D, for point p(x,y)  : x^xDeg * y^yDeg
00050     in 3D, for point p(x,y,z): x^xDeg * y^yDeg * z^zDeg
00051 */
00052 double computeMonomial(FieldContainer<double> & p, int xDeg, int yDeg=0, int zDeg=0) {
00053   double val = 1.0;
00054   int polydeg[3];
00055   polydeg[0] = xDeg; polydeg[1] = yDeg; polydeg[2] = zDeg;
00056   for (int i=0; i<p.dimension(0); i++) {
00057     val *= std::pow(p(i),polydeg[i]);
00058   }
00059   return val;
00060 }
00061 
00062 
00063 /*
00064   Computes integrals of monomials over a given reference cell.
00065 */
00066 double computeIntegral(shards::CellTopology & cellTopology, int cubDegree, int xDeg, int yDeg) {
00067  
00068   DefaultCubatureFactory<double>  cubFactory;                                         // create factory
00069   Teuchos::RCP<Cubature<double> > myCub = cubFactory.create(cellTopology, cubDegree); // create default cubature
00070 
00071   double val       = 0.0;
00072   int cubDim       = myCub->getDimension();
00073   int numCubPoints = myCub->getNumPoints();
00074 
00075   FieldContainer<double> point(cubDim);
00076   FieldContainer<double> cubPoints(numCubPoints, cubDim);
00077   FieldContainer<double> cubWeights(numCubPoints);
00078 
00079   myCub->getCubature(cubPoints, cubWeights);
00080 
00081   for (int i=0; i<numCubPoints; i++) {
00082     for (int j=0; j<cubDim; j++) {
00083       point(j) = cubPoints(i,j);
00084     }
00085     val += computeMonomial(point, xDeg, yDeg)*cubWeights(i);
00086   }
00087 
00088   return val;
00089 }
00090 
00091 
00092 int main(int argc, char *argv[]) {
00093 
00094   Teuchos::GlobalMPISession mpiSession(&argc, &argv);
00095 
00096   // This little trick lets us print to std::cout only if
00097   // a (dummy) command-line argument is provided.
00098   int iprint     = argc - 1;
00099   Teuchos::RCP<std::ostream> outStream;
00100   Teuchos::oblackholestream bhs; // outputs nothing
00101   if (iprint > 0)
00102     outStream = Teuchos::rcp(&std::cout, false);
00103   else
00104     outStream = Teuchos::rcp(&bhs, false);
00105 
00106   // Save the format state of the original std::cout.
00107   Teuchos::oblackholestream oldFormatState;
00108   oldFormatState.copyfmt(std::cout);
00109  
00110   *outStream \
00111   << "===============================================================================\n" \
00112   << "|                                                                             |\n" \
00113   << "|     Unit Test (CubatureDirect,CubatureTensor,DefaultCubatureFactory)        |\n" \
00114   << "|                                                                             |\n" \
00115   << "|     1) Computing integrals of monomials on reference cells in 2D            |\n" \
00116   << "|                                                                             |\n" \
00117   << "|  Questions? Contact  Pavel Bochev (pbboche@sandia.gov) or                   |\n" \
00118   << "|                      Denis Ridzal (dridzal@sandia.gov).                     |\n" \
00119   << "|                                                                             |\n" \
00120   << "|  Intrepid's website: http://trilinos.sandia.gov/packages/intrepid           |\n" \
00121   << "|  Trilinos website:   http://trilinos.sandia.gov                             |\n" \
00122   << "|                                                                             |\n" \
00123   << "===============================================================================\n"\
00124   << "| TEST 1: integrals of monomials in 2D                                        |\n"\
00125   << "===============================================================================\n";
00126 
00127   // internal variables:
00128   int                                      errorFlag = 0;
00129   int                                      polyCt = 0;
00130   int                                      offset = 0;
00131   Teuchos::Array< Teuchos::Array<double> > testInt;
00132   Teuchos::Array< Teuchos::Array<double> > analyticInt;
00133   Teuchos::Array<double>                   tmparray(1);
00134   double                                   reltol = 1.0e+03 * INTREPID_TOL;
00135   int                                      maxDeg[2];
00136   int                                      numPoly[2];
00137   maxDeg[0]  = INTREPID_CUBATURE_TRI_DEFAULT_MAX;
00138   maxDeg[1]  = INTREPID_CUBATURE_LINE_GAUSS_MAX;
00139   numPoly[0] = (INTREPID_CUBATURE_TRI_DEFAULT_MAX+1)*(INTREPID_CUBATURE_TRI_DEFAULT_MAX+2)/2;
00140   numPoly[1] = (INTREPID_CUBATURE_LINE_GAUSS_MAX+1)*(INTREPID_CUBATURE_LINE_GAUSS_MAX+2)/2;
00141 
00142   // get names of files with analytic values
00143   std::string basedir = "./data";
00144   std::stringstream namestream[2];
00145   std::string filename[2];
00146   namestream[0] << basedir << "/TRI_integrals" << ".dat";
00147   namestream[0] >> filename[0];
00148   namestream[1] << basedir << "/QUAD_integrals" << ".dat";
00149   namestream[1] >> filename[1];
00150 
00151   shards::CellTopology cellType[] = {shards::getCellTopologyData< shards::Triangle<> >(),
00152                                      shards::getCellTopologyData< shards::Quadrilateral<> >()};
00153 
00154   // compute and compare integrals
00155   try {
00156     for (int cellCt=0; cellCt < 2; cellCt++) {
00157       testInt.assign(numPoly[cellCt], tmparray);
00158       analyticInt.assign(numPoly[cellCt], tmparray);
00159       *outStream << "\nIntegrals of monomials on a reference " << cellType[cellCt].getBaseTopology()->name << ":\n";
00160       std::ifstream filecompare(&filename[cellCt][0]);
00161       // compute integrals
00162       for (int cubDeg=0; cubDeg <= maxDeg[cellCt]; cubDeg++) {
00163         polyCt = 0;
00164         testInt[cubDeg].resize((cubDeg+1)*(cubDeg+2)/2);
00165         for (int xDeg=0; xDeg <= cubDeg; xDeg++) {
00166           for (int yDeg=0; yDeg <= cubDeg-xDeg; yDeg++) {
00167             testInt[cubDeg][polyCt] = computeIntegral(cellType[cellCt], cubDeg, xDeg, yDeg);
00168             polyCt++; 
00169           }
00170         }
00171       }
00172       // get analytic values
00173       if (filecompare.is_open()) {
00174         getAnalytic(analyticInt, filecompare);
00175         // close file
00176         filecompare.close();
00177       }
00178       // perform comparison
00179       for (int cubDeg=0; cubDeg <= maxDeg[cellCt]; cubDeg++) {
00180         polyCt = 0;
00181         offset = 0;
00182         for (int xDeg=0; xDeg <= cubDeg; xDeg++) {
00183           for (int yDeg=0; yDeg <= cubDeg-xDeg; yDeg++) {
00184             double abstol = ( analyticInt[polyCt+offset][0] == 0.0 ? reltol : std::fabs(reltol*analyticInt[polyCt+offset][0]) );
00185             double absdiff = std::fabs(analyticInt[polyCt+offset][0] - testInt[cubDeg][polyCt]);
00186             *outStream << "Cubature order " << std::setw(2) << std::left << cubDeg << " integrating "
00187                        << "x^" << std::setw(2) << std::left << xDeg << " * y^" << std::setw(2) << yDeg << ":" << "   "
00188                        << std::scientific << std::setprecision(16) << testInt[cubDeg][polyCt] << "   " << analyticInt[polyCt+offset][0] << "   "
00189                        << std::setprecision(4) << absdiff << "   " << "<?" << "   " << abstol << "\n";
00190             if (absdiff > abstol) {
00191               errorFlag++;
00192               *outStream << std::right << std::setw(111) << "^^^^---FAILURE!\n";
00193             }
00194             polyCt++;
00195           }
00196           offset = offset + maxDeg[cellCt] - cubDeg;
00197         }
00198         *outStream << "\n";
00199       }
00200       *outStream << "\n";
00201     }  // end for cellCt
00202   }
00203   catch (std::logic_error err) {
00204     *outStream << err.what() << "\n";
00205     errorFlag = -1;
00206   };
00207 
00208 
00209   if (errorFlag != 0)
00210     std::cout << "End Result: TEST FAILED\n";
00211   else
00212     std::cout << "End Result: TEST PASSED\n";
00213 
00214   // reset format state of std::cout
00215   std::cout.copyfmt(oldFormatState);
00216 
00217   return errorFlag;
00218 }