Intrepid
/usr/src/RPM/BUILD/trilinos10-10.6.4/packages/intrepid/test/Discretization/Integration/test_04.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   Monomial evaluation.
00047     in 1D, for point p(x)    : x^xDeg
00048     in 2D, for point p(x,y)  : x^xDeg * y^yDeg
00049     in 3D, for point p(x,y,z): x^xDeg * y^yDeg * z^zDeg
00050 */
00051 double computeMonomial(FieldContainer<double> & p, int xDeg, int yDeg=0, int zDeg=0) {
00052   double val = 1.0;
00053   int polydeg[3];
00054   polydeg[0] = xDeg; polydeg[1] = yDeg; polydeg[2] = zDeg;
00055   for (int i=0; i<p.dimension(0); i++) {
00056     val *= std::pow(p(i),polydeg[i]);
00057   }
00058   return val;
00059 }
00060 
00061 
00062 /*
00063   Computes integrals of monomials over a given reference cell.
00064 */
00065 double computeIntegral(shards::CellTopology & cellTopology, int cubDegree, int xDeg, int yDeg, int zDeg) {
00066 
00067   DefaultCubatureFactory<double>  cubFactory;                                         // create factory
00068   Teuchos::RCP<Cubature<double> > myCub = cubFactory.create(cellTopology, cubDegree); // create default cubature
00069 
00070   double val       = 0.0;
00071   int cubDim       = myCub->getDimension();
00072   int numCubPoints = myCub->getNumPoints();
00073 
00074   FieldContainer<double> point(cubDim);
00075   FieldContainer<double> cubPoints(numCubPoints, cubDim);
00076   FieldContainer<double> cubWeights(numCubPoints);
00077 
00078   myCub->getCubature(cubPoints, cubWeights);
00079 
00080   for (int i=0; i<numCubPoints; i++) {
00081     for (int j=0; j<cubDim; j++) {
00082       point(j) = cubPoints(i,j);
00083     }
00084     val += computeMonomial(point, xDeg, yDeg, zDeg)*cubWeights(i);
00085   }
00086 
00087   return val;
00088 }
00089 
00090 
00091 int main(int argc, char *argv[]) {
00092 
00093   Teuchos::GlobalMPISession mpiSession(&argc, &argv);
00094 
00095   // This little trick lets us print to std::cout only if
00096   // a (dummy) command-line argument is provided.
00097   int iprint     = argc - 1;
00098   Teuchos::RCP<std::ostream> outStream;
00099   Teuchos::oblackholestream bhs; // outputs nothing
00100   if (iprint > 0)
00101     outStream = Teuchos::rcp(&std::cout, false);
00102   else
00103     outStream = Teuchos::rcp(&bhs, false);
00104 
00105   // Save the format state of the original std::cout.
00106   Teuchos::oblackholestream oldFormatState;
00107   oldFormatState.copyfmt(std::cout);
00108  
00109   *outStream \
00110   << "===============================================================================\n" \
00111   << "|                                                                             |\n" \
00112   << "|       Unit Test (CubatureDirect,CubatureTensor,DefaultCubatureFactory)      |\n" \
00113   << "|                                                                             |\n" \
00114   << "|     1) Computing integrals of monomials on reference cells in 3D            |\n" \
00115   << "|                - no BLAS, i.e. standard addition loops -                    |\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 3D (non-BLAS version)                     |\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+04 * INTREPID_TOL;
00135   int                                      maxDeg[3];
00136   int                                      maxOffset[3];
00137   int                                      numPoly[3];
00138   int                                      numAnalytic[3];
00139   // max polynomial degree tested, per cell type:
00140   maxDeg[0]                              = INTREPID_CUBATURE_TET_DEFAULT_MAX;
00141   maxDeg[1]                              = 20; // can be as large as INTREPID_CUBATURE_LINE_GAUSS_MAX, but runtime is excessive
00142   maxDeg[2]                              = INTREPID_CUBATURE_TRI_DEFAULT_MAX;
00143   // max polynomial degree recorded in analytic comparison files, per cell type:
00144   maxOffset[0]                           = INTREPID_CUBATURE_TET_DEFAULT_MAX;
00145   maxOffset[1]                           = INTREPID_CUBATURE_LINE_GAUSS_MAX;
00146   maxOffset[2]                           = INTREPID_CUBATURE_TRI_DEFAULT_MAX;
00147   for (int i=0; i<3; i++) {
00148     numPoly[i] = (maxDeg[i]+1)*(maxDeg[i]+2)*(maxDeg[i]+3)/6;
00149   }
00150   for (int i=0; i<3; i++) {
00151     numAnalytic[i] = (maxOffset[i]+1)*(maxOffset[i]+2)*(maxOffset[i]+3)/6;
00152   }
00153 
00154   // get names of files with analytic values
00155   std::string basedir = "./data";
00156   std::stringstream namestream[3];
00157   std::string filename[3];
00158   namestream[0] << basedir << "/TET_integrals" << ".dat";
00159   namestream[0] >> filename[0];
00160   namestream[1] << basedir << "/HEX_integrals" << ".dat";
00161   namestream[1] >> filename[1];
00162   namestream[2] << basedir << "/TRIPRISM_integrals" << ".dat";
00163   namestream[2] >> filename[2];
00164 
00165   // reference cells tested
00166   shards::CellTopology cellType[] = {shards::getCellTopologyData< shards::Tetrahedron<> >(),
00167                                      shards::getCellTopologyData< shards::Hexahedron<> >(),
00168                                      shards::getCellTopologyData< shards::Wedge<> >()};
00169   // format of data files with analytic values
00170   TypeOfExactData dataFormat[] = {INTREPID_UTILS_SCALAR, INTREPID_UTILS_FRACTION, INTREPID_UTILS_FRACTION};
00171 
00172   // compute and compare integrals
00173   try {
00174     for (int cellCt=0; cellCt < 3; cellCt++) {
00175       testInt.assign(numPoly[cellCt], tmparray);
00176       analyticInt.assign(numAnalytic[cellCt], tmparray);
00177       *outStream << "\nIntegrals of monomials on a reference " << cellType[cellCt].getBaseTopology()->name << ":\n";
00178       std::ifstream filecompare(&filename[cellCt][0]);
00179       // compute integrals
00180       for (int cubDeg=0; cubDeg <= maxDeg[cellCt]; cubDeg++) {
00181         polyCt = 0;
00182         testInt[cubDeg].resize((cubDeg+1)*(cubDeg+2)*(cubDeg+3)/6);
00183         for (int xDeg=0; xDeg <= cubDeg; xDeg++) {
00184           for (int yDeg=0; yDeg <= cubDeg-xDeg; yDeg++) {
00185             for (int zDeg=0; zDeg <= cubDeg-xDeg-yDeg; zDeg++) {
00186               testInt[cubDeg][polyCt] = computeIntegral(cellType[cellCt], cubDeg, xDeg, yDeg, zDeg);
00187               polyCt++; 
00188             }
00189           }
00190         }
00191       }
00192       // get analytic values
00193       if (filecompare.is_open()) {
00194         getAnalytic(analyticInt, filecompare, dataFormat[cellCt]);
00195         // close file
00196         filecompare.close();
00197       }
00198       // perform comparison
00199       for (int cubDeg=0; cubDeg <= maxDeg[cellCt]; cubDeg++) {
00200         polyCt = 0;
00201         offset = 0;
00202         int oldErrorFlag = errorFlag;
00203         for (int xDeg=0; xDeg <= cubDeg; xDeg++) {
00204           for (int yDeg=0; yDeg <= cubDeg-xDeg; yDeg++) {
00205             for (int zDeg=0; zDeg <= cubDeg-xDeg-yDeg; zDeg++) {
00206               double abstol = ( analyticInt[polyCt+offset][0] == 0.0 ? reltol : std::fabs(reltol*analyticInt[polyCt+offset][0]) );
00207               double absdiff = std::fabs(analyticInt[polyCt+offset][0] - testInt[cubDeg][polyCt]);
00208               if (absdiff > abstol) {
00209                 *outStream << "Cubature order " << std::setw(2) << std::left << cubDeg << " integrating "
00210                            << "x^" << std::setw(2) << std::left << xDeg << " * y^" << std::setw(2) << yDeg
00211                            << " * z^" << std::setw(2) << zDeg << ":" << "   "
00212                            << std::scientific << std::setprecision(16)
00213                            << testInt[cubDeg][polyCt] << "   " << analyticInt[polyCt+offset][0] << "   "
00214                            << std::setprecision(4) << absdiff << "   " << "<?" << "   " << abstol << "\n";
00215                 errorFlag++;
00216                 *outStream << std::right << std::setw(118) << "^^^^---FAILURE!\n";
00217               }
00218               polyCt++;
00219             }
00220             offset = offset + maxOffset[cellCt] - cubDeg;
00221           }
00222           offset = offset + (maxOffset[cellCt] - cubDeg)*(maxOffset[cellCt] - cubDeg + 1)/2;
00223         }
00224         *outStream << "Cubature order " << std::setw(2) << std::left << cubDeg;
00225         if (errorFlag == oldErrorFlag)
00226          *outStream << " passed.\n";
00227         else
00228          *outStream << " failed.\n";
00229       }
00230       *outStream << "\n";
00231     }  // end for cellCt
00232   }
00233   catch (std::logic_error err) {
00234     *outStream << err.what() << "\n";
00235     errorFlag = -1;
00236   };
00237 
00238 
00239   if (errorFlag != 0)
00240     std::cout << "End Result: TEST FAILED\n";
00241   else
00242     std::cout << "End Result: TEST PASSED\n";
00243 
00244   // reset format state of std::cout
00245   std::cout.copyfmt(oldFormatState);
00246 
00247   return errorFlag;
00248 }