ATLAS Offline Software
Loading...
Searching...
No Matches
eflowRecursiveGaussLegendreIntegrator< IntegrandType > Class Template Reference

Class to perform a generic recursive Gauss-Legendre Integration, see http://en.wikipedia.org/wiki/Gaussian_quadrature#Gauss.E2.80.93Legendre_quadrature Templated in the type of the integrand. More...

#include <eflowCellIntegrator.h>

Collaboration diagram for eflowRecursiveGaussLegendreIntegrator< IntegrandType >:

Public Member Functions

 eflowRecursiveGaussLegendreIntegrator (IntegrandType *integrand, double error)
virtual ~eflowRecursiveGaussLegendreIntegrator ()
virtual double integrate (const eflowRange &range)
double getError () const

Private Member Functions

double DoGaussLegendreIntegration (const eflowRange &range, int nOrder)
double RecurseIntegration (const eflowRange &range, int nSubRanges)

Private Attributes

IntegrandType * m_integrand
double m_error
int m_depth

Detailed Description

template<class IntegrandType>
class eflowRecursiveGaussLegendreIntegrator< IntegrandType >

Class to perform a generic recursive Gauss-Legendre Integration, see http://en.wikipedia.org/wiki/Gaussian_quadrature#Gauss.E2.80.93Legendre_quadrature Templated in the type of the integrand.

All the integrand class needs to have, is an evaluate() method, that takes a double as a parameter and returns a double. Templates are used instead of polymorphism here due to the better CPU performance.

Definition at line 34 of file eflowCellIntegrator.h.

Constructor & Destructor Documentation

◆ eflowRecursiveGaussLegendreIntegrator()

template<class IntegrandType>
eflowRecursiveGaussLegendreIntegrator< IntegrandType >::eflowRecursiveGaussLegendreIntegrator ( IntegrandType * integrand,
double error )
inline

Definition at line 36 of file eflowCellIntegrator.h.

36 :
39 }
Class to perform a generic recursive Gauss-Legendre Integration, see http://en.wikipedia....

◆ ~eflowRecursiveGaussLegendreIntegrator()

template<class IntegrandType>
virtual eflowRecursiveGaussLegendreIntegrator< IntegrandType >::~eflowRecursiveGaussLegendreIntegrator ( )
inlinevirtual

Definition at line 40 of file eflowCellIntegrator.h.

40{ }

Member Function Documentation

◆ DoGaussLegendreIntegration()

template<class IntegrandType>
double eflowRecursiveGaussLegendreIntegrator< IntegrandType >::DoGaussLegendreIntegration ( const eflowRange & range,
int nOrder )
inlineprivate

Definition at line 61 of file eflowCellIntegrator.h.

61 {
62 /* Perform nth order Gauss-Legendre quadrature, see
63 * http://en.wikipedia.org/wiki/Gaussian_quadrature#Gauss.E2.80.93Legendre_quadrature */
64
65 /* Array offset for legendre weights/roots */
66 int j = nOrder * (nOrder - 1) / 2;
67 const double* roots = legendreRoots + j;
68 const double* weights = legendreWeights + j;
69
70 double rangeCenter = range.getCenter();
71 double rangeHalfWidth = range.getWidth()/2;
72
73 double I = 0.0;
74 double x;
75 if(m_integrand){
76 for (int i = 0; i < nOrder; i++) {
78 I += weights[i] * m_integrand->evaluate(x);
79 }
80 }
81 else{
82 std::cerr << "eflowCellIntergrator::DoGaussLegendreIntegration ERROR : Invalid pointer to IntegrandType and cannot perform integration" << std::endl;
83 }
84 if (I < 0. || rangeHalfWidth < 0.) std::cerr << "eflowCellIntergrator::DoGaussLegendreIntegration WARNING: I = " << I << "\trange = " << range.print() << std::endl;
85
86 return I * rangeHalfWidth;
87 }

◆ getError()

template<class IntegrandType>
double eflowRecursiveGaussLegendreIntegrator< IntegrandType >::getError ( ) const
inline

Definition at line 57 of file eflowCellIntegrator.h.

57{return m_error;}

◆ integrate()

template<class IntegrandType>
virtual double eflowRecursiveGaussLegendreIntegrator< IntegrandType >::integrate ( const eflowRange & range)
inlinevirtual

Definition at line 42 of file eflowCellIntegrator.h.

42 {
43 /* Try a 5th and 6th order Gauss-Legendre quadrature */
46
47 /* If results agree within m_error threshold, return the mean... */
48 double Imean = (I5 + I6) / 2.0;
49 if (fabs(I5 - I6) / Imean <= m_error) {
50 return Imean;
51 } else {
52 /* ...else recursively part up the integration into n subRanges */
53 return RecurseIntegration(range, 2);
54 }
55 }
double DoGaussLegendreIntegration(const eflowRange &range, int nOrder)
double RecurseIntegration(const eflowRange &range, int nSubRanges)

◆ RecurseIntegration()

template<class IntegrandType>
double eflowRecursiveGaussLegendreIntegrator< IntegrandType >::RecurseIntegration ( const eflowRange & range,
int nSubRanges )
inlineprivate

Definition at line 89 of file eflowCellIntegrator.h.

89 {
90 m_depth++;
91 /* Part up integration range in n subRanges */
92 double subRangeWidth = range.getWidth() / nSubRanges;
93 eflowRange subRange(range.getMin(), range.getMin() + subRangeWidth);
94
95 /* Loop over subranges */
96 double Integral(0.0);
97 for (int i = 0; i < nSubRanges; i++) {
100 }
101 m_depth--;
102
103 return Integral;
104 }
virtual double integrate(const eflowRange &range)

Member Data Documentation

◆ m_depth

template<class IntegrandType>
int eflowRecursiveGaussLegendreIntegrator< IntegrandType >::m_depth
private

Definition at line 108 of file eflowCellIntegrator.h.

◆ m_error

template<class IntegrandType>
double eflowRecursiveGaussLegendreIntegrator< IntegrandType >::m_error
private

Definition at line 107 of file eflowCellIntegrator.h.

◆ m_integrand

template<class IntegrandType>
IntegrandType* eflowRecursiveGaussLegendreIntegrator< IntegrandType >::m_integrand
private

Definition at line 106 of file eflowCellIntegrator.h.


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