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