Logo ROOT  
Reference Guide
GaussLegendreIntegrator.h
Go to the documentation of this file.
1// @(#)root/mathcore:$Id$
2// Authors: David Gonzalez Maline 01/2008
3
4/**********************************************************************
5 * *
6 * Copyright (c) 2006 , LCG ROOT MathLib Team *
7 * *
8 * *
9 **********************************************************************/
10
11// Header file for GaussIntegrator
12//
13// Created by: David Gonzalez Maline : Wed Jan 16 2008
14//
15
16#ifndef ROOT_Math_GaussLegendreIntegrator
17#define ROOT_Math_GaussLegendreIntegrator
18
19
21
22
23namespace ROOT {
24namespace Math {
25
26//___________________________________________________________________________________________
27/**
28 User class for performing function integration.
29
30 It will use the Gauss-Legendre Method for function integration in a given interval.
31 This class is implemented from TF1::Integral().
32
33 @ingroup Integration
34
35 */
36
38public:
39
40 /** Basic contructor of GaussLegendreIntegrator.
41 \@param num Number of desired points to calculate the integration.
42 \@param eps Desired relative error.
43 */
44 GaussLegendreIntegrator(int num = 10 ,double eps=1e-12);
45
46 /** Default Destructor */
48
49 /** Set the number of points used in the calculation of the
50 integral */
51 void SetNumberPoints(int num);
52
53 /** Set the desired relative Error. */
54 virtual void SetRelTolerance (double);
55
56 /** This method is not implemented. */
57 virtual void SetAbsTolerance (double);
58
59
60 /** Returns the arrays x and w containing the abscissa and weight of
61 the Gauss-Legendre n-point quadrature formula.
62
63 Gauss-Legendre: W(x)=1 -1<x<1
64 (j+1)P_{j+1} = (2j+1)xP_j-jP_{j-1}
65 */
66 void GetWeightVectors(double *x, double *w) const;
67
68 int GetNumberPoints() const { return fNum; }
69
70 /**
71 return number of function evaluations in calculating the integral
72 This is equivalent to the number of points
73 */
74 int NEval() const { return fNum; }
75
76
77 /// get the option used for the integration
79
80 // set the options
81 virtual void SetOptions(const ROOT::Math::IntegratorOneDimOptions & opt);
82
83private:
84
85 /**
86 Integration surrugate method. Return integral of passed function in interval [a,b]
87 Reimplement method of GaussIntegrator using CalcGaussLegendreSamplingPoints
88 */
89 virtual double DoIntegral (double a, double b, const IGenFunction* func);
90
91 /**
92 Type: unsafe but fast interface filling the arrays x and w (static method)
93
94 Given the number of sampling points this routine fills the arrays x and w
95 of length num, containing the abscissa and weight of the Gauss-Legendre
96 n-point quadrature formula.
97
98 Gauss-Legendre: W(x)=1 -1<x<1
99 (j+1)P_{j+1} = (2j+1)xP_j-jP_{j-1}
100
101 num is the number of sampling points (>0)
102 x and w are arrays of size num
103 eps is the relative precision
104
105 If num<=0 or eps<=0 no action is done.
106
107 Reference: Numerical Recipes in C, Second Edition
108 */
110
111
112protected:
113 int fNum; // Number of points used in the stimation of the integral.
114 double* fX; // Abscisa of the points used.
115 double* fW; // Weights of the points used.
116
117};
118
119} // end namespace Math
120
121} // end namespace ROOT
122
123#endif /* ROOT_Math_GaussLegendreIntegrator */
#define b(i)
Definition: RSha256.hxx:100
#define e(i)
Definition: RSha256.hxx:103
User class for performing function integration.
User class for performing function integration.
virtual double DoIntegral(double a, double b, const IGenFunction *func)
Integration surrugate method.
virtual void SetAbsTolerance(double)
This method is not implemented.
virtual ~GaussLegendreIntegrator()
Default Destructor.
virtual void SetOptions(const ROOT::Math::IntegratorOneDimOptions &opt)
set the options (should be re-implemented by derived classes -if more options than tolerance exist
void GetWeightVectors(double *x, double *w) const
Returns the arrays x and w containing the abscissa and weight of the Gauss-Legendre n-point quadratur...
int NEval() const
return number of function evaluations in calculating the integral This is equivalent to the number of...
GaussLegendreIntegrator(int num=10, double eps=1e-12)
Basic contructor of GaussLegendreIntegrator.
void CalcGaussLegendreSamplingPoints()
Type: unsafe but fast interface filling the arrays x and w (static method)
virtual void SetRelTolerance(double)
Set the desired relative Error.
void SetNumberPoints(int num)
Set the number of points used in the calculation of the integral.
virtual ROOT::Math::IntegratorOneDimOptions Options() const
get the option used for the integration
Interface (abstract class) for generic functions objects of one-dimension Provides a method to evalua...
Definition: IFunction.h:135
Numerical one dimensional integration options.
Double_t x[n]
Definition: legend1.C:17
Namespace for new Math classes and functions.
tbb::task_arena is an alias of tbb::interface7::task_arena, which doesn't allow to forward declare tb...
Definition: StringConv.hxx:21
auto * a
Definition: textangle.C:12