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