ROOT logo
ROOT » MATH » MATHCORE » ROOT::Math::GaussIntegrator

class ROOT::Math::GaussIntegrator: public ROOT::Math::VirtualIntegratorOneDim


   User class for performing function integration.

   It will use the Gauss Method for function integration in a given interval.
   This class is implemented from TF1::Integral().

   @ingroup Integration


Function Members (Methods)

public:
virtual~GaussIntegrator()
voidAbsValue(bool flag)
virtual doubleError() const
ROOT::Math::GaussIntegratorGaussIntegrator(double relTol = 1.E-12)
ROOT::Math::GaussIntegratorGaussIntegrator(const ROOT::Math::GaussIntegrator&)
virtual doubleIntegral()
virtual doubleIntegral(const vector<double>& pts)
virtual doubleIntegral(double a, double b)
virtual doubleIntegralCauchy(double a, double b, double c)
virtual doubleIntegralLow(double b)
virtual doubleIntegralUp(double a)
virtual intROOT::Math::VirtualIntegrator::NEval() const
ROOT::Math::VirtualIntegratorOneDim&ROOT::Math::VirtualIntegratorOneDim::operator=(const ROOT::Math::VirtualIntegratorOneDim&)
virtual ROOT::Math::IntegratorOneDimOptionsOptions() const
virtual doubleResult() const
virtual voidSetAbsTolerance(double)
virtual voidSetFunction(const ROOT::Math::IGenFunction&)
virtual voidSetOptions(const ROOT::Math::IntegratorOneDimOptions& opt)
virtual voidSetRelTolerance(double)
virtual intStatus() const
virtual ROOT::Math::IntegrationOneDim::TypeROOT::Math::VirtualIntegratorOneDim::Type() const
private:
virtual doubleDoIntegral(double a, double b, const ROOT::Math::IGenFunction* func)

Data Members

protected:
doublefEpsilonRelative error.
const ROOT::Math::IGenFunction*fFunctionPointer to function used.
doublefLastErrorError from the last stimation.
doublefLastResultResult from the last stimation.
boolfUsedOnceBool value to check if the function was at least called once.
static boolfgAbsValueAbsValue used for the calculation of the integral

Class Charts

Inheritance Inherited Members Includes Libraries
Class Charts

Function documentation

virtual ~GaussIntegrator()
 Destructor 
GaussIntegrator(double relTol = 1.E-12)
 Default Constructor. 
void AbsValue(bool flag)
 Static function: set the fgAbsValue flag.
       By default TF1::Integral uses the original function value to compute the integral
       However, TF1::Moment, CentralMoment require to compute the integral
       using the absolute value of the function.

void SetRelTolerance(double )
 Implementing VirtualIntegrator Interface
 Set the desired relative Error. 
void SetAbsTolerance(double )
 This method is not implemented. 
double Result() const
 Returns the result of the last Integral calculation. 
double Error() const
 Return the estimate of the absolute Error of the last Integral calculation. 
int Status() const
 return the status of the last integration - 0 in case of success 
double Integral(double a, double b)
 Implementing VirtualIntegratorOneDim Interface

     Returns Integral of function between a and b.
     Based on original CERNLIB routine DGAUSS by Sigfried Kolbig
     converted to C++ by Rene Brun

     This function computes, to an attempted specified accuracy, the value
     of the integral.

    Method:
       For any interval [a,b] we define g8(a,b) and g16(a,b) to be the 8-point
       and 16-point Gaussian quadrature approximations to

I = #int^{b}_{a} f(x)dx
      and define

r(a,b) = #frac{#||{g_{16}(a,b)-g_{8}(a,b)}}{1+#||{g_{16}(a,b)}}
      Then,

G = #sum_{i=1}^{k}g_{16}(x_{i-1},x_{i})
      where, starting with x0 = A and finishing with xk = B,
      the subdivision points xi(i=1,2,...) are given by

x_{i} = x_{i-1} + #lambda(B-x_{i-1})

#lambda
      is equal to the first member of the
      sequence 1,1/2,1/4,... for which r(xi-1, xi) < EPS.
      If, at any stage in the process of subdivision, the ratio

q = #||{#frac{x_{i}-x_{i-1}}{B-A}}
      is so small that 1+0.005q is indistinguishable from 1 to
      machine accuracy, an error exit occurs with the function value
      set equal to zero.

   Accuracy:
      Unless there is severe cancellation of positive and negative values of
      f(x) over the interval [A,B], the relative error may be considered as
      specifying a bound on the <I>relative</I> error of I in the case
      |I|&gt;1, and a bound on the absolute error in the case |I|&lt;1. More
      precisely, if k is the number of sub-intervals contributing to the
      approximation (see Method), and if

I_{abs} = #int^{B}_{A} #||{f(x)}dx
      then the relation

#frac{#||{G-I}}{I_{abs}+k} < EPS
      will nearly always be true, provided the routine terminates without
      printing an error message. For functions f having no singularities in
      the closed interval [A,B] the accuracy will usually be much higher than
      this.

    Error handling:
      The requested accuracy cannot be obtained (see Method).
      The function value is set equal to zero.

    Note 1:
      Values of the function f(x) at the interval end-points A and B are not
      required. The subprogram may therefore be used when these values are
      undefined

double Integral()
 Returns Integral of function on an infinite interval.
      This function computes, to an attempted specified accuracy, the value of the integral:

I = #int^{#infinity}_{-#infinity} f(x)dx
      Usage:
        In any arithmetic expression, this function has the approximate value
        of the integral I.

      The integral is mapped onto [0,1] using a transformation then integral computation is surrogated to DoIntegral.

double IntegralUp(double a)
 Returns Integral of function on an upper semi-infinite interval.
      This function computes, to an attempted specified accuracy, the value of the integral:

I = #int^{#infinity}_{A} f(x)dx
      Usage:
        In any arithmetic expression, this function has the approximate value
        of the integral I.
        - A: lower end-point of integration interval.

      The integral is mapped onto [0,1] using a transformation then integral computation is surrogated to DoIntegral.

double IntegralLow(double b)
 Returns Integral of function on a lower semi-infinite interval.
       This function computes, to an attempted specified accuracy, the value of the integral:

I = #int^{B}_{#infinity} f(x)dx
      Usage:
         In any arithmetic expression, this function has the approximate value
         of the integral I.
         - B: upper end-point of integration interval.

      The integral is mapped onto [0,1] using a transformation then integral computation is surrogated to DoIntegral.

void SetFunction(const ROOT::Math::IGenFunction& )
 Set integration function (flag control if function must be copied inside).
       \@param f Function to be used in the calculations.

double Integral(const vector<double>& pts)
 This method is not implemented. 
double IntegralCauchy(double a, double b, double c)
 This method is not implemented. 
void SetOptions(const ROOT::Math::IntegratorOneDimOptions& opt)
 set the options
double DoIntegral(double a, double b, const ROOT::Math::IGenFunction* func)
      Integration surrugate method. Return integral of passed function in  interval [a,b]
      Derived class (like GaussLegendreIntegrator)  can re-implement this method to modify to use
      an improved algorithm