Logo ROOT   6.08/07
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 
20 #ifndef ROOT_Math_GaussIntegrator
21 #include "Math/GaussIntegrator.h"
22 #endif
23 
24 
25 namespace ROOT {
26 namespace Math {
27 
28 //___________________________________________________________________________________________
29 /**
30  User class for performing function integration.
31 
32  It will use the Gauss-Legendre Method for function integration in a given interval.
33  This class is implemented from TF1::Integral().
34 
35  @ingroup Integration
36 
37  */
38 
40 public:
41 
42  /** Basic contructor of GaussLegendreIntegrator.
43  \@param num Number of desired points to calculate the integration.
44  \@param eps Desired relative error.
45  */
46  GaussLegendreIntegrator(int num = 10 ,double eps=1e-12);
47 
48  /** Default Destructor */
49  virtual ~GaussLegendreIntegrator();
50 
51  /** Set the number of points used in the calculation of the
52  integral */
53  void SetNumberPoints(int num);
54 
55  /** Set the desired relative Error. */
56  virtual void SetRelTolerance (double);
57 
58  /** This method is not implemented. */
59  virtual void SetAbsTolerance (double);
60 
61 
62  /** Returns the arrays x and w containing the abscissa and weight of
63  the Gauss-Legendre n-point quadrature formula.
64 
65  Gauss-Legendre: W(x)=1 -1<x<1
66  (j+1)P_{j+1} = (2j+1)xP_j-jP_{j-1}
67  */
68  void GetWeightVectors(double *x, double *w) const;
69 
70  int GetNumberPoints() const { return fNum; }
71 
72  /**
73  return number of function evaluations in calculating the integral
74  This is equivalent to the number of points
75  */
76  int NEval() const { return fNum; }
77 
78 
79  /// get the option used for the integration
81 
82  // set the options
83  virtual void SetOptions(const ROOT::Math::IntegratorOneDimOptions & opt);
84 
85 private:
86 
87  /**
88  Integration surrugate method. Return integral of passed function in interval [a,b]
89  Reimplement method of GaussIntegrator using CalcGaussLegendreSamplingPoints
90  */
91  virtual double DoIntegral (double a, double b, const IGenFunction* func);
92 
93  /**
94  Type: unsafe but fast interface filling the arrays x and w (static method)
95 
96  Given the number of sampling points this routine fills the arrays x and w
97  of length num, containing the abscissa and weight of the Gauss-Legendre
98  n-point quadrature formula.
99 
100  Gauss-Legendre: W(x)=1 -1<x<1
101  (j+1)P_{j+1} = (2j+1)xP_j-jP_{j-1}
102 
103  num is the number of sampling points (>0)
104  x and w are arrays of size num
105  eps is the relative precision
106 
107  If num<=0 or eps<=0 no action is done.
108 
109  Reference: Numerical Recipes in C, Second Edition
110  */
112 
113 
114 protected:
115  int fNum; // Number of points used in the stimation of the integral.
116  double* fX; // Abscisa of the points used.
117  double* fW; // Weights of the points used.
118 
119 };
120 
121 } // end namespace Math
122 
123 } // end namespace ROOT
124 
125 #endif /* ROOT_Math_GaussLegendreIntegrator */
Interface (abstract class) for generic functions objects of one-dimension Provides a method to evalua...
Definition: IFunction.h:133
This namespace contains pre-defined functions to be used in conjuction with TExecutor::Map and TExecu...
Definition: StringConv.hxx:21
virtual ROOT::Math::IntegratorOneDimOptions Options() const
get the option used for the integration
GaussLegendreIntegrator(int num=10, double eps=1e-12)
Basic contructor of GaussLegendreIntegrator.
TArc * a
Definition: textangle.C:12
void CalcGaussLegendreSamplingPoints()
Type: unsafe but fast interface filling the arrays x and w (static method)
User class for performing function integration.
int NEval() const
return number of function evaluations in calculating the integral This is equivalent to the number of...
virtual double DoIntegral(double a, double b, const IGenFunction *func)
Integration surrugate method.
Double_t x[n]
Definition: legend1.C:17
virtual void SetOptions(const ROOT::Math::IntegratorOneDimOptions &opt)
set the options (should be re-implemented by derived classes -if more options than tolerance exist ...
virtual ~GaussLegendreIntegrator()
Default Destructor.
virtual void SetRelTolerance(double)
Set the desired relative Error.
Numerical one dimensional integration options.
User class for performing function integration.
void SetNumberPoints(int num)
Set the number of points used in the calculation of the integral.
double func(double *x, double *p)
Definition: stressTF1.cxx:213
you should not use this method at all Int_t Int_t Double_t Double_t Double_t e
Definition: TRolke.cxx:630
Namespace for new Math classes and functions.
virtual void SetAbsTolerance(double)
This method is not implemented.
you should not use this method at all Int_t Int_t Double_t Double_t Double_t Int_t Double_t Double_t Double_t Double_t b
Definition: TRolke.cxx:630
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...