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_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
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, shards::CellTopology & cellTopology, int cubDegree) {
00068
00069 DefaultCubatureFactory<double> cubFactory;
00070 Teuchos::RCP<Cubature<double> > myCub = cubFactory.create(cellTopology, cubDegree);
00071
00072 int cubDim = myCub->getDimension();
00073 int numCubPoints = myCub->getNumPoints();
00074 int numPolys = (cubDegree+1)*(cubDegree+2)*(cubDegree+3)/6;
00075
00076 FieldContainer<double> point(cubDim);
00077 FieldContainer<double> cubPoints(numCubPoints, cubDim);
00078 FieldContainer<double> cubWeights(numCubPoints);
00079 FieldContainer<double> functValues(numCubPoints, numPolys);
00080
00081 myCub->getCubature(cubPoints, cubWeights);
00082
00083 int polyCt = 0;
00084 for (int xDeg=0; xDeg <= cubDegree; xDeg++) {
00085 for (int yDeg=0; yDeg <= cubDegree-xDeg; yDeg++) {
00086 for (int zDeg=0; zDeg <= cubDegree-xDeg-yDeg; zDeg++) {
00087 for (int i=0; i<numCubPoints; i++) {
00088 for (int j=0; j<cubDim; j++) {
00089 point(j) = cubPoints(i,j);
00090 }
00091 functValues(i,polyCt) = computeMonomial(point, xDeg, yDeg, zDeg);
00092 }
00093 polyCt++;
00094 }
00095 }
00096 }
00097
00098 Teuchos::BLAS<int, double> myblas;
00099 int inc = 1;
00100 double alpha = 1.0;
00101 double beta = 0.0;
00102 myblas.GEMV(Teuchos::NO_TRANS, numPolys, numCubPoints, alpha, &functValues[0], numPolys,
00103 &cubWeights[0], inc, beta, &testIntFixDeg[0], inc);
00104 }
00105
00106
00107 int main(int argc, char *argv[]) {
00108
00109 Teuchos::GlobalMPISession mpiSession(&argc, &argv);
00110
00111
00112
00113 int iprint = argc - 1;
00114 Teuchos::RCP<std::ostream> outStream;
00115 Teuchos::oblackholestream bhs;
00116 if (iprint > 0)
00117 outStream = Teuchos::rcp(&std::cout, false);
00118 else
00119 outStream = Teuchos::rcp(&bhs, false);
00120
00121
00122 Teuchos::oblackholestream oldFormatState;
00123 oldFormatState.copyfmt(std::cout);
00124
00125 *outStream \
00126 << "===============================================================================\n" \
00127 << "| |\n" \
00128 << "| Unit Test (CubatureDirect,CubatureTensor,DefaultCubatureFactory) |\n" \
00129 << "| |\n" \
00130 << "| 1) Computing integrals of monomials on reference cells in 3D |\n" \
00131 << "| - using Level 2 BLAS - |\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 3D (Level 2 BLAS version) |\n"\
00141 << "===============================================================================\n";
00142
00143
00144 int errorFlag = 0;
00145 int polyCt = 0;
00146 int offset = 0;
00147 Teuchos::Array< Teuchos::Array<double> > testInt;
00148 Teuchos::Array< Teuchos::Array<double> > analyticInt;
00149 Teuchos::Array<double> tmparray(1);
00150 double reltol = 1.0e+04 * INTREPID_TOL;
00151 int maxDeg[3];
00152 int maxOffset[3];
00153 int numPoly[3];
00154 int numAnalytic[3];
00155
00156 maxDeg[0] = INTREPID_CUBATURE_TET_DEFAULT_MAX;
00157 maxDeg[1] = 20;
00158 maxDeg[2] = INTREPID_CUBATURE_TRI_DEFAULT_MAX;
00159
00160 maxOffset[0] = INTREPID_CUBATURE_TET_DEFAULT_MAX;
00161 maxOffset[1] = INTREPID_CUBATURE_LINE_GAUSS_MAX;
00162 maxOffset[2] = INTREPID_CUBATURE_TRI_DEFAULT_MAX;
00163 for (int i=0; i<3; i++) {
00164 numPoly[i] = (maxDeg[i]+1)*(maxDeg[i]+2)*(maxDeg[i]+3)/6;
00165 }
00166 for (int i=0; i<3; i++) {
00167 numAnalytic[i] = (maxOffset[i]+1)*(maxOffset[i]+2)*(maxOffset[i]+3)/6;
00168 }
00169
00170
00171
00172 std::string basedir = "./data";
00173 std::stringstream namestream[3];
00174 std::string filename[3];
00175 namestream[0] << basedir << "/TET_integrals" << ".dat";
00176 namestream[0] >> filename[0];
00177 namestream[1] << basedir << "/HEX_integrals" << ".dat";
00178 namestream[1] >> filename[1];
00179 namestream[2] << basedir << "/TRIPRISM_integrals" << ".dat";
00180 namestream[2] >> filename[2];
00181
00182
00183 shards::CellTopology cellType[] = {shards::getCellTopologyData< shards::Tetrahedron<> >(),
00184 shards::getCellTopologyData< shards::Hexahedron<> >(),
00185 shards::getCellTopologyData< shards::Wedge<> >()};
00186
00187 TypeOfExactData dataFormat[] = {INTREPID_UTILS_SCALAR, INTREPID_UTILS_FRACTION, INTREPID_UTILS_FRACTION};
00188
00189
00190 try {
00191 for (int cellCt=0; cellCt < 3; cellCt++) {
00192 testInt.assign(numPoly[cellCt], tmparray);
00193 analyticInt.assign(numAnalytic[cellCt], tmparray);
00194
00195 *outStream << "\nIntegrals of monomials on a reference " << cellType[cellCt].getBaseTopology()->name << ":\n";
00196 std::ifstream filecompare(&filename[cellCt][0]);
00197
00198 for (int cubDeg=0; cubDeg <= maxDeg[cellCt]; cubDeg++) {
00199 int numMonomials = (cubDeg+1)*(cubDeg+2)*(cubDeg+3)/6;
00200 testInt[cubDeg].resize(numMonomials);
00201 computeIntegral(testInt[cubDeg], cellType[cellCt], cubDeg);
00202 }
00203
00204 if (filecompare.is_open()) {
00205 getAnalytic(analyticInt, filecompare, dataFormat[cellCt]);
00206
00207 filecompare.close();
00208 }
00209
00210 for (int cubDeg=0; cubDeg <= maxDeg[cellCt]; cubDeg++) {
00211 polyCt = 0;
00212 offset = 0;
00213 int oldErrorFlag = errorFlag;
00214 for (int xDeg=0; xDeg <= cubDeg; xDeg++) {
00215 for (int yDeg=0; yDeg <= cubDeg-xDeg; yDeg++) {
00216 for (int zDeg=0; zDeg <= cubDeg-xDeg-yDeg; zDeg++) {
00217 double abstol = ( analyticInt[polyCt+offset][0] == 0.0 ? reltol : std::fabs(reltol*analyticInt[polyCt+offset][0]) );
00218 double absdiff = std::fabs(analyticInt[polyCt+offset][0] - testInt[cubDeg][polyCt]);
00219 if (absdiff > abstol) {
00220 *outStream << "Cubature order " << std::setw(2) << std::left << cubDeg << " integrating "
00221 << "x^" << std::setw(2) << std::left << xDeg << " * y^" << std::setw(2) << yDeg
00222 << " * z^" << std::setw(2) << zDeg << ":" << " "
00223 << std::scientific << std::setprecision(16)
00224 << testInt[cubDeg][polyCt] << " " << analyticInt[polyCt+offset][0] << " "
00225 << std::setprecision(4) << absdiff << " " << "<?" << " " << abstol << "\n";
00226 errorFlag++;
00227 *outStream << std::right << std::setw(118) << "^^^^---FAILURE!\n";
00228 }
00229 polyCt++;
00230 }
00231 offset = offset + maxOffset[cellCt] - cubDeg;
00232 }
00233 offset = offset + (maxOffset[cellCt] - cubDeg)*(maxOffset[cellCt] - cubDeg + 1)/2;
00234 }
00235 *outStream << "Cubature order " << std::setw(2) << std::left << cubDeg;
00236 if (errorFlag == oldErrorFlag)
00237 *outStream << " passed.\n";
00238 else
00239 *outStream << " failed.\n";
00240 }
00241 *outStream << "\n";
00242 }
00243 }
00244 catch (std::logic_error err) {
00245 *outStream << err.what() << "\n";
00246 errorFlag = -1;
00247 };
00248
00249
00250 if (errorFlag != 0)
00251 std::cout << "End Result: TEST FAILED\n";
00252 else
00253 std::cout << "End Result: TEST PASSED\n";
00254
00255
00256 std::cout.copyfmt(oldFormatState);
00257
00258 return errorFlag;
00259 }