Logo ROOT  
Reference Guide
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
Loading...
Searching...
No Matches
mathmoreIntegration.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_math
3/// \notebook -nodraw
4/// Example on the usage of the adaptive 1D integration algorithm of MathMore.
5///
6/// It calculates the numerically cumulative integral of a distribution (like in this case the BreitWigner)
7/// to execute the macro type it (you need to compile with AClic)
8///
9/// ~~~{.cpp}
10/// root[0] .x mathmoreIntegration.C+
11/// ~~~
12///
13/// This tutorial requires having libMathMore built with ROOT.
14///
15/// To build mathmore you need to have a version of GSL >= 1.8 installed in your system
16/// The ROOT configure will automatically find GSL if the script gsl-config (from GSL) is in your PATH,.
17/// otherwise you need to configure root with the options --gsl-incdir and --gsl-libdir.
18///
19/// \macro_image
20/// \macro_output
21/// \macro_code
22///
23/// \authors M. Slawinska, L. Moneta
24
25#include "TMath.h"
26#include "TH1.h"
27#include "TCanvas.h"
28#include "TLegend.h"
29
30/*#include "TLabel.h"*/
31#include "Math/Functor.h"
33#include "Math/IFunction.h"
34#include "Math/Integrator.h"
35#include <iostream>
36
37#include "TStopwatch.h"
38#include "TF1.h"
39
40#include <limits>
41
42//!calculates exact integral of Breit Wigner distribution
43//!and compares with existing methods
44
45int nc = 0;
46double exactIntegral( double a, double b) {
47
48 return (TMath::ATan(2*b)- TMath::ATan(2*a))/ TMath::Pi();
49}
50double func( double x){
51 nc++;
52 return TMath::BreitWigner(x);
53}
54
55// TF1 requires the function to have the ( )( double *, double *) signature
56double func2(const double *x, const double * = nullptr){
57 nc++;
58 return TMath::BreitWigner(x[0]);
59}
60
61
62
63
64void testIntegPerf(double x1, double x2, int n = 100000){
65
66
67 std::cout << "\n\n***************************************************************\n";
68 std::cout << "Test integration performances in interval [ " << x1 << " , " << x2 << " ]\n\n";
69
71
72 double dx = (x2-x1)/double(n);
73
74 //ROOT::Math::Functor1D<ROOT::Math::IGenFunction> f1(& TMath::BreitWigner);
76
77 timer.Start();
79 double s1 = 0.0;
80 nc = 0;
81 for (int i = 0; i < n; ++i) {
82 double x = x1 + dx*i;
83 s1+= ig.Integral(x1,x);
84 }
85 timer.Stop();
86 std::cout << "Time using ROOT::Math::Integrator :\t" << timer.RealTime() << std::endl;
87 std::cout << "Number of function calls = " << nc/n << std::endl;
88 int pr = std::cout.precision(18); std::cout << s1 << std::endl; std::cout.precision(pr);
89
90
91
92 //TF1 *fBW = new TF1("fBW","TMath::BreitWigner(x)",x1, x2); // this is faster but cannot measure number of function calls
93 TF1 *fBW = new TF1("fBW",func2,x1, x2,0);
94
95 timer.Start();
96 nc = 0;
97 double s2 = 0;
98 for (int i = 0; i < n; ++i) {
99 double x = x1 + dx*i;
100 s2+= fBW->Integral(x1,x );
101 }
102 timer.Stop();
103 std::cout << "Time using TF1::Integral :\t\t\t" << timer.RealTime() << std::endl;
104 std::cout << "Number of function calls = " << nc/n << std::endl;
105 pr = std::cout.precision(18); std::cout << s1 << std::endl; std::cout.precision(pr);
106
107
108}
109
110void DrawCumulative(double x1, double x2, int n = 100){
111
112 std::cout << "\n\n***************************************************************\n";
113 std::cout << "Drawing cumulatives of BreitWigner in interval [ " << x1 << " , " << x2 << " ]\n\n";
114
115
116 double dx = (x2-x1)/double(n);
117
118 TH1D *cum0 = new TH1D("cum0", "", n, x1, x2); //exact cumulative
119 for (int i = 1; i <= n; ++i) {
120 double x = x1 + dx*i;
121 cum0->SetBinContent(i, exactIntegral(x1, x));
122
123 }
124
125 // alternative method using ROOT::Math::Functor class
127
128
130
131 TH1D *cum1 = new TH1D("cum1", "", n, x1, x2);
132 for (int i = 1; i <= n; ++i) {
133 double x = x1 + dx*i;
134 cum1->SetBinContent(i, ig.Integral(x1,x));
135 }
136
137
138 TF1 *fBW = new TF1("fBW","TMath::BreitWigner(x, 0, 1)",x1, x2);
139
140
141 TH1D *cum2 = new TH1D("cum2", "", n, x1, x2);
142 for (int i = 1; i <= n; ++i) {
143 double x = x1 + dx*i;
144 cum2->SetBinContent(i, fBW->Integral(x1,x));
145 }
146
147 TH1D *cum10 = new TH1D("cum10", "", n, x1, x2); //difference between 1 and exact
148 TH1D *cum20 = new TH1D("cum23", "", n, x1, x2); //difference between 2 and exact
149 for (int i = 1; i <= n; ++i) {
150 double delta = cum1->GetBinContent(i) - cum0->GetBinContent(i);
151 double delta2 = cum2->GetBinContent(i) - cum0->GetBinContent(i);
152 //std::cout << " diff for " << x << " is " << delta << " " << cum1->GetBinContent(i) << std::endl;
153 cum10->SetBinContent(i, delta );
154 cum10->SetBinError(i, std::numeric_limits<double>::epsilon() * cum1->GetBinContent(i) );
155 cum20->SetBinContent(i, delta2 );
156 }
157
158
159 TCanvas *c1 = new TCanvas("c1","Integration example",20,10,800,500);
160 c1->Divide(2,1);
161 c1->Draw();
162
163 cum0->SetLineColor(kBlack);
164 cum0->SetTitle("BreitWigner - the cumulative");
165 cum0->SetStats(false);
166 cum1->SetLineStyle(kDashed);
167 cum2->SetLineStyle(kDotted);
168 cum1->SetLineColor(kBlue);
169 cum2->SetLineColor(kRed);
170 c1->cd(1);
171 cum0->DrawCopy("h");
172 cum1->DrawCopy("same");
173 //cum2->DrawCopy("same");
174 cum2->DrawCopy("same");
175
176 c1->cd(2);
177 cum10->SetTitle("Difference");
178 cum10->SetStats(false);
179 cum10->SetLineColor(kBlue);
180 cum10->Draw("e0");
181 cum20->SetLineColor(kRed);
182 cum20->Draw("hsame");
183
184 TLegend * l = new TLegend(0.11, 0.8, 0.7 ,0.89);
185 l->AddEntry(cum10, "GSL integration - analytical ");
186 l->AddEntry(cum20, "TF1::Integral - analytical ");
187 l->Draw();
188
189
190 c1->Update();
191 std::cout << "\n***************************************************************\n";
192
193
194}
195
196
197
198void mathmoreIntegration(double a = -2, double b = 2)
199{
201 testIntegPerf(a, b);
202}
203
#define b(i)
Definition RSha256.hxx:100
#define a(i)
Definition RSha256.hxx:99
#define s1(x)
Definition RSha256.hxx:91
@ kRed
Definition Rtypes.h:66
@ kBlack
Definition Rtypes.h:65
@ kBlue
Definition Rtypes.h:66
@ kDashed
Definition TAttLine.h:52
@ kDotted
Definition TAttLine.h:52
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
Option_t Option_t TPoint TPoint const char x2
Option_t Option_t TPoint TPoint const char x1
Functor1D class for one-dimensional functions.
Definition Functor.h:95
User Class for performing numerical integration of a function in one dimension.
Definition Integrator.h:98
Template class to wrap any C++ callable object which takes one argument i.e.
The Canvas class.
Definition TCanvas.h:23
1-Dim function class
Definition TF1.h:234
1-D histogram with a double per channel (see TH1 documentation)
Definition TH1.h:698
This class displays a legend box (TPaveText) containing several legend entries.
Definition TLegend.h:23
virtual void Draw(Option_t *option="")
Default Draw method for all objects.
Definition TObject.cxx:293
Stopwatch class.
Definition TStopwatch.h:28
@ kADAPTIVE
to be used for general functions without singularities
return c1
Definition legend1.C:41
Double_t x[n]
Definition legend1.C:17
const Int_t n
Definition legend1.C:16
TF1 * f1
Definition legend1.C:11
Double_t ATan(Double_t)
Returns the principal value of the arc tangent of x, expressed in radians.
Definition TMath.h:644
Double_t BreitWigner(Double_t x, Double_t mean=0, Double_t gamma=1)
Calculates a Breit Wigner function with mean and gamma.
Definition TMath.cxx:442
constexpr Double_t Pi()
Definition TMath.h:37
TLine l
Definition textangle.C:4