53 #include "Teuchos_oblackholestream.hpp"
54 #include "Teuchos_RCP.hpp"
55 #include "Teuchos_BLAS.hpp"
56 #include "Teuchos_GlobalMPISession.hpp"
58 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,
int zDeg) {
83 Teuchos::RCP<Cubature<double> > myCub = cubFactory.
create(cellTopology, cubDegree);
86 int cubDim = myCub->getDimension();
87 int numCubPoints = myCub->getNumPoints();
94 myCub->getCubature(cubPoints, cubWeights);
96 for (
int i=0; i<numCubPoints; i++) {
97 for (
int j=0; j<cubDim; j++) {
98 point(j) = cubPoints(i,j);
100 functValues(i) = computeMonomial(point, xDeg, yDeg, zDeg);
103 Teuchos::BLAS<int, double> myblas;
105 val = myblas.DOT(numCubPoints, &functValues[0], inc, &cubWeights[0], inc);
112 int main(
int argc,
char *argv[]) {
114 Teuchos::GlobalMPISession mpiSession(&argc, &argv);
118 int iprint = argc - 1;
119 Teuchos::RCP<std::ostream> outStream;
120 Teuchos::oblackholestream bhs;
122 outStream = Teuchos::rcp(&std::cout,
false);
124 outStream = Teuchos::rcp(&bhs,
false);
127 Teuchos::oblackholestream oldFormatState;
128 oldFormatState.copyfmt(std::cout);
131 <<
"===============================================================================\n" \
133 <<
"| Unit Test (CubatureDirect,CubatureTensor,DefaultCubatureFactory) |\n" \
135 <<
"| 1) Computing integrals of monomials on reference cells in 3D |\n" \
136 <<
"| - using Level 1 BLAS - |\n" \
138 <<
"| Questions? Contact Pavel Bochev (pbboche@sandia.gov) or |\n" \
139 <<
"| Denis Ridzal (dridzal@sandia.gov). |\n" \
141 <<
"| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \
142 <<
"| Trilinos website: http://trilinos.sandia.gov |\n" \
144 <<
"===============================================================================\n"\
145 <<
"| TEST 1: integrals of monomials in 3D (Level 1 BLAS version) |\n"\
146 <<
"===============================================================================\n";
152 Teuchos::Array< Teuchos::Array<double> > testInt;
153 Teuchos::Array< Teuchos::Array<double> > analyticInt;
154 Teuchos::Array<double> tmparray(1);
155 double reltol = 1.0e+04 * INTREPID_TOL;
170 for (
int i=0; i<4; i++) {
171 numPoly[i] = (maxDeg[i]+1)*(maxDeg[i]+2)*(maxDeg[i]+3)/6;
173 for (
int i=0; i<4; i++) {
174 numAnalytic[i] = (maxOffset[i]+1)*(maxOffset[i]+2)*(maxOffset[i]+3)/6;
178 std::string basedir =
"./data";
179 std::stringstream namestream[4];
180 std::string filename[4];
181 namestream[0] << basedir <<
"/TET_integrals" <<
".dat";
182 namestream[0] >> filename[0];
183 namestream[1] << basedir <<
"/HEX_integrals" <<
".dat";
184 namestream[1] >> filename[1];
185 namestream[2] << basedir <<
"/TRIPRISM_integrals" <<
".dat";
186 namestream[2] >> filename[2];
187 namestream[3] << basedir <<
"/PYR_integrals" <<
".dat";
188 namestream[3] >> filename[3];
191 shards::CellTopology cellType[] = {shards::getCellTopologyData< shards::Tetrahedron<> >(),
192 shards::getCellTopologyData< shards::Hexahedron<> >(),
193 shards::getCellTopologyData< shards::Wedge<> >(),
194 shards::getCellTopologyData< shards::Pyramid<> >() };
196 TypeOfExactData dataFormat[] = {INTREPID_UTILS_SCALAR, INTREPID_UTILS_FRACTION, INTREPID_UTILS_FRACTION, INTREPID_UTILS_FRACTION};
200 for (
int cellCt=0; cellCt < 4; cellCt++) {
201 testInt.assign(numPoly[cellCt], tmparray);
202 analyticInt.assign(numAnalytic[cellCt], tmparray);
204 *outStream <<
"\nIntegrals of monomials on a reference " << cellType[cellCt].getBaseCellTopologyData()->name <<
":\n";
205 std::ifstream filecompare(&filename[cellCt][0]);
207 for (
int cubDeg=0; cubDeg <= maxDeg[cellCt]; cubDeg++) {
209 testInt[cubDeg].resize((cubDeg+1)*(cubDeg+2)*(cubDeg+3)/6);
210 for (
int xDeg=0; xDeg <= cubDeg; xDeg++) {
211 for (
int yDeg=0; yDeg <= cubDeg-xDeg; yDeg++) {
212 for (
int zDeg=0; zDeg <= cubDeg-xDeg-yDeg; zDeg++) {
213 testInt[cubDeg][polyCt] = computeIntegral(cellType[cellCt], cubDeg, xDeg, yDeg, zDeg);
220 if (filecompare.is_open()) {
221 getAnalytic(analyticInt, filecompare, dataFormat[cellCt]);
226 for (
int cubDeg=0; cubDeg <= maxDeg[cellCt]; cubDeg++) {
229 int oldErrorFlag = errorFlag;
230 for (
int xDeg=0; xDeg <= cubDeg; xDeg++) {
231 for (
int yDeg=0; yDeg <= cubDeg-xDeg; yDeg++) {
232 for (
int zDeg=0; zDeg <= cubDeg-xDeg-yDeg; zDeg++) {
233 double abstol = ( analyticInt[polyCt+offset][0] == 0.0 ? reltol : std::fabs(reltol*analyticInt[polyCt+offset][0]) );
234 double absdiff = std::fabs(analyticInt[polyCt+offset][0] - testInt[cubDeg][polyCt]);
236 *outStream <<
"Cubature order " << std::setw(2) << std::left << cubDeg <<
" integrating "
237 <<
"x^" << std::setw(2) << std::left << xDeg <<
" * y^" << std::setw(2) << yDeg
238 <<
" * z^" << std::setw(2) << zDeg <<
":" <<
" "
239 << std::scientific << std::setprecision(16)
240 << testInt[cubDeg][polyCt] <<
" " << analyticInt[polyCt+offset][0] <<
" "
241 << std::setprecision(4) << absdiff <<
" " <<
"<?" <<
" " << abstol <<
"\n";
242 if (absdiff > abstol) {
244 *outStream << std::right << std::setw(118) <<
"^^^^---FAILURE!\n";
248 offset = offset + maxOffset[cellCt] - cubDeg;
250 offset = offset + (maxOffset[cellCt] - cubDeg)*(maxOffset[cellCt] - cubDeg + 1)/2;
252 *outStream <<
"Cubature order " << std::setw(2) << std::left << cubDeg;
253 if (errorFlag == oldErrorFlag)
254 *outStream <<
" passed.\n";
256 *outStream <<
" failed.\n";
261 catch (
const std::logic_error & err) {
262 *outStream << err.what() <<
"\n";
268 std::cout <<
"End Result: TEST FAILED\n";
270 std::cout <<
"End Result: TEST PASSED\n";
273 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.