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