53 #include "Teuchos_oblackholestream.hpp"
54 #include "Teuchos_RCP.hpp"
55 #include "Teuchos_GlobalMPISession.hpp"
57 using namespace Intrepid;
69 polydeg[0] = xDeg; polydeg[1] = yDeg; polydeg[2] = zDeg;
71 val *= std::pow(p(i),polydeg[i]);
80 double computeIntegral(
int cubDegree,
int polyDegree) {
83 shards::CellTopology line(shards::getCellTopologyData< shards::Line<> >());
84 Teuchos::RCP<Cubature<double> > lineCub = cubFactory.
create(line, cubDegree);
87 int cubDim = lineCub->getDimension();
89 int numCubPoints = lineCub->getNumPoints();
95 lineCub->getCubature(cubPoints, cubWeights);
97 for (
int i=0; i<numCubPoints; i++) {
98 for (
int j=0; j<cubDim; j++) {
99 point(j) = cubPoints(i,j);
101 val += computeMonomial(point, polyDegree)*cubWeights(i);
108 int main(
int argc,
char *argv[]) {
110 Teuchos::GlobalMPISession mpiSession(&argc, &argv);
114 int iprint = argc - 1;
115 Teuchos::RCP<std::ostream> outStream;
116 Teuchos::oblackholestream bhs;
118 outStream = Teuchos::rcp(&std::cout,
false);
120 outStream = Teuchos::rcp(&bhs,
false);
123 Teuchos::oblackholestream oldFormatState;
124 oldFormatState.copyfmt(std::cout);
127 <<
"===============================================================================\n" \
129 <<
"| Unit Test (CubatureDirectLineGauss) |\n" \
131 <<
"| 1) Computing integrals of monomials on reference cells in 1D |\n" \
133 <<
"| Questions? Contact Pavel Bochev (pbboche@sandia.gov) or |\n" \
134 <<
"| Denis Ridzal (dridzal@sandia.gov). |\n" \
136 <<
"| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \
137 <<
"| Trilinos website: http://trilinos.sandia.gov |\n" \
139 <<
"===============================================================================\n"\
140 <<
"| TEST 1: integrals of monomials in 1D |\n"\
141 <<
"===============================================================================\n";
145 Teuchos::Array< Teuchos::Array<double> > testInt;
146 Teuchos::Array< Teuchos::Array<double> > analyticInt;
147 Teuchos::Array<double> tmparray(1);
148 double reltol = 1.0e+01 * INTREPID_TOL;
153 std::string basedir =
"./data";
154 std::stringstream namestream;
155 std::string filename;
156 namestream << basedir <<
"/EDGE_integrals" <<
".dat";
157 namestream >> filename;
158 std::ifstream filecompare(&filename[0]);
160 *outStream <<
"\nIntegrals of monomials on a reference line (edge):\n";
166 testInt[cubDeg].resize(cubDeg+1);
167 for (
int polyDeg=0; polyDeg <= cubDeg; polyDeg++) {
168 testInt[cubDeg][polyDeg] = computeIntegral(cubDeg, polyDeg);
172 if (filecompare.is_open()) {
173 getAnalytic(analyticInt, filecompare);
179 for (
int polyDeg=0; polyDeg <= cubDeg; polyDeg++) {
180 double abstol = ( analyticInt[polyDeg][0] == 0.0 ? reltol : std::fabs(reltol*analyticInt[polyDeg][0]) );
181 double absdiff = std::fabs(analyticInt[polyDeg][0] - testInt[cubDeg][polyDeg]);
182 *outStream <<
"Cubature order " << std::setw(2) << std::left << cubDeg <<
" integrating "
183 <<
"x^" << std::setw(2) << std::left << polyDeg <<
":" <<
" "
184 << std::scientific << std::setprecision(16) << testInt[cubDeg][polyDeg] <<
" " << analyticInt[polyDeg][0] <<
" "
185 << std::setprecision(4) << absdiff <<
" " <<
"<?" <<
" " << abstol <<
"\n";
186 if (absdiff > abstol) {
188 *outStream << std::right << std::setw(104) <<
"^^^^---FAILURE!\n";
194 catch (std::logic_error err) {
195 *outStream << err.what() <<
"\n";
201 std::cout <<
"End Result: TEST FAILED\n";
203 std::cout <<
"End Result: TEST PASSED\n";
206 std::cout.copyfmt(oldFormatState);
int dimension(const int whichDim) const
Returns the specified dimension.
#define INTREPID_CUBATURE_LINE_GAUSS_MAX
The maximum degree of the polynomial that can be integrated exactly by a direct line rule of the Gaus...
Header file for the abstract base class Intrepid::DefaultCubatureFactory.
A factory class that generates specific instances of cubatures.
Teuchos::RCP< Cubature< Scalar, ArrayPoint, ArrayWeight > > create(const shards::CellTopology &cellTopology, const std::vector< int > °ree)
Factory method.