00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016 #include "TMath.h"
00017 #include "TH1.h"
00018 #include "TCanvas.h"
00019 #include "TLegend.h"
00020
00021
00022 #include "Math/Functor.h"
00023 #include "Math/WrappedFunction.h"
00024 #include "Math/IFunction.h"
00025 #include "Math/Integrator.h"
00026 #include <iostream>
00027
00028 #include "TStopwatch.h"
00029 #include "TF1.h"
00030
00031 #include <limits>
00032
00033
00034
00035
00036 int nc = 0;
00037 double exactIntegral( double a, double b) {
00038
00039 return (TMath::ATan(2*b)- TMath::ATan(2*a))/ TMath::Pi();
00040 }
00041 double func( double x){
00042 nc++;
00043 return TMath::BreitWigner(x);
00044 }
00045
00046 double func2(const double *x, const double * = 0){
00047 nc++;
00048 return TMath::BreitWigner(x[0]);
00049 }
00050
00051
00052
00053
00054 void testIntegPerf(double x1, double x2, int n = 100000){
00055
00056
00057 std::cout << "\n\n***************************************************************\n";
00058 std::cout << "Test integration performances in interval [ " << x1 << " , " << x2 << " ]\n\n";
00059
00060 TStopwatch timer;
00061
00062 double dx = (x2-x1)/double(n);
00063
00064
00065 ROOT::Math::WrappedFunction<> f1(func);
00066
00067 timer.Start();
00068 ROOT::Math::Integrator ig(f1 );
00069 double s1 = 0.0;
00070 nc = 0;
00071 for (int i = 0; i < n; ++i) {
00072 double x = x1 + dx*i;
00073 s1+= ig.Integral(x1,x);
00074 }
00075 timer.Stop();
00076 std::cout << "Time using ROOT::Math::Integrator :\t" << timer.RealTime() << std::endl;
00077 std::cout << "Number of function calls = " << nc/n << std::endl;
00078 int pr = std::cout.precision(18); std::cout << s1 << std::endl; std::cout.precision(pr);
00079
00080
00081
00082
00083 TF1 *fBW = new TF1("fBW",func2,x1, x2,0);
00084
00085 timer.Start();
00086 nc = 0;
00087 double s2 = 0;
00088 for (int i = 0; i < n; ++i) {
00089 double x = x1 + dx*i;
00090 s2+= fBW->Integral(x1,x );
00091 }
00092 timer.Stop();
00093 std::cout << "Time using TF1::Integral :\t\t\t" << timer.RealTime() << std::endl;
00094 std::cout << "Number of function calls = " << nc/n << std::endl;
00095 pr = std::cout.precision(18); std::cout << s1 << std::endl; std::cout.precision(pr);
00096
00097
00098 }
00099
00100 void DrawCumulative(double x1, double x2, int n = 100){
00101
00102 std::cout << "\n\n***************************************************************\n";
00103 std::cout << "Drawing cumulatives of BreitWigner in interval [ " << x1 << " , " << x2 << " ]\n\n";
00104
00105
00106 double dx = (x2-x1)/double(n);
00107
00108 TH1D *cum0 = new TH1D("cum0", "", n, x1, x2);
00109 for (int i = 1; i <= n; ++i) {
00110 double x = x1 + dx*i;
00111 cum0->SetBinContent(i, exactIntegral(x1, x));
00112
00113 }
00114
00115
00116 ROOT::Math::Functor1D f1(& func);
00117
00118
00119 ROOT::Math::Integrator ig(f1, ROOT::Math::IntegrationOneDim::kADAPTIVE,1.E-12,1.E-12);
00120
00121 TH1D *cum1 = new TH1D("cum1", "", n, x1, x2);
00122 for (int i = 1; i <= n; ++i) {
00123 double x = x1 + dx*i;
00124 cum1->SetBinContent(i, ig.Integral(x1,x));
00125 }
00126
00127
00128 TF1 *fBW = new TF1("fBW","TMath::BreitWigner(x, 0, 1)",x1, x2);
00129
00130
00131 TH1D *cum2 = new TH1D("cum2", "", n, x1, x2);
00132 for (int i = 1; i <= n; ++i) {
00133 double x = x1 + dx*i;
00134 cum2->SetBinContent(i, fBW->Integral(x1,x));
00135 }
00136
00137 TH1D *cum10 = new TH1D("cum10", "", n, x1, x2);
00138 TH1D *cum20 = new TH1D("cum23", "", n, x1, x2);
00139 for (int i = 1; i <= n; ++i) {
00140 double delta = cum1->GetBinContent(i) - cum0->GetBinContent(i);
00141 double delta2 = cum2->GetBinContent(i) - cum0->GetBinContent(i);
00142
00143 cum10->SetBinContent(i, delta );
00144 cum10->SetBinError(i, std::numeric_limits<double>::epsilon() * cum1->GetBinContent(i) );
00145 cum20->SetBinContent(i, delta2 );
00146 }
00147
00148
00149 TCanvas *c1 = new TCanvas("c1","Integration example",20,10,800,500);
00150 c1->Divide(2,1);
00151 c1->Draw();
00152
00153 cum0->SetLineColor(kBlack);
00154 cum0->SetTitle("BreitWigner - the cumulative");
00155 cum0->SetStats(0);
00156 cum1->SetLineStyle(2);
00157 cum2->SetLineStyle(3);
00158 cum1->SetLineColor(kBlue);
00159 cum2->SetLineColor(kRed);
00160 c1->cd(1);
00161 cum0->DrawCopy("h");
00162 cum1->DrawCopy("same");
00163
00164 cum2->DrawCopy("same");
00165
00166 c1->cd(2);
00167 cum10->SetTitle("Difference");
00168 cum10->SetStats(0);
00169 cum10->SetLineColor(kBlue);
00170 cum10->Draw("e0");
00171 cum20->SetLineColor(kRed);
00172 cum20->Draw("hsame");
00173
00174 TLegend * l = new TLegend(0.11, 0.8, 0.7 ,0.89);
00175 l->AddEntry(cum10, "GSL integration - analytical ");
00176 l->AddEntry(cum20, "TF1::Integral - analytical ");
00177 l->Draw();
00178
00179
00180 c1->Update();
00181 std::cout << "\n***************************************************************\n";
00182
00183
00184 }
00185
00186
00187
00188 void mathmoreIntegration(double a = -2, double b = 2)
00189 {
00190 #if defined(__CINT__) && !defined(__MAKECINT__)
00191 cout << "WARNING: This tutorial can run only using ACliC, you must run it by doing: " << endl;
00192 cout << "\t .x $ROOTSYS/tutorials/math/mathmoreIntegration.C+" << endl;
00193 return;
00194 #endif
00195
00196 DrawCumulative(a, b);
00197 testIntegPerf(a, b);
00198 }
00199