Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
vectorizedFit.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_fit
3 /// \notebook /// Tutorial for creating a Vectorized TF1 function using a formunla expression and
4 /// use it for fitting an histogram
5 ///
6 /// To create a vectorized function (if ROOT has been compiled with support for vectorization)
7 /// is very easy. One needs to create the TF1 object with the option "VEC" or call the method
8 /// TF1::SetVectorized
9 ///
10 /// \macro_image
11 /// \macro output
12 /// \macro_code
13 ///
14 /// \author Lorenzo Moneta
15 
16 void vectorizedFit() {
17 
18  gStyle->SetOptFit(111111);
19 
20 
21  ROOT::Math::MinimizerOptions::SetDefaultMinimizer("Minuit2");
22 
23  int nbins = 40000;
24  auto h1 = new TH1D("h1","h1",nbins,-3,3);
25  h1->FillRandom("gaus",nbins*50);
26  auto c1 = new TCanvas("Fit","Fit",800,1000);
27  c1->Divide(1,2);
28  c1->cd(1);
29  TStopwatch w;
30 
31  std::cout << "Doing Serial Gaussian Fit " << std::endl;
32  auto f1 = new TF1("f1","gaus");
33  f1->SetNpx(nbins*10);
34  w.Start();
35  h1->Fit(f1);
36  h1->Fit(f1,"L+");
37  w.Print();
38 
39  std::cout << "Doing Vectorized Gaussian Fit " << std::endl;
40  auto f2 = new TF1("f2","gaus",-3,3,"VEC");
41  // alternativly you can also use the TF1::SetVectorized function
42  //f2->SetVectorized(true);
43  w.Start();
44  h1->Fit(f2);
45  h1->Fit(f2,"L+");
46  w.Print();
47  // rebin histograms and scale it back to the function
48  h1->Rebin(nbins/100);
49  h1->Scale(100./nbins);
50  ((TF1 *)h1->GetListOfFunctions()->At(0))->SetTitle("Chi2 Fit");
51  ((TF1 *)h1->GetListOfFunctions()->At(1))->SetTitle("Likelihood Fit");
52  ((TF1 *)h1->GetListOfFunctions()->At(1))->SetLineColor(kBlue);
53  //c1->cd(1)->BuildLegend();
54 
55  /// Do a polynomail fit now
56  c1->cd(2);
57  auto f3 = new TF1("f3","[A]*x^2+[B]*x+[C]",0,10);
58  f3->SetParameters(0.5,3,2);
59  f3->SetNpx(nbins*10);
60  // generate the events
61  auto h2 = new TH1D("h2","h2",nbins,0,10);
62  h2->FillRandom("f3",10*nbins);
63  std::cout << "Doing Serial Polynomial Fit " << std::endl;
64  f3->SetParameters(2,2,2);
65  w.Start();
66  h2->Fit(f3);
67  h2->Fit(f3,"L+");
68  w.Print();
69 
70  std::cout << "Doing Vectorized Polynomial Fit " << std::endl;
71  auto f4 = new TF1("f4","[A]*x*x+[B]*x+[C]",0,10);
72  f4->SetVectorized(true);
73  f4->SetParameters(2,2,2);
74  w.Start();
75  h2->Fit(f4);
76  h2->Fit(f4,"L+");
77  w.Print();
78 
79  // rebin histograms and scale it back to the function
80  h2->Rebin(nbins/100);
81  h2->Scale(100./nbins);
82  ((TF1 *)h2->GetListOfFunctions()->At(0))->SetTitle("Chi2 Fit");
83  ((TF1 *)h2->GetListOfFunctions()->At(1))->SetTitle("Likelihood Fit");
84  ((TF1 *)h2->GetListOfFunctions()->At(1))->SetLineColor(kBlue);
85  //c1->cd(2)->BuildLegend();
86 
87 }