Intrepid
/usr/src/RPM/BUILD/trilinos10-10.6.4/packages/intrepid/test/Discretization/Integration/test_10.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_CubaturePolylib.hpp"
00038 #include "Intrepid_Utils.hpp"
00039 #include "Teuchos_oblackholestream.hpp"
00040 #include "Teuchos_RCP.hpp"
00041 #include "Teuchos_GlobalMPISession.hpp"
00042 
00043 #define INTREPID_CUBATURE_LINE_MAX 61
00044 
00045 using namespace Intrepid;
00046 
00047 
00048 /*
00049   Monomial evaluation.
00050     in 1D, for point p(x)    : x^xDeg
00051     in 2D, for point p(x,y)  : x^xDeg * y^yDeg
00052     in 3D, for point p(x,y,z): x^xDeg * y^yDeg * z^zDeg
00053 */
00054 double computeMonomial(FieldContainer<double> & p, int xDeg, int yDeg=0, int zDeg=0) {
00055   double val = 1.0;
00056   int polydeg[3];
00057   polydeg[0] = xDeg; polydeg[1] = yDeg; polydeg[2] = zDeg;
00058   for (int i=0; i<p.dimension(0); i++) {
00059     val *= std::pow(p(i),polydeg[i]);
00060   }
00061   return val;
00062 }
00063 
00064 
00065 /*
00066   Computes integrals of monomials over a given reference cell.
00067 */
00068 double computeIntegral(int cubDegree, int polyDegree, EIntrepidPLPoly poly_type) {
00069 
00070   CubaturePolylib<double> lineCub(cubDegree, poly_type);
00071   double val = 0.0;
00072 
00073   int cubDim =  lineCub.getDimension();
00074 
00075   int numCubPoints = lineCub.getNumPoints();
00076 
00077   FieldContainer<double> point(cubDim);
00078   FieldContainer<double> cubPoints(numCubPoints, cubDim);
00079   FieldContainer<double> cubWeights(numCubPoints);
00080 
00081   lineCub.getCubature(cubPoints, cubWeights);
00082 
00083   for (int i=0; i<numCubPoints; i++) {
00084     for (int j=0; j<cubDim; j++) {
00085       point(j) = cubPoints(i,j);
00086     }
00087     val += computeMonomial(point, polyDegree)*cubWeights(i);
00088   }
00089 
00090   return val;
00091 }
00092 
00093 
00094 int main(int argc, char *argv[]) {
00095   
00096   Teuchos::GlobalMPISession mpiSession(&argc, &argv);
00097 
00098   // This little trick lets us print to std::cout only if
00099   // a (dummy) command-line argument is provided.
00100   int iprint     = argc - 1;
00101   Teuchos::RCP<std::ostream> outStream;
00102   Teuchos::oblackholestream bhs; // outputs nothing
00103   if (iprint > 0)
00104     outStream = Teuchos::rcp(&std::cout, false);
00105   else
00106     outStream = Teuchos::rcp(&bhs, false);
00107 
00108   // Save the format state of the original std::cout.
00109   Teuchos::oblackholestream oldFormatState;
00110   oldFormatState.copyfmt(std::cout);
00111  
00112   *outStream \
00113   << "===============================================================================\n" \
00114   << "|                                                                             |\n" \
00115   << "|                         Unit Test (CubaturePolylib)                         |\n" \
00116   << "|                                                                             |\n" \
00117   << "|     1) Computing integrals of monomials on reference cells in 1D            |\n" \
00118   << "|                                                                             |\n" \
00119   << "|  Questions? Contact  Pavel Bochev (pbboche@sandia.gov) or                   |\n" \
00120   << "|                      Denis Ridzal (dridzal@sandia.gov).                     |\n" \
00121   << "|                                                                             |\n" \
00122   << "|  Intrepid's website: http://trilinos.sandia.gov/packages/intrepid           |\n" \
00123   << "|  Trilinos website:   http://trilinos.sandia.gov                             |\n" \
00124   << "|                                                                             |\n" \
00125   << "===============================================================================\n"\
00126   << "| TEST 1: integrals of monomials in 1D                                        |\n"\
00127   << "===============================================================================\n";
00128 
00129   // internal variables:
00130   int                                      errorFlag = 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   testInt.assign(INTREPID_CUBATURE_LINE_MAX+1, tmparray);
00136   analyticInt.assign(INTREPID_CUBATURE_LINE_MAX+1, tmparray);
00137 
00138   // open file with analytic values
00139   std::string basedir = "./data";
00140   std::stringstream namestream;
00141   std::string filename;
00142   namestream <<  basedir << "/EDGE_integrals" << ".dat";
00143   namestream >> filename;
00144   std::ifstream filecompare(&filename[0]);
00145 
00146   *outStream << "\nIntegrals of monomials on a reference line (edge):\n";
00147 
00148   // compute and compare integrals
00149   try {
00150     for (EIntrepidPLPoly poly_type=PL_GAUSS; poly_type <= PL_GAUSS_LOBATTO; poly_type++) {
00151       // compute integrals
00152       for (int cubDeg=0; cubDeg <= INTREPID_CUBATURE_LINE_MAX; cubDeg++) {
00153         testInt[cubDeg].resize(cubDeg+1);
00154         for (int polyDeg=0; polyDeg <= cubDeg; polyDeg++) {
00155           testInt[cubDeg][polyDeg] = computeIntegral(cubDeg, polyDeg, poly_type);
00156         }
00157       }
00158       // get analytic values
00159       if (filecompare.is_open()) {
00160         getAnalytic(analyticInt, filecompare);
00161         // close file
00162         filecompare.close();
00163       }
00164       // perform comparison
00165       for (int cubDeg=0; cubDeg <= INTREPID_CUBATURE_LINE_MAX; cubDeg++) {
00166         for (int polyDeg=0; polyDeg <= cubDeg; polyDeg++) {
00167           double abstol = ( analyticInt[polyDeg][0] == 0.0 ? reltol : std::fabs(reltol*analyticInt[polyDeg][0]) );
00168           double absdiff = std::fabs(analyticInt[polyDeg][0] - testInt[cubDeg][polyDeg]);
00169           *outStream << "Cubature order " << std::setw(2) << std::left << cubDeg << " integrating "
00170                      << "x^" << std::setw(2) << std::left << polyDeg <<  ":" << "   "
00171                      << std::scientific << std::setprecision(16) << testInt[cubDeg][polyDeg] << "   " << analyticInt[polyDeg][0] << "   "
00172                      << std::setprecision(4) << absdiff << "   " << "<?" << "   " << abstol << "\n";
00173           if (absdiff > abstol) {
00174             errorFlag++;
00175             *outStream << std::right << std::setw(104) << "^^^^---FAILURE!\n";
00176           }
00177         }
00178         *outStream << "\n";
00179       } // end for cubDeg
00180     } // end for poly_type
00181   }
00182   catch (std::logic_error err) {
00183     *outStream << err.what() << "\n";
00184     errorFlag = -1;
00185   };
00186 
00187 
00188   if (errorFlag != 0)
00189     std::cout << "End Result: TEST FAILED\n";
00190   else
00191     std::cout << "End Result: TEST PASSED\n";
00192 
00193   // reset format state of std::cout
00194   std::cout.copyfmt(oldFormatState);
00195 
00196   return errorFlag;
00197 }