54 #include "Teuchos_oblackholestream.hpp" 
   55 #include "Teuchos_RCP.hpp" 
   56 #include "Teuchos_GlobalMPISession.hpp" 
   58 using namespace Intrepid;
 
   63 long double evalQuad(std::vector<int> power, 
int dimension, 
int order, 
 
   64                      std::vector<EIntrepidBurkardt> rule, 
 
   65                      std::vector<EIntrepidGrowth> growth) {
 
   68   int size = lineCub.getNumPoints();
 
   71   lineCub.getCubature(cubPoints,cubWeights);
 
   79   int l1      = growthRule1D(order,growth[0],rule[0]);
 
   80   int l2      = growthRule1D(order,growth[1],rule[1]);
 
   84   for (
int i=0; i<l1; i++) {
 
   85     locMid = i*l1+mid2; Qi = 0.0;
 
   87       Qi = cubWeights(locMid)*powl(cubPoints(locMid,1),power[1]); cnt++;
 
   88       for (
int j=1; j<=mid2; j++) {
 
   89         Qi += cubWeights(locMid-j)*powl(cubPoints(locMid-j,1),power[1]) 
 
   90           +cubWeights(locMid+j)*powl(cubPoints(locMid+j,1),power[1]); cnt += 2;
 
   94       for (
int j=0; j<mid2; j++) {
 
   95         Qi += cubWeights(locMid-j)*powl(cubPoints(locMid-j,1),power[1]) 
 
   96           +cubWeights(locMid+j+1)*powl(cubPoints(locMid+j+1,1),power[1]); cnt += 2;
 
   99     Qi *= powl(cubPoints(locMid,0),power[0]);
 
  126 long double evalInt(
int dimension, std::vector<int> power, std::vector<EIntrepidBurkardt> rule) {
 
  129   for (
int i=0; i<dimension; i++) {
 
  130     if (rule[i]==BURK_CLENSHAWCURTIS||rule[i]==BURK_FEJER2||
 
  131         rule[i]==BURK_LEGENDRE||rule[i]==BURK_PATTERSON || 
 
  132         rule[i]==BURK_TRAPEZOIDAL) {
 
  136         I *= 2.0/((
long double)power[i]+1.0);
 
  138     else if (rule[i]==BURK_LAGUERRE) {
 
  139       I *= tgammal((
long double)(power[i]+1));
 
  141     else if (rule[i]==BURK_CHEBYSHEV1) {
 
  142       long double bot, top;
 
  145         for (
int j=2;j<=power[i];j+=2) {
 
  146           top *= (
long double)(j-1);
 
  147           bot *= (
long double)j;
 
  155     else if (rule[i]==BURK_CHEBYSHEV2) {
 
  156       long double bot, top;
 
  159       for (
int j=2;j<=power[i];j+=2) {
 
  160         top *= (
long double)(j-1);
 
  161         bot *= (
long double)j;
 
  163       bot *= (
long double)(power[i]+2);
 
  170     else if (rule[i]==BURK_HERMITE||rule[i]==BURK_GENZKEISTER) {
 
  175         long double value = 1.0;
 
  176         if ((power[i]-1)>=1) {
 
  177           int n_copy = power[i]-1;
 
  179             value  *= (
long double)n_copy;
 
  183         I *= value*sqrt(M_PI)/powl(2.0,(
long double)power[i]/2.0);
 
  190 int main(
int argc, 
char *argv[]) {
 
  192   Teuchos::GlobalMPISession mpiSession(&argc, &argv);
 
  196   int iprint     = argc - 1;
 
  197   Teuchos::RCP<std::ostream> outStream;
 
  198   Teuchos::oblackholestream bhs; 
 
  200     outStream = Teuchos::rcp(&std::cout, 
false);
 
  202     outStream = Teuchos::rcp(&bhs, 
false);
 
  205   Teuchos::oblackholestream oldFormatState;
 
  206   oldFormatState.copyfmt(std::cout);
 
  209   << 
"===============================================================================\n" \
 
  211   << 
"|                         Unit Test (CubatureTensorSorted)                    |\n" \
 
  213   << 
"|     1) Computing integrals of monomials in 2D                               |\n" \
 
  215   << 
"|  Questions? Contact  Drew Kouri (dpkouri@sandia.gov) or                     |\n" \
 
  216   << 
"|                      Denis Ridzal (dridzal@sandia.gov).                     |\n" \
 
  218   << 
"|  Intrepid's website: http://trilinos.sandia.gov/packages/intrepid           |\n" \
 
  219   << 
"|  Trilinos website:   http://trilinos.sandia.gov                             |\n" \
 
  221   << 
"===============================================================================\n"\
 
  222   << 
"| TEST 14: integrals of monomials in 2D - Isotropic with growth rules         |\n"\
 
  223   << 
"===============================================================================\n";
 
  229   long double reltol      = 1.0e+06*INTREPID_TOL;
 
  232   long double analyticInt = 0;
 
  233   long double testInt     = 0;
 
  236   std::vector<int> power(2,0);
 
  237   std::vector<EIntrepidBurkardt> rule1(2,BURK_CLENSHAWCURTIS);
 
  239   std::vector<EIntrepidGrowth> growth(2,GROWTH_DEFAULT);
 
  241   *outStream << 
"\nIntegrals of monomials on a reference line (edge):\n";
 
  244     for (EIntrepidBurkardt rule=BURK_CHEBYSHEV1; rule <= BURK_LAGUERRE; rule++) {   
 
  246       rule1[0] = rule; rule1[1] = rule;
 
  247       if (rule!=BURK_PATTERSON&&rule!=BURK_GENZKEISTER&&rule!=BURK_TRAPEZOIDAL) {
 
  248         *outStream << 
"Testing " << EIntrepidBurkardtToString(rule) << 
"\n";
 
  249         for (
int i=1; i <= maxOrder; i++) {
 
  250           l1 = growthRule1D(i,growth[0],rule);
 
  251           l2 = growthRule1D(i,growth[1],rule);
 
  252           if (rule==BURK_CHEBYSHEV1||rule==BURK_CHEBYSHEV2||rule==BURK_LEGENDRE||
 
  253               rule==BURK_LAGUERRE||rule==BURK_HERMITE) {
 
  257           else if (rule==BURK_CLENSHAWCURTIS||rule==BURK_FEJER2) {
 
  263           for (
int j=0; j <= maxDegx; j++) {
 
  265             for (
int k=0; k <= maxDegy; k++) {
 
  267               analyticInt = evalInt(dimension, power, rule1);
 
  268               testInt     = evalQuad(power,dimension,order,rule1,growth);
 
  270               long double abstol  = (analyticInt == 0.0 ? reltol : std::fabs(reltol*analyticInt) );
 
  271               long double absdiff = std::fabs(analyticInt - testInt);
 
  272               *outStream << 
"Cubature order (" << std::setw(2) << std::left 
 
  273                          << l1 << 
", " << std::setw(2) << std::left << l2
 
  275                          << 
"x^" << std::setw(2) << std::left << j << 
"y^"  
  276                          << std::setw(2) << std::left 
 
  277                          << k <<  
":" << 
"   " << std::scientific 
 
  278                          << std::setprecision(16) << testInt 
 
  279                          << 
"   " << analyticInt << 
"   "  
  280                          << std::setprecision(4) << absdiff << 
"   "  
  281                          << 
"<?" << 
"   " << abstol << 
"\n";
 
  282               if (absdiff > abstol) {
 
  284                 *outStream << std::right << std::setw(104) 
 
  285                            << 
"^^^^---FAILURE!\n";
 
  293       else if (rule==BURK_PATTERSON) {
 
  294         *outStream << 
"Testing " << EIntrepidBurkardtToString(rule) << 
"\n";
 
  295         for (
int i=0; i < 3; i++) {       
 
  296           l1 = growthRule1D(i,growth[0],rule);
 
  297           l2 = growthRule1D(i,growth[1],rule);
 
  303             maxDegx = (int)(1.5*(
double)l1-0.5);
 
  304             maxDegy = (int)(1.5*(
double)l2-0.5);
 
  308           for (
int j=0; j <= maxDegx; j++) {        
 
  310             for (
int k=0; k <= maxDegy; k++) {        
 
  312               analyticInt = evalInt(dimension, power, rule1);
 
  313               testInt     = evalQuad(power,dimension,order,rule1,growth);
 
  315               long double abstol  = (analyticInt == 0.0 ? reltol : std::fabs(reltol*analyticInt) );
 
  316               long double absdiff = std::fabs(analyticInt - testInt);
 
  317               *outStream << 
"Cubature order (" << std::setw(2) << std::left 
 
  318                          << l1 << 
", " << std::setw(2) << std::left << l2 
 
  320                          << 
"x^" << std::setw(2) << std::left << j << 
"y^"  
  321                          << std::setw(2) << std::left 
 
  322                          << k <<  
":" << 
"   " << std::scientific 
 
  323                          << std::setprecision(16) << testInt 
 
  324                          << 
"   " << analyticInt << 
"   "  
  325                          << std::setprecision(4) << absdiff << 
"   "  
  326                          << 
"<?" << 
"   " << abstol << 
"\n";
 
  327               if (absdiff > abstol) {
 
  329                 *outStream << std::right << std::setw(104) 
 
  330                            << 
"^^^^---FAILURE!\n";
 
  338       else if (rule==BURK_GENZKEISTER) {
 
  339         *outStream << 
"Testing " << EIntrepidBurkardtToString(rule) << 
"\n";
 
  340         for (
int i=0; i < 3; i++) {
 
  341           l1 = growthRule1D(i,growth[0],rule);
 
  342           l2 = growthRule1D(i,growth[1],rule);
 
  348             maxDegx = (int)(1.5*(
double)l1-0.5);
 
  349             maxDegy = (int)(1.5*(
double)l2-0.5);
 
  353           for (
int j=0; j <= maxDegx; j++) {        
 
  355             for (
int k=0; k <= maxDegy; k++) {        
 
  357               analyticInt = evalInt(dimension, power, rule1);
 
  358               testInt     = evalQuad(power,dimension,order,rule1,growth);
 
  360               long double abstol  = (analyticInt == 0.0 ? reltol : std::fabs(reltol*analyticInt) );
 
  361               long double absdiff = std::fabs(analyticInt - testInt);
 
  362               *outStream << 
"Cubature order (" << std::setw(2) << std::left 
 
  363                          << l1 << 
", " << std::setw(2) << std::left << l2
 
  365                          << 
"x^" << std::setw(2) << std::left << j << 
"y^"  
  366                          << std::setw(2) << std::left 
 
  367                          << k <<  
":" << 
"   " << std::scientific 
 
  368                          << std::setprecision(16) << testInt 
 
  369                          << 
"   " << analyticInt << 
"   "  
  370                          << std::setprecision(4) << absdiff << 
"   "  
  371                          << 
"<?" << 
"   " << abstol << 
"\n";
 
  372               if (absdiff > abstol) {
 
  374                 *outStream << std::right << std::setw(104) 
 
  375                            << 
"^^^^---FAILURE!\n";
 
  385   catch (std::logic_error err) {
 
  386     *outStream << err.what() << 
"\n";
 
  392     std::cout << 
"End Result: TEST FAILED\n";
 
  394     std::cout << 
"End Result: TEST PASSED\n";
 
  397   std::cout.copyfmt(oldFormatState);
 
Implementation of a templated lexicographical container for a multi-indexed scalar quantity...
Header file for the Intrepid::CubatureTensorSorted class. 
Utilizes 1D cubature (integration) rules contained in the library sandia_rules (John Burkardt...