queso-0.51.1
Protected Attributes | List of all members
QUESO::GaussianHermite1DQuadrature Class Reference

Class for Hermite-Gauss quadrature rule for one-dimensional functions. More...

#include <1DQuadrature.h>

Inheritance diagram for QUESO::GaussianHermite1DQuadrature:
Inheritance graph
[legend]
Collaboration diagram for QUESO::GaussianHermite1DQuadrature:
Collaboration graph
[legend]

Public Member Functions

Constructor/Destructor methods
 GaussianHermite1DQuadrature (double mean, double stddev, unsigned int order)
 Default constructor. More...
 
 ~GaussianHermite1DQuadrature ()
 Destructor. More...
 
Mathematical methods
void dumbRoutine () const
 A bogus method. More...
 
- Public Member Functions inherited from QUESO::Base1DQuadrature
 Base1DQuadrature (double minDomainValue, double maxDomainValue, unsigned int order)
 Default constructor. More...
 
virtual ~Base1DQuadrature ()
 Virtual destructor. More...
 
double minDomainValue () const
 Returns the minimum value of the domain of the (one-dimensional) function. More...
 
double maxDomainValue () const
 Returns the maximum value of the domain of the (one-dimensional) function. More...
 
unsigned int order () const
 Returns the order of the quadrature rule. More...
 
const std::vector< double > & positions () const
 Array of the positions for the numerical integration. More...
 
const std::vector< double > & weights () const
 Array of the weights used in the numerical integration. More...
 

Protected Attributes

double m_mean
 
double m_stddev
 
- Protected Attributes inherited from QUESO::Base1DQuadrature
double m_minDomainValue
 
double m_maxDomainValue
 
unsigned int m_order
 
std::vector< double > m_positions
 
std::vector< double > m_weights
 

Detailed Description

Class for Hermite-Gauss quadrature rule for one-dimensional functions.

Hermite-Gauss quadrature, also called Hermite quadrature, is a Gaussian quadrature over the interval $(-\infty,\infty)$ with weighting function $ W(x)=e^{-x^2}$.
The abscissas for quadrature order $ n $ are given by the roots $ x_i $ of the Hermite polynomials $ H_n(x)$, which occur symmetrically about 0.
The abscissas and weights can be computed analytically for small $ n $:

$ n $$ x_i $$ w_i $
2 $\pm \frac{1}{2}\sqrt{2} $ $ \frac{1}{2}\sqrt{\pi} $
3 $ 0 $ $ \frac{2}{3}\sqrt{\pi} $
$\pm \frac{1}{2}\sqrt{6} $ $ \frac{1}{6}\sqrt{\pi} $
4 $\pm \sqrt{\frac{3-\sqrt{6}}{2}} $ $ \frac{\sqrt{\pi}}{4(3-\sqrt{6})} $
$\pm \sqrt{\frac{3-\sqrt{6}}{2}} $ $ \frac{\sqrt{\pi}}{4(3+\sqrt{6})} $
See Also
Weisstein, Eric W. "Hermite-Gauss Quadrature." From MathWorld–A Wolfram Web Resource. http://mathworld.wolfram.com/Hermite-GaussQuadrature.html.

Definition at line 219 of file 1DQuadrature.h.

Constructor & Destructor Documentation

QUESO::GaussianHermite1DQuadrature::GaussianHermite1DQuadrature ( double  mean,
double  stddev,
unsigned int  order 
)

Default constructor.

Constructs a Gaussian-Hermite quadrature of order order. Valid values for the order of the quadrature rule are: 1-9, 19.

Todo:
: Prepare the code to include both parameters mean and stddev.

Definition at line 446 of file 1DQuadrature.C.

References QUESO::Base1DQuadrature::m_order, QUESO::Base1DQuadrature::m_positions, QUESO::Base1DQuadrature::m_weights, UQ_FATAL_TEST_MACRO, and QUESO::UQ_UNAVAILABLE_RANK.

