|
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 // 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_CubatureSparse.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 //Teuchos::RCP<std::ostream> outStream; 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 double computeIntegral(int cubDegree, int xDeg, int yDeg) { 00082 00083 CubatureSparse<double,2> myCub(cubDegree); 00084 00085 double val = 0.0; 00086 int cubDim = myCub.getDimension(); 00087 int numCubPoints = myCub.getNumPoints(); 00088 00089 FieldContainer<double> point(cubDim); 00090 FieldContainer<double> cubPoints(numCubPoints, cubDim); 00091 FieldContainer<double> cubWeights(numCubPoints); 00092 00093 myCub.getCubature(cubPoints, cubWeights); 00094 00095 for (int i=0; i<numCubPoints; i++) { 00096 for (int j=0; j<cubDim; j++) { 00097 point(j) = cubPoints(i,j); 00098 } 00099 val += computeMonomial(point, xDeg, yDeg)*cubWeights(i); 00100 } 00101 00102 return val; 00103 } 00104 00105 00106 int main(int argc, char *argv[]) { 00107 00108 Teuchos::GlobalMPISession mpiSession(&argc, &argv); 00109 00110 // This little trick lets us print to std::cout only if 00111 // a (dummy) command-line argument is provided. 00112 int iprint = argc - 1; 00113 Teuchos::RCP<std::ostream> outStream; 00114 Teuchos::oblackholestream bhs; // outputs nothing 00115 if (iprint > 0) 00116 outStream = Teuchos::rcp(&std::cout, false); 00117 else 00118 outStream = Teuchos::rcp(&bhs, false); 00119 00120 // Save the format state of the original std::cout. 00121 Teuchos::oblackholestream oldFormatState; 00122 oldFormatState.copyfmt(std::cout); 00123 00124 *outStream \ 00125 << "===============================================================================\n" \ 00126 << "| |\n" \ 00127 << "| Unit Test (CubatureSparse) |\n" \ 00128 << "| |\n" \ 00129 << "| 1) Computing integrals of monomials on reference cells in 2D |\n" \ 00130 << "| |\n" \ 00131 << "| Questions? Contact Pavel Bochev (pbboche@sandia.gov), |\n" \ 00132 << "| Denis Ridzal (dridzal@sandia.gov) or |\n" \ 00133 << "| Matthew Keegan (mskeega@sandia.gov) |\n" \ 00134 << "| |\n" \ 00135 << "| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \ 00136 << "| Trilinos website: http://trilinos.sandia.gov |\n" \ 00137 << "| |\n" \ 00138 << "===============================================================================\n"\ 00139 << "| TEST 1: integrals of monomials in 2D for Sparse Grid Construction |\n"\ 00140 << "===============================================================================\n"; 00141 00142 // internal variables: 00143 int errorFlag = 0; 00144 int polyCt = 0; 00145 int offset = 0; 00146 Teuchos::Array< Teuchos::Array<double> > testInt; 00147 Teuchos::Array< Teuchos::Array<double> > analyticInt; 00148 Teuchos::Array<double> tmparray(1); 00149 double reltol = 1.0e+03 * INTREPID_TOL; 00150 int maxDeg = 30; // can be as large as INTREPID_CUBATURE_SPARSE2D_GAUSS_MAX, but runtime is excessive 00151 int maxOffset = INTREPID_CUBATURE_LINE_GAUSS_MAX; 00152 int numPoly = (maxDeg+1)*(maxDeg+2)/2; 00153 int numAnalytic = (maxOffset+1)*(maxOffset+2)/2; 00154 testInt.assign(numPoly, tmparray); 00155 analyticInt.assign(numAnalytic, tmparray); 00156 00157 // get names of files with analytic values 00158 std::string basedir = "./data"; 00159 std::stringstream namestream; 00160 std::string filename; 00161 namestream << basedir << "/QUAD_integrals" << ".dat"; 00162 namestream >> filename; 00163 00164 // compute and compare integrals 00165 try { 00166 *outStream << "\nIntegrals of monomials:\n"; 00167 00168 std::ifstream filecompare(&filename[0]); 00169 // compute integrals 00170 for (int cubDeg=0; cubDeg <= maxDeg; cubDeg++) { 00171 polyCt = 0; 00172 testInt[cubDeg].resize((cubDeg+1)*(cubDeg+2)/2); 00173 for (int xDeg=0; xDeg <= cubDeg; xDeg++) { 00174 for (int yDeg=0; yDeg <= cubDeg-xDeg; yDeg++) { 00175 testInt[cubDeg][polyCt] = computeIntegral(cubDeg, xDeg, yDeg); 00176 polyCt++; 00177 } 00178 } 00179 } 00180 00181 // get analytic values 00182 if (filecompare.is_open()) { 00183 getAnalytic(analyticInt, filecompare); 00184 // close file 00185 filecompare.close(); 00186 } 00187 // perform comparison 00188 for (int cubDeg=0; cubDeg <= maxDeg; cubDeg++) { 00189 polyCt = 0; 00190 offset = 0; 00191 for (int xDeg=0; xDeg <= cubDeg; xDeg++) { 00192 for (int yDeg=0; yDeg <= cubDeg-xDeg; yDeg++) { 00193 double abstol = ( analyticInt[polyCt+offset][0] == 0.0 ? reltol : std::fabs(reltol*analyticInt[polyCt+offset][0]) ); 00194 double absdiff = std::fabs(analyticInt[polyCt+offset][0] - testInt[cubDeg][polyCt]); 00195 *outStream << "Cubature order " << std::setw(2) << std::left << cubDeg << " integrating " 00196 << "x^" << std::setw(2) << std::left << xDeg << " * y^" << std::setw(2) << yDeg << ":" << " " 00197 << std::scientific << std::setprecision(16) << testInt[cubDeg][polyCt] << " " << analyticInt[polyCt+offset][0] << " " 00198 << std::setprecision(4) << absdiff << " " << "<?" << " " << abstol << "\n"; 00199 if (absdiff > abstol) { 00200 errorFlag++; 00201 *outStream << std::right << std::setw(111) << "^^^^---FAILURE!\n"; 00202 } 00203 polyCt++; 00204 } 00205 offset = offset + maxOffset - cubDeg; 00206 } 00207 *outStream << "\n"; 00208 } 00209 *outStream << "\n"; 00210 } 00211 catch (std::logic_error err) { 00212 *outStream << err.what() << "\n"; 00213 errorFlag = -1; 00214 }; 00215 00216 00217 if (errorFlag != 0) 00218 std::cout << "End Result: TEST FAILED\n"; 00219 else 00220 std::cout << "End Result: TEST PASSED\n"; 00221 00222 // reset format state of std::cout 00223 std::cout.copyfmt(oldFormatState); 00224 00225 return errorFlag; 00226 }
1.7.6.1