Class for performing numerical integration of a function in one dimension.
It uses the numerical integration algorithms of GSL, which reimplements the
algorithms used in the QUADPACK, a numerical integration package written in Fortran.
Various types of adaptive and non-adaptive integration are supported. These include
integration over infinite and semi-infinite ranges and singular integrals.
The integration type is selected using the Integration::type enumeration
in the class constructor.
The default type is adaptive integration with singularity
(ADAPTIVESINGULAR or QAGS in the QUADPACK convention) applying a Gauss-Kronrod 21-point integration rule.
In the case of ADAPTIVE type, the integration rule can also be specified via the
Integration::GKRule. The default rule is 31 points.
In the case of integration over infinite and semi-infinite ranges, the type used is always
ADAPTIVESINGULAR applying a transformation from the original interval into (0,1).
The ADAPTIVESINGULAR type is the most sophicticated type. When performances are
important, it is then recommened to use the NONADAPTIVE type in case of smooth functions or
ADAPTIVE with a lower Gauss-Kronrod rule.
For detailed description on GSL integration algorithms see the
<A HREF="http://www.gnu.org/software/gsl/manual/html_node/Numerical-Integration.html">GSL Manual</A>.
@ingroup Integration
virtual | ~GSLIntegrator() |
virtual double | Error() const |
ROOT::Math::IntegrationOneDim::Type | GetType() const |
const char* | GetTypeName() const |
ROOT::Math::GSLIntegrator | GSLIntegrator(double absTol = 1.0000000000000001E-9, double relTol = 9.9999999999999995E-7, size_t size = 1000) |
ROOT::Math::GSLIntegrator | GSLIntegrator(const ROOT::Math::Integration::Type type, double absTol = 1.0000000000000001E-9, double relTol = 9.9999999999999995E-7, size_t size = 1000) |
ROOT::Math::GSLIntegrator | GSLIntegrator(const ROOT::Math::Integration::Type type, const ROOT::Math::Integration::GKRule rule, double absTol = 1.0000000000000001E-9, double relTol = 9.9999999999999995E-7, size_t size = 1000) |
ROOT::Math::GSLIntegrator | GSLIntegrator(const char* type, int rule, double absTol, double relTol, size_t size) |
virtual double | Integral() |
double | Integral(const ROOT::Math::IGenFunction& f) |
virtual double | Integral(const vector<double>& pts) |
double | Integral(const ROOT::Math::IGenFunction& f, const vector<double>& pts) |
virtual double | Integral(double a, double b) |
double | Integral(ROOT::Math::GSLFuncPointer f, void* p) |
double | Integral(const ROOT::Math::IGenFunction& f, double a, double b) |
double | Integral(ROOT::Math::GSLFuncPointer f, void* p, const vector<double>& pts) |
double | Integral(ROOT::Math::GSLFuncPointer f, void* p, double a, double b) |
virtual double | IntegralCauchy(double a, double b, double c) |
double | IntegralCauchy(const ROOT::Math::IGenFunction& f, double a, double b, double c) |
virtual double | IntegralLow(double b) |
double | IntegralLow(const ROOT::Math::IGenFunction& f, double b) |
double | IntegralLow(ROOT::Math::GSLFuncPointer f, void* p, double b) |
virtual double | IntegralUp(double a) |
double | IntegralUp(const ROOT::Math::IGenFunction& f, double a) |
double | IntegralUp(ROOT::Math::GSLFuncPointer f, void* p, double a) |
virtual int | NEval() const |
virtual ROOT::Math::IntegratorOneDimOptions | Options() const |
virtual double | Result() const |
virtual void | SetAbsTolerance(double absTolerance) |
virtual void | SetFunction(const ROOT::Math::IGenFunction& f) |
void | SetFunction(ROOT::Math::GSLFuncPointer f, void* p = 0) |
void | SetIntegrationRule(ROOT::Math::Integration::GKRule) |
virtual void | SetOptions(const ROOT::Math::IntegratorOneDimOptions& opt) |
virtual void | SetRelTolerance(double relTolerance) |
virtual int | Status() const |
virtual ROOT::Math::IntegrationOneDim::Type | ROOT::Math::VirtualIntegratorOneDim::Type() const |
ROOT::Math::VirtualIntegrator | ROOT::Math::VirtualIntegrator::VirtualIntegrator() |
ROOT::Math::VirtualIntegrator | ROOT::Math::VirtualIntegrator::VirtualIntegrator(const ROOT::Math::VirtualIntegrator&) |
ROOT::Math::VirtualIntegratorOneDim | ROOT::Math::VirtualIntegratorOneDim::VirtualIntegratorOneDim() |
ROOT::Math::VirtualIntegratorOneDim | ROOT::Math::VirtualIntegratorOneDim::VirtualIntegratorOneDim(const ROOT::Math::VirtualIntegratorOneDim&) |
bool | CheckFunction() |
ROOT::Math::GSLIntegrator | GSLIntegrator(const ROOT::Math::GSLIntegrator&) |
ROOT::Math::GSLIntegrator& | operator=(const ROOT::Math::GSLIntegrator&) |
constructor of GSL Integrator. In the case of Adaptive integration the Gauss-Krond rule of 31 points is used @param type type of integration. The possible types are defined in the Integration::Type enumeration @param absTol desired absolute Error @param relTol desired relative Error @param size maximum number of sub-intervals
generic constructor for GSL Integrator @param type type of integration. The possible types are defined in the Integration::Type enumeration @param rule Gauss-Kronrod rule. It is used only for ADAPTIVE::Integration types. The possible rules are defined in the Integration::GKRule enumeration @param absTol desired absolute Error @param relTol desired relative Error @param size maximum number of sub-intervals
constructor of GSL Integrator. In the case of Adaptive integration the Gauss-Krond rule of 31 points is used This is used by the plug-in manager (need a char * instead of enumerations) @param type type of integration. The possible types are defined in the Integration::Type enumeration @param rule Gauss-Kronrod rule (from 1 to 6) @param absTol desired absolute Error @param relTol desired relative Error @param size maximum number of sub-intervals
Set function from a GSL pointer function type
methods using IGenFunction evaluate the Integral of a function f over the defined interval (a,b) @param f integration function. The function type must implement the mathlib::IGenFunction interface @param a lower value of the integration interval @param b upper value of the integration interval
evaluate the Integral of a function f over the infinite interval (-inf,+inf) @param f integration function. The function type must implement the mathlib::IGenFunction interface
evaluate the Cauchy principal value of the integral of a previously defined function f over the defined interval (a,b) with a singularity at c @param a lower interval value @param b lower interval value @param c singular value of f
evaluate the Cauchy principal value of the integral of a function f over the defined interval (a,b) with a singularity at c @param f integration function. The function type must implement the mathlib::IGenFunction interface @param a lower interval value @param b lower interval value @param c singular value of f
evaluate the Integral of a function f over the semi-infinite interval (a,+inf) @param f integration function. The function type must implement the mathlib::IGenFunction interface @param a lower value of the integration interval
evaluate the Integral of a function f over the over the semi-infinite interval (-inf,b) @param f integration function. The function type must implement the mathlib::IGenFunction interface @param b upper value of the integration interval
evaluate the Integral of a function f with known singular points over the defined Integral (a,b) @param f integration function. The function type must implement the mathlib::IGenFunction interface @param pts vector containing both the function singular points and the lower/upper edges of the interval. The vector must have as first element the lower edge of the integration Integral ( \a a) and last element the upper value.
evaluate using cached function evaluate the Integral over the defined interval (a,b) using the function previously set with GSLIntegrator::SetFunction method @param a lower value of the integration interval @param b upper value of the integration interval
evaluate the Integral over the infinite interval (-inf,+inf) using the function previously set with GSLIntegrator::SetFunction method.
evaluate the Integral of a function f over the semi-infinite interval (a,+inf) using the function previously set with GSLIntegrator::SetFunction method. @param a lower value of the integration interval
evaluate the Integral of a function f over the over the semi-infinite interval (-inf,b) using the function previously set with GSLIntegrator::SetFunction method. @param b upper value of the integration interval
evaluate the Integral over the defined interval (a,b) using the function previously set with GSLIntegrator::SetFunction method. The function has known singular points. @param pts vector containing both the function singular points and the lower/upper edges of the interval. The vector must have as first element the lower edge of the integration Integral ( \a a) and last element the upper value.
evaluate using free function pointer (same GSL signature) signature for function pointers used by GSL typedef double ( * GSLFuncPointer ) ( double, void * ); evaluate the Integral of of a function f over the defined interval (a,b) passing a free function pointer The integration function must be a free function and have a signature consistent with GSL functions: <em>double my_function ( double x, void * p ) { ...... } </em> This method is the most efficient since no internal adapter to GSL function is created. @param f pointer to the integration function @param p pointer to the Parameters of the function @param a lower value of the integration interval @param b upper value of the integration interval
evaluate the Integral of a function f over the infinite interval (-inf,+inf) passing a free function pointer
evaluate the Integral of a function f over the semi-infinite interval (a,+inf) passing a free function pointer
evaluate the Integral of a function f over the over the semi-infinite interval (-inf,b) passing a free function pointer
evaluate the Integral of a function f with knows singular points over the over a defined interval passing a free function pointer
return number of function evaluations in calculating the integral
{ return fNEval; }
setter for control Parameters (getters are not needed so far ) set the desired relative Error
set the integration rule (Gauss-Kronrod rule). The possible rules are defined in the Integration::GKRule enumeration. The integration rule can be modified only for ADAPTIVE type integrations