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