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_BLAS.hpp"
00042 #include "Teuchos_GlobalMPISession.hpp"
00043
00044 using namespace Intrepid;
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 void computeIntegral(Teuchos::Array<double>& testIntFixDeg, int cubDegree) {
00068
00069 CubatureSparse<double,3> myCub(cubDegree);
00070
00071 int cubDim = myCub.getDimension();
00072 int numCubPoints = myCub.getNumPoints();
00073 int numPolys = (cubDegree+1)*(cubDegree+2)*(cubDegree+3)/6;
00074
00075 FieldContainer<double> point(cubDim);
00076 FieldContainer<double> cubPoints(numCubPoints, cubDim);
00077 FieldContainer<double> cubWeights(numCubPoints);
00078 FieldContainer<double> functValues(numCubPoints, numPolys);
00079
00080 myCub.getCubature(cubPoints, cubWeights);
00081
00082 int polyCt = 0;
00083 for (int xDeg=0; xDeg <= cubDegree; xDeg++) {
00084 for (int yDeg=0; yDeg <= cubDegree-xDeg; yDeg++) {
00085 for (int zDeg=0; zDeg <= cubDegree-xDeg-yDeg; zDeg++) {
00086 for (int i=0; i<numCubPoints; i++) {
00087 for (int j=0; j<cubDim; j++) {
00088 point(j) = cubPoints(i,j);
00089 }
00090 functValues(i,polyCt) = computeMonomial(point, xDeg, yDeg, zDeg);
00091 }
00092 polyCt++;
00093 }
00094 }
00095 }
00096
00097 Teuchos::BLAS<int, double> myblas;
00098 int inc = 1;
00099 double alpha = 1.0;
00100 double beta = 0.0;
00101 myblas.GEMV(Teuchos::NO_TRANS, numPolys, numCubPoints, alpha, &functValues[0], numPolys,
00102 &cubWeights[0], inc, beta, &testIntFixDeg[0], inc);
00103 }
00104
00105
00106 int main(int argc, char *argv[]) {
00107
00108 Teuchos::GlobalMPISession mpiSession(&argc, &argv);
00109
00110
00111
00112 int iprint = argc - 1;
00113 Teuchos::RCP<std::ostream> outStream;
00114 Teuchos::oblackholestream bhs;
00115 if (iprint > 0)
00116 outStream = Teuchos::rcp(&std::cout, false);
00117 else
00118 outStream = Teuchos::rcp(&bhs, false);
00119
00120
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 3D |\n" \
00130 << "| - using Level 2 BLAS - |\n" \
00131 << "| |\n" \
00132 << "| Questions? Contact Pavel Bochev (pbboche@sandia.gov), |\n" \
00133 << "| Denis Ridzal (dridzal@sandia.gov) or |\n" \
00134 << "| Matthew Keegan (mskeega@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 3D (Level 2 BLAS version) using Sparse |\n"\
00141 << "| Grid Construction |\n"\
00142 << "===============================================================================\n";
00143
00144
00145 int errorFlag = 0;
00146 int polyCt = 0;
00147 int offset = 0;
00148 Teuchos::Array< Teuchos::Array<double> > testInt;
00149 Teuchos::Array< Teuchos::Array<double> > analyticInt;
00150 Teuchos::Array<double> tmparray(1);
00151 double reltol = 1.0e+04 * INTREPID_TOL;
00152 int maxDeg = 20;
00153 int maxOffset = INTREPID_CUBATURE_LINE_GAUSS_MAX;
00154 int numPoly = (maxDeg+1)*(maxDeg+2)*(maxDeg+3)/6;
00155 int numAnalytic = (maxOffset+1)*(maxOffset+2)*(maxOffset+3)/6;
00156 testInt.assign(numPoly, tmparray);
00157 analyticInt.assign(numAnalytic, tmparray);
00158
00159
00160 std::string basedir = "./data";
00161 std::stringstream namestream;
00162 std::string filename;
00163 namestream << basedir << "/HEX_integrals" << ".dat";
00164 namestream >> filename;
00165
00166
00167 TypeOfExactData dataFormat = INTREPID_UTILS_FRACTION;
00168
00169
00170 try {
00171 *outStream << "\nIntegrals of monomials:\n";
00172 std::ifstream filecompare(&filename[0]);
00173
00174 for (int cubDeg=0; cubDeg <= maxDeg; cubDeg++) {
00175 int numMonomials = (cubDeg+1)*(cubDeg+2)*(cubDeg+3)/6;
00176 testInt[cubDeg].resize(numMonomials);
00177 computeIntegral(testInt[cubDeg], cubDeg);
00178 }
00179
00180 if (filecompare.is_open()) {
00181 getAnalytic(analyticInt, filecompare, dataFormat);
00182
00183 filecompare.close();
00184 }
00185
00186 for (int cubDeg=0; cubDeg <= maxDeg; cubDeg++) {
00187 polyCt = 0;
00188 offset = 0;
00189 int oldErrorFlag = errorFlag;
00190 for (int xDeg=0; xDeg <= cubDeg; xDeg++) {
00191 for (int yDeg=0; yDeg <= cubDeg-xDeg; yDeg++) {
00192 for (int zDeg=0; zDeg <= cubDeg-xDeg-yDeg; zDeg++) {
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 if (absdiff > abstol) {
00196 *outStream << "Cubature order " << std::setw(2) << std::left << cubDeg << " integrating "
00197 << "x^" << std::setw(2) << std::left << xDeg << " * y^" << std::setw(2) << yDeg
00198 << " * z^" << std::setw(2) << zDeg << ":" << " "
00199 << std::scientific << std::setprecision(16)
00200 << testInt[cubDeg][polyCt] << " " << analyticInt[polyCt+offset][0] << " "
00201 << std::setprecision(4) << absdiff << " " << "<?" << " " << abstol << "\n";
00202 errorFlag++;
00203 *outStream << std::right << std::setw(118) << "^^^^---FAILURE!\n";
00204 }
00205 polyCt++;
00206 }
00207 offset = offset + maxOffset - cubDeg;
00208 }
00209 offset = offset + (maxOffset - cubDeg)*(maxOffset - cubDeg + 1)/2;
00210 }
00211 *outStream << "Cubature order " << std::setw(2) << std::left << cubDeg;
00212 if (errorFlag == oldErrorFlag)
00213 *outStream << " passed.\n";
00214 else
00215 *outStream << " failed.\n";
00216 }
00217 *outStream << "\n";
00218 }
00219 catch (std::logic_error err) {
00220 *outStream << err.what() << "\n";
00221 errorFlag = -1;
00222 };
00223
00224
00225 if (errorFlag != 0)
00226 std::cout << "End Result: TEST FAILED\n";
00227 else
00228 std::cout << "End Result: TEST PASSED\n";
00229
00230
00231 std::cout.copyfmt(oldFormatState);
00232 return errorFlag;
00233 }