53 #include "Teuchos_oblackholestream.hpp"
54 #include "Teuchos_RCP.hpp"
55 #include "Teuchos_BLAS.hpp"
56 #include "Teuchos_GlobalMPISession.hpp"
58 using namespace Intrepid;
70 polydeg[0] = xDeg; polydeg[1] = yDeg; polydeg[2] = zDeg;
72 val *= std::pow(p(i),polydeg[i]);
81 void computeIntegral(Teuchos::Array<double>& testIntFixDeg, shards::CellTopology & cellTopology,
int cubDegree) {
84 Teuchos::RCP<Cubature<double> > myCub = cubFactory.
create(cellTopology, cubDegree);
86 int cubDim = myCub->getDimension();
87 int numCubPoints = myCub->getNumPoints();
88 int numPolys = (cubDegree+1)*(cubDegree+2)*(cubDegree+3)/6;
95 myCub->getCubature(cubPoints, cubWeights);
98 for (
int xDeg=0; xDeg <= cubDegree; xDeg++) {
99 for (
int yDeg=0; yDeg <= cubDegree-xDeg; yDeg++) {
100 for (
int zDeg=0; zDeg <= cubDegree-xDeg-yDeg; zDeg++) {
101 for (
int i=0; i<numCubPoints; i++) {
102 for (
int j=0; j<cubDim; j++) {
103 point(j) = cubPoints(i,j);
105 functValues(i,polyCt) = computeMonomial(point, xDeg, yDeg, zDeg);
112 Teuchos::BLAS<int, double> myblas;
116 myblas.GEMV(Teuchos::NO_TRANS, numPolys, numCubPoints, alpha, &functValues[0], numPolys,
117 &cubWeights[0], inc, beta, &testIntFixDeg[0], inc);
121 int main(
int argc,
char *argv[]) {
123 Teuchos::GlobalMPISession mpiSession(&argc, &argv);
127 int iprint = argc - 1;
128 Teuchos::RCP<std::ostream> outStream;
129 Teuchos::oblackholestream bhs;
131 outStream = Teuchos::rcp(&std::cout,
false);
133 outStream = Teuchos::rcp(&bhs,
false);
136 Teuchos::oblackholestream oldFormatState;
137 oldFormatState.copyfmt(std::cout);
140 <<
"===============================================================================\n" \
142 <<
"| Unit Test (CubatureDirect,CubatureTensor,DefaultCubatureFactory) |\n" \
144 <<
"| 1) Computing integrals of monomials on reference cells in 3D |\n" \
145 <<
"| - using Level 2 BLAS - |\n" \
147 <<
"| Questions? Contact Pavel Bochev (pbboche@sandia.gov) or |\n" \
148 <<
"| Denis Ridzal (dridzal@sandia.gov). |\n" \
150 <<
"| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \
151 <<
"| Trilinos website: http://trilinos.sandia.gov |\n" \
153 <<
"===============================================================================\n"\
154 <<
"| TEST 1: integrals of monomials in 3D (Level 2 BLAS version) |\n"\
155 <<
"===============================================================================\n";
161 Teuchos::Array< Teuchos::Array<double> > testInt;
162 Teuchos::Array< Teuchos::Array<double> > analyticInt;
163 Teuchos::Array<double> tmparray(1);
164 double reltol = 1.0e+04 * INTREPID_TOL;
179 for (
int i=0; i<4; i++) {
180 numPoly[i] = (maxDeg[i]+1)*(maxDeg[i]+2)*(maxDeg[i]+3)/6;
182 for (
int i=0; i<4; i++) {
183 numAnalytic[i] = (maxOffset[i]+1)*(maxOffset[i]+2)*(maxOffset[i]+3)/6;
188 std::string basedir =
"./data";
189 std::stringstream namestream[4];
190 std::string filename[4];
191 namestream[0] << basedir <<
"/TET_integrals" <<
".dat";
192 namestream[0] >> filename[0];
193 namestream[1] << basedir <<
"/HEX_integrals" <<
".dat";
194 namestream[1] >> filename[1];
195 namestream[2] << basedir <<
"/TRIPRISM_integrals" <<
".dat";
196 namestream[2] >> filename[2];
197 namestream[3] << basedir <<
"/PYR_integrals" <<
".dat";
198 namestream[3] >> filename[3];
201 shards::CellTopology cellType[] = {shards::getCellTopologyData< shards::Tetrahedron<> >(),
202 shards::getCellTopologyData< shards::Hexahedron<> >(),
203 shards::getCellTopologyData< shards::Wedge<> >(),
204 shards::getCellTopologyData< shards::Pyramid<> >() };
206 TypeOfExactData dataFormat[] = {INTREPID_UTILS_SCALAR, INTREPID_UTILS_FRACTION, INTREPID_UTILS_FRACTION, INTREPID_UTILS_FRACTION};
210 for (
int cellCt=0; cellCt < 4; cellCt++) {
211 testInt.assign(numPoly[cellCt], tmparray);
212 analyticInt.assign(numAnalytic[cellCt], tmparray);
214 *outStream <<
"\nIntegrals of monomials on a reference " << cellType[cellCt].getBaseCellTopologyData()->name <<
":\n";
215 std::ifstream filecompare(&filename[cellCt][0]);
217 for (
int cubDeg=0; cubDeg <= maxDeg[cellCt]; cubDeg++) {
218 int numMonomials = (cubDeg+1)*(cubDeg+2)*(cubDeg+3)/6;
219 testInt[cubDeg].resize(numMonomials);
220 computeIntegral(testInt[cubDeg], cellType[cellCt], cubDeg);
223 if (filecompare.is_open()) {
224 getAnalytic(analyticInt, filecompare, dataFormat[cellCt]);
229 for (
int cubDeg=0; cubDeg <= maxDeg[cellCt]; cubDeg++) {
232 int oldErrorFlag = errorFlag;
233 for (
int xDeg=0; xDeg <= cubDeg; xDeg++) {
234 for (
int yDeg=0; yDeg <= cubDeg-xDeg; yDeg++) {
235 for (
int zDeg=0; zDeg <= cubDeg-xDeg-yDeg; zDeg++) {
236 double abstol = ( analyticInt[polyCt+offset][0] == 0.0 ? reltol : std::fabs(reltol*analyticInt[polyCt+offset][0]) );
237 double absdiff = std::fabs(analyticInt[polyCt+offset][0] - testInt[cubDeg][polyCt]);
238 if (absdiff > abstol) {
239 *outStream <<
"Cubature order " << std::setw(2) << std::left << cubDeg <<
" integrating "
240 <<
"x^" << std::setw(2) << std::left << xDeg <<
" * y^" << std::setw(2) << yDeg
241 <<
" * z^" << std::setw(2) << zDeg <<
":" <<
" "
242 << std::scientific << std::setprecision(16)
243 << testInt[cubDeg][polyCt] <<
" " << analyticInt[polyCt+offset][0] <<
" "
244 << std::setprecision(4) << absdiff <<
" " <<
"<?" <<
" " << abstol <<
"\n";
246 *outStream << std::right << std::setw(118) <<
"^^^^---FAILURE!\n";
250 offset = offset + maxOffset[cellCt] - cubDeg;
252 offset = offset + (maxOffset[cellCt] - cubDeg)*(maxOffset[cellCt] - cubDeg + 1)/2;
254 *outStream <<
"Cubature order " << std::setw(2) << std::left << cubDeg;
255 if (errorFlag == oldErrorFlag)
256 *outStream <<
" passed.\n";
258 *outStream <<
" failed.\n";
263 catch (std::logic_error err) {
264 *outStream << err.what() <<
"\n";
270 std::cout <<
"End Result: TEST FAILED\n";
272 std::cout <<
"End Result: TEST PASSED\n";
275 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...
#define INTREPID_CUBATURE_LINE_GAUSSJACOBI20_MAX
The maximum degree of the polynomial that can be integrated exactly by a direct line rule of the Gaus...
#define INTREPID_CUBATURE_TET_DEFAULT_MAX
The maximum degree of the polynomial that can be integrated exactly by a direct tetrahedron rule of t...
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.