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(shards::CellTopology & cellTopology,
int cubDegree,
int xDeg,
int yDeg) {
83 Teuchos::RCP<Cubature<double> > myCub = cubFactory.
create(cellTopology, cubDegree);
86 int cubDim = myCub->getDimension();
87 int numCubPoints = myCub->getNumPoints();
93 myCub->getCubature(cubPoints, cubWeights);
95 for (
int i=0; i<numCubPoints; i++) {
96 for (
int j=0; j<cubDim; j++) {
97 point(j) = cubPoints(i,j);
99 val += computeMonomial(point, xDeg, yDeg)*cubWeights(i);
106 int main(
int argc,
char *argv[]) {
108 Teuchos::GlobalMPISession mpiSession(&argc, &argv);
112 int iprint = argc - 1;
113 Teuchos::RCP<std::ostream> outStream;
114 Teuchos::oblackholestream bhs;
116 outStream = Teuchos::rcp(&std::cout,
false);
118 outStream = Teuchos::rcp(&bhs,
false);
121 Teuchos::oblackholestream oldFormatState;
122 oldFormatState.copyfmt(std::cout);
125 <<
"===============================================================================\n" \
127 <<
"| Unit Test (CubatureDirect,CubatureTensor,DefaultCubatureFactory) |\n" \
129 <<
"| 1) Computing integrals of monomials on reference cells in 2D |\n" \
131 <<
"| Questions? Contact Pavel Bochev (pbboche@sandia.gov) or |\n" \
132 <<
"| Denis Ridzal (dridzal@sandia.gov). |\n" \
134 <<
"| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \
135 <<
"| Trilinos website: http://trilinos.sandia.gov |\n" \
137 <<
"===============================================================================\n"\
138 <<
"| TEST 1: integrals of monomials in 2D |\n"\
139 <<
"===============================================================================\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+03 * INTREPID_TOL;
157 std::string basedir =
"./data";
158 std::stringstream namestream[2];
159 std::string filename[2];
160 namestream[0] << basedir <<
"/TRI_integrals" <<
".dat";
161 namestream[0] >> filename[0];
162 namestream[1] << basedir <<
"/QUAD_integrals" <<
".dat";
163 namestream[1] >> filename[1];
165 shards::CellTopology cellType[] = {shards::getCellTopologyData< shards::Triangle<> >(),
166 shards::getCellTopologyData< shards::Quadrilateral<> >()};
170 for (
int cellCt=0; cellCt < 2; cellCt++) {
171 testInt.assign(numPoly[cellCt], tmparray);
172 analyticInt.assign(numPoly[cellCt], tmparray);
173 *outStream <<
"\nIntegrals of monomials on a reference " << cellType[cellCt].getBaseCellTopologyData()->name <<
":\n";
174 std::ifstream filecompare(&filename[cellCt][0]);
176 for (
int cubDeg=0; cubDeg <= maxDeg[cellCt]; cubDeg++) {
178 testInt[cubDeg].resize((cubDeg+1)*(cubDeg+2)/2);
179 for (
int xDeg=0; xDeg <= cubDeg; xDeg++) {
180 for (
int yDeg=0; yDeg <= cubDeg-xDeg; yDeg++) {
181 testInt[cubDeg][polyCt] = computeIntegral(cellType[cellCt], cubDeg, xDeg, yDeg);
187 if (filecompare.is_open()) {
188 getAnalytic(analyticInt, filecompare);
193 for (
int cubDeg=0; cubDeg <= maxDeg[cellCt]; cubDeg++) {
196 for (
int xDeg=0; xDeg <= cubDeg; xDeg++) {
197 for (
int yDeg=0; yDeg <= cubDeg-xDeg; yDeg++) {
198 double abstol = ( analyticInt[polyCt+offset][0] == 0.0 ? reltol : std::fabs(reltol*analyticInt[polyCt+offset][0]) );
199 double absdiff = std::fabs(analyticInt[polyCt+offset][0] - testInt[cubDeg][polyCt]);
200 *outStream <<
"Cubature order " << std::setw(2) << std::left << cubDeg <<
" integrating "
201 <<
"x^" << std::setw(2) << std::left << xDeg <<
" * y^" << std::setw(2) << yDeg <<
":" <<
" "
202 << std::scientific << std::setprecision(16) << testInt[cubDeg][polyCt] <<
" " << analyticInt[polyCt+offset][0] <<
" "
203 << std::setprecision(4) << absdiff <<
" " <<
"<?" <<
" " << abstol <<
"\n";
204 if (absdiff > abstol) {
206 *outStream << std::right << std::setw(111) <<
"^^^^---FAILURE!\n";
210 offset = offset + maxDeg[cellCt] - cubDeg;
217 catch (std::logic_error err) {
218 *outStream << err.what() <<
"\n";
224 std::cout <<
"End Result: TEST FAILED\n";
226 std::cout <<
"End Result: TEST PASSED\n";
229 std::cout.copyfmt(oldFormatState);
int dimension(const int whichDim) const
Returns the specified dimension.
#define INTREPID_CUBATURE_TRI_DEFAULT_MAX
The maximum degree of the polynomial that can be integrated exactly by a direct triangle rule of the ...
#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.