53 #include "Teuchos_oblackholestream.hpp" 
   54 #include "Teuchos_RCP.hpp" 
   55 #include "Teuchos_GlobalMPISession.hpp" 
   57 using namespace Intrepid;
 
   68   polydeg[0] = xDeg; polydeg[1] = yDeg; polydeg[2] = zDeg;
 
   70     val *= std::pow(p(i),polydeg[i]);
 
   79 double computeIntegral(shards::CellTopology & cellTopology, 
int cubDegree, 
int xDeg, 
int yDeg, 
int zDeg) {
 
   82   Teuchos::RCP<Cubature<double> > myCub = cubFactory.
create(cellTopology, cubDegree); 
 
   85   int cubDim       = myCub->getDimension();
 
   86   int numCubPoints = myCub->getNumPoints();
 
   92   myCub->getCubature(cubPoints, cubWeights);
 
   94   for (
int i=0; i<numCubPoints; i++) {
 
   95     for (
int j=0; j<cubDim; j++) {
 
   96       point(j) = cubPoints(i,j);
 
   98     val += computeMonomial(point, xDeg, yDeg, zDeg)*cubWeights(i);
 
  105 int main(
int argc, 
char *argv[]) {
 
  107   Teuchos::GlobalMPISession mpiSession(&argc, &argv);
 
  111   int iprint     = argc - 1;
 
  112   Teuchos::RCP<std::ostream> outStream;
 
  113   Teuchos::oblackholestream bhs; 
 
  115     outStream = Teuchos::rcp(&std::cout, 
false);
 
  117     outStream = Teuchos::rcp(&bhs, 
false);
 
  120   Teuchos::oblackholestream oldFormatState;
 
  121   oldFormatState.copyfmt(std::cout);
 
  124   << 
"===============================================================================\n" \
 
  126   << 
"|       Unit Test (CubatureDirect,CubatureTensor,DefaultCubatureFactory)      |\n" \
 
  128   << 
"|     1) Computing integrals of monomials on reference cells in 3D            |\n" \
 
  129   << 
"|                - no BLAS, i.e. standard addition loops -                    |\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 3D (non-BLAS version)                     |\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+04 * INTREPID_TOL;
 
  163   for (
int i=0; i<4; i++) {
 
  164     numPoly[i] = (maxDeg[i]+1)*(maxDeg[i]+2)*(maxDeg[i]+3)/6;
 
  166   for (
int i=0; i<4; i++) {
 
  167     numAnalytic[i] = (maxOffset[i]+1)*(maxOffset[i]+2)*(maxOffset[i]+3)/6;
 
  171   std::string basedir = 
"./data";
 
  172   std::stringstream namestream[4];
 
  173   std::string filename[4];
 
  174   namestream[0] << basedir << 
"/TET_integrals" << 
".dat";
 
  175   namestream[0] >> filename[0];
 
  176   namestream[1] << basedir << 
"/HEX_integrals" << 
".dat";
 
  177   namestream[1] >> filename[1];
 
  178   namestream[2] << basedir << 
"/TRIPRISM_integrals" << 
".dat";
 
  179   namestream[2] >> filename[2];
 
  180   namestream[3] << basedir << 
"/PYR_integrals" << 
".dat";
 
  181   namestream[3] >> filename[3];
 
  184   shards::CellTopology cellType[] = {shards::getCellTopologyData< shards::Tetrahedron<> >(),
 
  185                                      shards::getCellTopologyData< shards::Hexahedron<> >(),
 
  186                                      shards::getCellTopologyData< shards::Wedge<> >(),
 
  187                                      shards::getCellTopologyData< shards::Pyramid<> >() };
 
  189   TypeOfExactData dataFormat[] = {INTREPID_UTILS_SCALAR, INTREPID_UTILS_FRACTION, INTREPID_UTILS_FRACTION, INTREPID_UTILS_FRACTION};
 
  193     for (
int cellCt=0; cellCt < 4; cellCt++) {
 
  194       testInt.assign(numPoly[cellCt], tmparray);
 
  195       analyticInt.assign(numAnalytic[cellCt], tmparray);
 
  196       *outStream << 
"\nIntegrals of monomials on a reference " << cellType[cellCt].getBaseCellTopologyData()->name << 
":\n";
 
  197       std::ifstream filecompare(&filename[cellCt][0]);
 
  199       for (
int cubDeg=0; cubDeg <= maxDeg[cellCt]; cubDeg++) {
 
  201         testInt[cubDeg].resize((cubDeg+1)*(cubDeg+2)*(cubDeg+3)/6);
 
  202         for (
int xDeg=0; xDeg <= cubDeg; xDeg++) {
 
  203           for (
int yDeg=0; yDeg <= cubDeg-xDeg; yDeg++) {
 
  204             for (
int zDeg=0; zDeg <= cubDeg-xDeg-yDeg; zDeg++) {
 
  205               testInt[cubDeg][polyCt] = computeIntegral(cellType[cellCt], cubDeg, xDeg, yDeg, zDeg);
 
  212       if (filecompare.is_open()) {
 
  213         getAnalytic(analyticInt, filecompare, dataFormat[cellCt]);
 
  218       for (
int cubDeg=0; cubDeg <= maxDeg[cellCt]; cubDeg++) {
 
  221         int oldErrorFlag = errorFlag;
 
  222         for (
int xDeg=0; xDeg <= cubDeg; xDeg++) {
 
  223           for (
int yDeg=0; yDeg <= cubDeg-xDeg; yDeg++) {
 
  224             for (
int zDeg=0; zDeg <= cubDeg-xDeg-yDeg; zDeg++) {
 
  225               double abstol = ( analyticInt[polyCt+offset][0] == 0.0 ? reltol : std::fabs(reltol*analyticInt[polyCt+offset][0]) );
 
  226               double absdiff = std::fabs(analyticInt[polyCt+offset][0] - testInt[cubDeg][polyCt]);
 
  227               if (absdiff > abstol) {
 
  228                 *outStream << 
"Cubature order " << std::setw(2) << std::left << cubDeg << 
" integrating " 
  229                            << 
"x^" << std::setw(2) << std::left << xDeg << 
" * y^" << std::setw(2) << yDeg
 
  230                            << 
" * z^" << std::setw(2) << zDeg << 
":" << 
"   " 
  231                            << std::scientific << std::setprecision(16)
 
  232                            << testInt[cubDeg][polyCt] << 
"   " << analyticInt[polyCt+offset][0] << 
"   " 
  233                            << std::setprecision(4) << absdiff << 
"   " << 
"<?" << 
"   " << abstol << 
"\n";
 
  235                 *outStream << std::right << std::setw(118) << 
"^^^^---FAILURE!\n";
 
  239             offset = offset + maxOffset[cellCt] - cubDeg;
 
  241           offset = offset + (maxOffset[cellCt] - cubDeg)*(maxOffset[cellCt] - cubDeg + 1)/2;
 
  243         *outStream << 
"Cubature order " << std::setw(2) << std::left << cubDeg;
 
  244         if (errorFlag == oldErrorFlag)
 
  245          *outStream << 
" passed.\n";
 
  247          *outStream << 
" failed.\n";
 
  252   catch (std::logic_error err) {
 
  253     *outStream << err.what() << 
"\n";
 
  259     std::cout << 
"End Result: TEST FAILED\n";
 
  261     std::cout << 
"End Result: TEST PASSED\n";
 
  264   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.