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_GlobalMPISession.hpp"
00042
00043 using namespace Intrepid;
00044
00045
00046
00047
00048
00049
00050
00051
00052 double computeMonomial(FieldContainer<double> & p, int xDeg, int yDeg=0, int zDeg=0) {
00053 double val = 1.0;
00054 int polydeg[3];
00055 polydeg[0] = xDeg; polydeg[1] = yDeg; polydeg[2] = zDeg;
00056 for (int i=0; i<p.dimension(0); i++) {
00057 val *= std::pow(p(i),polydeg[i]);
00058 }
00059 return val;
00060 }
00061
00062
00063
00064
00065
00066 double computeIntegral(int cubDegree, int polyDegree) {
00067
00068 DefaultCubatureFactory<double> cubFactory;
00069 shards::CellTopology line(shards::getCellTopologyData< shards::Line<> >());
00070 Teuchos::RCP<Cubature<double> > lineCub = cubFactory.create(line, cubDegree);
00071 double val = 0.0;
00072
00073 int cubDim = lineCub->getDimension();
00074
00075 int numCubPoints = lineCub->getNumPoints();
00076
00077 FieldContainer<double> point(cubDim);
00078 FieldContainer<double> cubPoints(numCubPoints, cubDim);
00079 FieldContainer<double> cubWeights(numCubPoints);
00080
00081 lineCub->getCubature(cubPoints, cubWeights);
00082
00083 for (int i=0; i<numCubPoints; i++) {
00084 for (int j=0; j<cubDim; j++) {
00085 point(j) = cubPoints(i,j);
00086 }
00087 val += computeMonomial(point, polyDegree)*cubWeights(i);
00088 }
00089
00090 return val;
00091 }
00092
00093
00094 int main(int argc, char *argv[]) {
00095
00096 Teuchos::GlobalMPISession mpiSession(&argc, &argv);
00097
00098
00099
00100 int iprint = argc - 1;
00101 Teuchos::RCP<std::ostream> outStream;
00102 Teuchos::oblackholestream bhs;
00103 if (iprint > 0)
00104 outStream = Teuchos::rcp(&std::cout, false);
00105 else
00106 outStream = Teuchos::rcp(&bhs, false);
00107
00108
00109 Teuchos::oblackholestream oldFormatState;
00110 oldFormatState.copyfmt(std::cout);
00111
00112 *outStream \
00113 << "===============================================================================\n" \
00114 << "| |\n" \
00115 << "| Unit Test (CubatureDirectLineGauss) |\n" \
00116 << "| |\n" \
00117 << "| 1) Computing integrals of monomials on reference cells in 1D |\n" \
00118 << "| |\n" \
00119 << "| Questions? Contact Pavel Bochev (pbboche@sandia.gov) or |\n" \
00120 << "| Denis Ridzal (dridzal@sandia.gov). |\n" \
00121 << "| |\n" \
00122 << "| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \
00123 << "| Trilinos website: http://trilinos.sandia.gov |\n" \
00124 << "| |\n" \
00125 << "===============================================================================\n"\
00126 << "| TEST 1: integrals of monomials in 1D |\n"\
00127 << "===============================================================================\n";
00128
00129
00130 int errorFlag = 0;
00131 Teuchos::Array< Teuchos::Array<double> > testInt;
00132 Teuchos::Array< Teuchos::Array<double> > analyticInt;
00133 Teuchos::Array<double> tmparray(1);
00134 double reltol = 1.0e+01 * INTREPID_TOL;
00135 testInt.assign(INTREPID_CUBATURE_LINE_GAUSS_MAX+1, tmparray);
00136 analyticInt.assign(INTREPID_CUBATURE_LINE_GAUSS_MAX+1, tmparray);
00137
00138
00139 std::string basedir = "./data";
00140 std::stringstream namestream;
00141 std::string filename;
00142 namestream << basedir << "/EDGE_integrals" << ".dat";
00143 namestream >> filename;
00144 std::ifstream filecompare(&filename[0]);
00145
00146 *outStream << "\nIntegrals of monomials on a reference line (edge):\n";
00147
00148
00149 try {
00150
00151 for (int cubDeg=0; cubDeg <= INTREPID_CUBATURE_LINE_GAUSS_MAX; cubDeg++) {
00152 testInt[cubDeg].resize(cubDeg+1);
00153 for (int polyDeg=0; polyDeg <= cubDeg; polyDeg++) {
00154 testInt[cubDeg][polyDeg] = computeIntegral(cubDeg, polyDeg);
00155 }
00156 }
00157
00158 if (filecompare.is_open()) {
00159 getAnalytic(analyticInt, filecompare);
00160
00161 filecompare.close();
00162 }
00163
00164 for (int cubDeg=0; cubDeg <= INTREPID_CUBATURE_LINE_GAUSS_MAX; cubDeg++) {
00165 for (int polyDeg=0; polyDeg <= cubDeg; polyDeg++) {
00166 double abstol = ( analyticInt[polyDeg][0] == 0.0 ? reltol : std::fabs(reltol*analyticInt[polyDeg][0]) );
00167 double absdiff = std::fabs(analyticInt[polyDeg][0] - testInt[cubDeg][polyDeg]);
00168 *outStream << "Cubature order " << std::setw(2) << std::left << cubDeg << " integrating "
00169 << "x^" << std::setw(2) << std::left << polyDeg << ":" << " "
00170 << std::scientific << std::setprecision(16) << testInt[cubDeg][polyDeg] << " " << analyticInt[polyDeg][0] << " "
00171 << std::setprecision(4) << absdiff << " " << "<?" << " " << abstol << "\n";
00172 if (absdiff > abstol) {
00173 errorFlag++;
00174 *outStream << std::right << std::setw(104) << "^^^^---FAILURE!\n";
00175 }
00176 }
00177 *outStream << "\n";
00178 }
00179 }
00180 catch (std::logic_error err) {
00181 *outStream << err.what() << "\n";
00182 errorFlag = -1;
00183 };
00184
00185
00186 if (errorFlag != 0)
00187 std::cout << "End Result: TEST FAILED\n";
00188 else
00189 std::cout << "End Result: TEST PASSED\n";
00190
00191
00192 std::cout.copyfmt(oldFormatState);
00193
00194 return errorFlag;
00195 }