450  :
451  Base1DQuadrature(-INFINITY,INFINITY,order),
452  m_mean (mean),
453  m_stddev(stddev)
454 {
455  // FIX ME: prepare code for mean != 0 and stddev != 1
456  m_positions.resize(m_order+1,0.); // Yes, '+1'
457  m_weights.resize (m_order+1,0.); // Yes, '+1'
458 
459  // http://www.efunda.com/math/num_integration/findgausshermite.cfm
460  switch (m_order) { // eep2011
461  case 1:
462  m_weights [0] = sqrt(M_PI)/2.;
463  m_weights [1] = sqrt(M_PI)/2.;
464 
465  m_positions[0] = -1./sqrt(2.);
466  m_positions[1] = 1./sqrt(2.);
467  break;
468 
469  case 2:
470  m_weights [0] = sqrt(M_PI)/6.;
471  m_weights [1] = 2.*sqrt(M_PI)/3.;
472  m_weights [2] = sqrt(M_PI)/6.;
473 
474  m_positions[0] = -sqrt(1.5);
475  m_positions[1] = 0.;
476  m_positions[2] = sqrt(1.5);
477  break;
478 
479  case 3:
480  m_weights [0] = sqrt(M_PI)/4./(3.-sqrt(6.));
481  m_weights [1] = sqrt(M_PI)/4./(3.+sqrt(6.));
482  m_weights [2] = sqrt(M_PI)/4./(3.+sqrt(6.));
483  m_weights [3] = sqrt(M_PI)/4./(3.-sqrt(6.));
484 
485  m_positions[0] = -sqrt(1.5+sqrt(1.5));
486  m_positions[1] = -sqrt(1.5-sqrt(1.5));
487  m_positions[2] = sqrt(1.5-sqrt(1.5));
488  m_positions[3] = sqrt(1.5+sqrt(1.5));
489  break;
490 
491  case 4:
492  m_weights [0] = 0.019953242049;
493  m_weights [1] = 0.393619323152;
494  m_weights [2] = 0.945308720483;
495  m_weights [3] = 0.393619323152;
496  m_weights [4] = 0.019953242059;
497 
498  m_positions[0] = -sqrt(2.5+sqrt(2.5));
499  m_positions[1] = -sqrt(2.5-sqrt(2.5));
500  m_positions[2] = 0.;
501  m_positions[3] = sqrt(2.5-sqrt(2.5));
502  m_positions[4] = sqrt(2.5+sqrt(2.5));
503  break;
504 
505  case 5:
506  m_weights [0] = 0.00453000990551;
507  m_weights [1] = 0.157067320323;
508  m_weights [2] = 0.724629595224;
509  m_weights [3] = 0.724629595224;
510  m_weights [4] = 0.157067320323;
511  m_weights [5] = 0.00453000990551;
512 
513  m_positions [0] = -2.35060497367;
514  m_positions [1] = -1.33584907401;
515  m_positions [2] = -0.436077411928;
516  m_positions [3] = 0.436077411928;
517  m_positions [4] = 1.33584907401;
518  m_positions [5] = 2.35060497367;
519  break;
520 
521  case 6:
522  m_weights [0] = 0.0009717812451;
523  m_weights [1] = 0.0545155828191;
524  m_weights [2] = 0.42560725261;
525  m_weights [3] = 0.810264617557;
526  m_weights [4] = 0.42560725261;
527  m_weights [5] = 0.0545155828191;
528  m_weights [6] = 0.0009717812451;
529 
530  m_positions [0] = -2.65196135684;
531  m_positions [1] = -1.67355162877;
532  m_positions [2] = -0.816287882859;
533  m_positions [3] = 0.;
534  m_positions [4] = 0.816287882859;
535  m_positions [5] = 1.67355162877;
536  m_positions [6] = 2.65196135684;
537  break;
538 
539  case 7:
540  m_weights [0] = 0.000199604072211;
541  m_weights [1] = 0.0170779830074;
542  m_weights [2] = 0.207802325815;
543  m_weights [3] = 0.661147012558;
544  m_weights [4] = 0.661147012558;
545  m_weights [5] = 0.207802325815;
546  m_weights [6] = 0.0170779830074;
547  m_weights [7] = 0.000199604072211;
548 
549  m_positions [0] = -2.93063742026;
550  m_positions [1] = -1.9816567567;
551  m_positions [2] = -1.15719371245;
552  m_positions [3] = -0.381186990207;
553  m_positions [4] = 0.381186990207;
554  m_positions [5] = 1.15719371245;
555  m_positions [6] = 1.9816567567;
556  m_positions [7] = 2.93063742026;
557  break;
558 
559  case 8:
560  m_weights [0] = 3.96069772633e-5;
561  m_weights [1] = 0.00494362427554;
562  m_weights [2] = 0.0884745273944;
563  m_weights [3] = 0.432651559003;
564  m_weights [4] = 0.720235215606;
565  m_weights [5] = 0.432651559003;
566  m_weights [6] = 0.0884745273944;
567  m_weights [7] = 0.00494362427554;
568  m_weights [8] = 3.96069772633e-5;
569 
570  m_positions [0] = -3.19099320178;
571  m_positions [1] = -2.26658058453;
572  m_positions [2] = -1.46855328922;
573  m_positions [3] = -0.723551018753;
574  m_positions [4] = 0.;
575  m_positions [5] = 0.723551018753;
576  m_positions [6] = 1.46855328922;
577  m_positions [7] = 2.26658058453;
578  m_positions [8] = 3.19099320178;
579  break;
580 
581  case 9:
582  m_weights [0] = 7.64043285523e-6;
583  m_weights [1] = 0.00134364574678;
584  m_weights [2] = 0.0338743944555;
585  m_weights [3] = 0.240138611082;
586  m_weights [4] = 0.610862633735;
587  m_weights [5] = 0.610862633735;
588  m_weights [6] = 0.240138611082;
589  m_weights [7] = 0.0338743944555;
590  m_weights [8] = 0.00134364574678;
591  m_weights [9] = 7.64043285523e-6;
592 
593  m_positions [0] = -3.43615911884;
594  m_positions [1] = -2.53273167423;
595  m_positions [2] = -1.7566836493;
596  m_positions [3] = -1.03661082979;
597  m_positions [4] = -0.342901327224;
598  m_positions [5] = 0.342901327224;
599  m_positions [6] = 1.03661082979;
600  m_positions [7] = 1.7566836493;
601  m_positions [8] = 2.53273167423;
602  m_positions [9] = 3.43615911884;
603  break;
604 
605  case 19:
606  m_weights [0] = 2.22939364554e-13;
607  m_weights [1] = 4.39934099226e-10;
608  m_weights [2] = 1.08606937077e-7;
609  m_weights [3] = 7.8025564785e-6;
610  m_weights [4] = 0.000228338636017;
611  m_weights [5] = 0.00324377334224;
612  m_weights [6] = 0.0248105208875;
613  m_weights [7] = 0.10901720602;
614  m_weights [8] = 0.286675505363;
615  m_weights [9] = 0.462243669601;
616  m_weights [10] = 0.462243669601;
617  m_weights [11] = 0.286675505363;
618  m_weights [12] = 0.10901720602;
619  m_weights [13] = 0.0248105208875;
620  m_weights [14] = 0.00324377334224;
621  m_weights [15] = 0.000228338636017;
622  m_weights [16] = 7.8025564785e-6;
623  m_weights [17] = 1.08606937077e-7;
624  m_weights [18] = 4.39934099226e-10;
625  m_weights [19] = 2.22939364554e-13;
626 
627  m_positions [0] = -5.38748089001;
628  m_positions [1] = -4.60368244955;
629  m_positions [2] = -3.94476404012;
630  m_positions [3] = -3.34785456738;
631  m_positions [4] = -2.78880605843;
632  m_positions [5] = -2.25497400209;
633  m_positions [6] = -1.73853771212;
634  m_positions [7] = -1.2340762154;
635  m_positions [8] = -0.737473728545;
636  m_positions [9] = -0.245340708301;
637  m_positions[10] = 0.245340708301;
638  m_positions[11] = 0.737473728545;
639  m_positions[12] = 1.2340762154;
640  m_positions[13] = 1.73853771212;
641  m_positions[14] = 2.25497400209;
642  m_positions[15] = 2.78880605843;
643  m_positions[16] = 3.34785456738;
644  m_positions[17] = 3.94476404012;
645  m_positions[18] = 4.60368244955;
646  m_positions[19] = 5.38748089001;
647  break;
648 
649  default:
650  std::cerr << "In GaussianHermite1DQuadrature::constructor()"
651  << ": m_order = " << m_order
652  << std::endl;
653  UQ_FATAL_TEST_MACRO(true,
655  "GaussianHermite1DQuadrature::constructor()",
656  "order not supported");
657  break;
658  }
659  for (unsigned int j = 0; j < (m_order+1); ++j) {
660  m_weights [j] *= sqrt(2.);
661  m_positions[j] *= sqrt(2.);
662  }
663 }
Base1DQuadrature(double minDomainValue, double maxDomainValue, unsigned int order)
Default constructor.
Definition: 1DQuadrature.C:32
std::vector< double > m_weights
Definition: 1DQuadrature.h:86
const int UQ_UNAVAILABLE_RANK
Definition: Defines.h:74
std::vector< double > m_positions
Definition: 1DQuadrature.h:85
unsigned int order() const
Returns the order of the quadrature rule.
Definition: 1DQuadrature.C:70
#define UQ_FATAL_TEST_MACRO(test, givenRank, where, what)
Definition: Defines.h:223
QUESO::GaussianHermite1DQuadrature::~GaussianHermite1DQuadrature ( )

Destructor.

Definition at line 665 of file 1DQuadrature.C.

666 {
667 }

Member Function Documentation

void QUESO::GaussianHermite1DQuadrature::dumbRoutine ( ) const
virtual

A bogus method.

Implements QUESO::Base1DQuadrature.

Definition at line 670 of file 1DQuadrature.C.

671 {
672  return;
673 }

Member Data Documentation

double QUESO::GaussianHermite1DQuadrature::m_mean
protected

Definition at line 246 of file 1DQuadrature.h.

double QUESO::GaussianHermite1DQuadrature::m_stddev
protected

Definition at line 247 of file 1DQuadrature.h.


The documentation for this class was generated from the following files:

Generated on Thu Apr 23 2015 19:26:17 for queso-0.51.1 by  doxygen 1.8.5