00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00037 #include "Intrepid_CubatureSparse.hpp"
00038 #include "Intrepid_Utils.hpp"
00039 #include "Teuchos_oblackholestream.hpp"
00040 #include "Teuchos_RCP.hpp"
00041 #include "Teuchos_GlobalMPISession.hpp"
00042
00043 using namespace Intrepid;
00044
00045
00046
00047
00048
00049
00050
00051
00052
00053 double computeMonomial(FieldContainer<double> & p, int xDeg, int yDeg=0, int zDeg=0) {
00054 double val = 1.0;
00055 int polydeg[3];
00056 polydeg[0] = xDeg; polydeg[1] = yDeg; polydeg[2] = zDeg;
00057 for (int i=0; i<p.dimension(0); i++) {
00058 val *= std::pow(p(i),polydeg[i]);
00059 }
00060 return val;
00061 }
00062
00063
00064
00065
00066
00067 double computeIntegral(int cubDegree, int xDeg, int yDeg) {
00068
00069 CubatureSparse<double,2> myCub(cubDegree);
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
00079 myCub.getCubature(cubPoints, cubWeights);
00080
00081 for (int i=0; i<numCubPoints; i++) {
00082 for (int j=0; j<cubDim; j++) {
00083 point(j) = cubPoints(i,j);
00084 }
00085 val += computeMonomial(point, xDeg, yDeg)*cubWeights(i);
00086 }
00087
00088 return val;
00089 }
00090
00091
00092 int main(int argc, char *argv[]) {
00093
00094 Teuchos::GlobalMPISession mpiSession(&argc, &argv);
00095
00096
00097
00098 int iprint = argc - 1;
00099 Teuchos::RCP<std::ostream> outStream;
00100 Teuchos::oblackholestream bhs;
00101 if (iprint > 0)
00102 outStream = Teuchos::rcp(&std::cout, false);
00103 else
00104 outStream = Teuchos::rcp(&bhs, false);
00105
00106
00107 Teuchos::oblackholestream oldFormatState;
00108 oldFormatState.copyfmt(std::cout);
00109
00110 *outStream \
00111 << "===============================================================================\n" \
00112 << "| |\n" \
00113 << "| Unit Test (CubatureSparse) |\n" \
00114 << "| |\n" \
00115 << "| 1) Computing integrals of monomials on reference cells in 2D |\n" \
00116 << "| |\n" \
00117 << "| Questions? Contact Pavel Bochev (pbboche@sandia.gov), |\n" \
00118 << "| Denis Ridzal (dridzal@sandia.gov) or |\n" \
00119 << "| Matthew Keegan (mskeega@sandia.gov) |\n" \
00120 << "| |\n" \
00121 << "| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \
00122 << "| Trilinos website: http://trilinos.sandia.gov |\n" \
00123 << "| |\n" \
00124 << "===============================================================================\n"\
00125 << "| TEST 1: integrals of monomials in 2D for Sparse Grid Construction |\n"\
00126 << "===============================================================================\n";
00127
00128
00129 int errorFlag = 0;
00130 int polyCt = 0;
00131 int offset = 0;
00132 Teuchos::Array< Teuchos::Array<double> > testInt;
00133 Teuchos::Array< Teuchos::Array<double> > analyticInt;
00134 Teuchos::Array<double> tmparray(1);
00135 double reltol = 1.0e+03 * INTREPID_TOL;
00136 int maxDeg = 30;
00137 int maxOffset = INTREPID_CUBATURE_LINE_GAUSS_MAX;
00138 int numPoly = (maxDeg+1)*(maxDeg+2)/2;
00139 int numAnalytic = (maxOffset+1)*(maxOffset+2)/2;
00140 testInt.assign(numPoly, tmparray);
00141 analyticInt.assign(numAnalytic, tmparray);
00142
00143
00144 std::string basedir = "./data";
00145 std::stringstream namestream;
00146 std::string filename;
00147 namestream << basedir << "/QUAD_integrals" << ".dat";
00148 namestream >> filename;
00149
00150
00151 try {
00152 *outStream << "\nIntegrals of monomials:\n";
00153
00154 std::ifstream filecompare(&filename[0]);
00155
00156 for (int cubDeg=0; cubDeg <= maxDeg; cubDeg++) {
00157 polyCt = 0;
00158 testInt[cubDeg].resize((cubDeg+1)*(cubDeg+2)/2);
00159 for (int xDeg=0; xDeg <= cubDeg; xDeg++) {
00160 for (int yDeg=0; yDeg <= cubDeg-xDeg; yDeg++) {
00161 testInt[cubDeg][polyCt] = computeIntegral(cubDeg, xDeg, yDeg);
00162 polyCt++;
00163 }
00164 }
00165 }
00166
00167
00168 if (filecompare.is_open()) {
00169 getAnalytic(analyticInt, filecompare);
00170
00171 filecompare.close();
00172 }
00173
00174 for (int cubDeg=0; cubDeg <= maxDeg; cubDeg++) {
00175 polyCt = 0;
00176 offset = 0;
00177 for (int xDeg=0; xDeg <= cubDeg; xDeg++) {
00178 for (int yDeg=0; yDeg <= cubDeg-xDeg; yDeg++) {
00179 double abstol = ( analyticInt[polyCt+offset][0] == 0.0 ? reltol : std::fabs(reltol*analyticInt[polyCt+offset][0]) );
00180 double absdiff = std::fabs(analyticInt[polyCt+offset][0] - testInt[cubDeg][polyCt]);
00181 *outStream << "Cubature order " << std::setw(2) << std::left << cubDeg << " integrating "
00182 << "x^" << std::setw(2) << std::left << xDeg << " * y^" << std::setw(2) << yDeg << ":" << " "
00183 << std::scientific << std::setprecision(16) << testInt[cubDeg][polyCt] << " " << analyticInt[polyCt+offset][0] << " "
00184 << std::setprecision(4) << absdiff << " " << "<?" << " " << abstol << "\n";
00185 if (absdiff > abstol) {
00186 errorFlag++;
00187 *outStream << std::right << std::setw(111) << "^^^^---FAILURE!\n";
00188 }
00189 polyCt++;
00190 }
00191 offset = offset + maxOffset - cubDeg;
00192 }
00193 *outStream << "\n";
00194 }
00195 *outStream << "\n";
00196 }
00197 catch (std::logic_error err) {
00198 *outStream << err.what() << "\n";
00199 errorFlag = -1;
00200 };
00201
00202
00203 if (errorFlag != 0)
00204 std::cout << "End Result: TEST FAILED\n";
00205 else
00206 std::cout << "End Result: TEST PASSED\n";
00207
00208
00209 std::cout.copyfmt(oldFormatState);
00210
00211 return errorFlag;
00212 }