Logo ROOT   6.10/09
Reference Guide
mp_H1_lambdas.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_multicore
3 /// \notebook -nodraw
4 /// Lambdas used to check and fit the result of the H1 analysis.
5 /// Used by mp104_processH1.C, mp105_processEntryList.C and roottest/root/multicore/tProcessExecutorH1Test.C
6 ///
7 /// \macro_code
8 ///
9 /// \author Gerardo Ganis
10 
11 // This function is used to check the result of the H1 analysis
12 auto checkH1 = [](TList *out) {
13 
14  // Make sure the output list is there
15  if (!out) {
16  std::cout << "checkH1 >>> Test failure: output list not found\n";
17  return -1;
18  }
19 
20  // Check the 'hdmd' histo
21  TH1F *hdmd = dynamic_cast<TH1F *>(out->FindObject("hdmd"));
22  if (!hdmd) {
23  std::cout << "checkH1 >>> Test failure: 'hdmd' histo not found\n";
24  return -1;
25  }
26  if ((Int_t)(hdmd->GetEntries()) != 7525) {
27  std::cout << "checkH1 >>> Test failure: 'hdmd' histo: wrong number"
28  " of entries ("
29  << (Int_t)(hdmd->GetEntries()) << ": expected 7525) \n";
30  return -1;
31  }
32  if (TMath::Abs((hdmd->GetMean() - 0.15512023) / 0.15512023) > 0.001) {
33  std::cout << "checkH1 >>> Test failure: 'hdmd' histo: wrong mean (" << hdmd->GetMean()
34  << ": expected 0.15512023) \n";
35  return -1;
36  }
37 
38  TH2F *h2 = dynamic_cast<TH2F *>(out->FindObject("h2"));
39  if (!h2) {
40  std::cout << "checkH1 >>> Test failure: 'h2' histo not found\n";
41  return -1;
42  }
43  if ((Int_t)(h2->GetEntries()) != 7525) {
44  std::cout << "checkH1 >>> Test failure: 'h2' histo: wrong number"
45  " of entries ("
46  << (Int_t)(h2->GetEntries()) << ": expected 7525) \n";
47  return -1;
48  }
49  if (TMath::Abs((h2->GetMean() - 0.15245688) / 0.15245688) > 0.001) {
50  std::cout << "checkH1 >>> Test failure: 'h2' histo: wrong mean (" << h2->GetMean() << ": expected 0.15245688) \n";
51  return -1;
52  }
53 
54  // Done
55  return 0;
56 };
57 
58 // This function is used to fit the result of the analysis with graphics
59 auto doFit = [](TList *out, const char *lfn = 0) -> Int_t {
60 
61  RedirectHandle_t redH;
62  if (lfn) gSystem->RedirectOutput(lfn, "a", &redH);
63 
64  auto hdmd = dynamic_cast<TH1F *>(out->FindObject("hdmd"));
65  auto h2 = dynamic_cast<TH2F *>(out->FindObject("h2"));
66 
67  // function called at the end of the event loop
68  if (hdmd == 0 || h2 == 0) {
69  std::cout << "doFit: hdmd = " << hdmd << " , h2 = " << h2 << "\n";
70  return -1;
71  if (lfn) gSystem->RedirectOutput(0, 0, &redH);
72  }
73 
74  // create the canvas for the h1analysis fit
75  gStyle->SetOptFit();
76  TCanvas *c1 = new TCanvas("c1", "h1analysis analysis", 10, 10, 800, 600);
77  c1->SetBottomMargin(0.15);
78  hdmd->GetXaxis()->SetTitle("m_{K#pi#pi} - m_{K#pi}[GeV/c^{2}]");
79  hdmd->GetXaxis()->SetTitleOffset(1.4);
80 
81  // fit histogram hdmd with function f5 using the log-likelihood option
82  if (gROOT->GetListOfFunctions()->FindObject("f5")) delete gROOT->GetFunction("f5");
83 
84  auto fdm5 = [](Double_t *xx, Double_t *par) -> Double_t {
85  const Double_t dxbin = (0.17 - 0.13) / 40; // Bin-width
86  Double_t x = xx[0];
87  if (x <= 0.13957) return 0;
88  Double_t xp3 = (x - par[3]) * (x - par[3]);
89  Double_t res = dxbin * (par[0] * TMath::Power(x - 0.13957, par[1]) +
90  par[2] / 2.5066 / par[4] * TMath::Exp(-xp3 / 2 / par[4] / par[4]));
91  return res;
92  };
93 
94  TF1 *f5 = new TF1("f5", fdm5, 0.139, 0.17, 5);
95  f5->SetParameters(1000000, .25, 2000, .1454, .001);
96  hdmd->Fit("f5", "lr");
97 
98  // Check the result of the fit
99  Double_t ref_f5[4] = {959915.0, 0.351114, 1185.03, 0.145569};
100  for (int i : {0, 1, 2, 3}) {
101  if ((TMath::Abs((f5->GetParameters())[i] - ref_f5[i]) / ref_f5[i]) > 0.001) {
102  std::cout << "\n >>> Test failure: fit to 'f5': parameter '" << f5->GetParName(i) << "' has wrong value ("
103  << (f5->GetParameters())[i] << ": expected" << ref_f5[i] << ") \n";
104  if (lfn) gSystem->RedirectOutput(0, 0, &redH);
105  return -1;
106  }
107  }
108 
109  // create the canvas for tau d0
110  gStyle->SetOptFit(0);
111  gStyle->SetOptStat(1100);
112  TCanvas *c2 = new TCanvas("c2", "tauD0", 100, 100, 800, 600);
113  c2->SetGrid();
114  c2->SetBottomMargin(0.15);
115 
116  // Project slices of 2-d histogram h2 along X , then fit each slice
117  // with function f2 and make a histogram for each fit parameter
118  // Note that the generated histograms are added to the list of objects
119  // in the current directory.
120  if (gROOT->GetListOfFunctions()->FindObject("f2")) delete gROOT->GetFunction("f2");
121 
122  auto fdm2 = [](Double_t *xx, Double_t *par) -> Double_t {
123  const Double_t dxbin = (0.17 - 0.13) / 40; // Bin-width
124  const Double_t sigma = 0.0012;
125  Double_t x = xx[0];
126  if (x <= 0.13957) return 0;
127  Double_t xp3 = (x - 0.1454) * (x - 0.1454);
128  Double_t res = dxbin * (par[0] * TMath::Power(x - 0.13957, 0.25) +
129  par[1] / 2.5066 / sigma * TMath::Exp(-xp3 / 2 / sigma / sigma));
130  return res;
131  };
132 
133  TF1 *f2 = new TF1("f2", fdm2, 0.139, 0.17, 2);
134  f2->SetParameters(10000, 10);
135 
136  // Restrict to three bins in this example
137  std::cout << "doFit: restricting fit to two bins only in this example...\n";
138 
139  h2->FitSlicesX(f2, 10, 20, 10, "g5 l");
140 
141  // Check the result of the fit
142  Double_t ref_f2[2] = {52432.2, 105.481};
143  for (int i : {0, 1}) {
144  if ((TMath::Abs((f2->GetParameters())[i] - ref_f2[i]) / ref_f2[i]) > 0.001) {
145  std::cout << "\n >>> Test failure: fit to 'f2': parameter '" << f2->GetParName(i) << "' has wrong value ("
146  << (f2->GetParameters())[i] << ": expected" << ref_f2[i] << ") \n";
147  if (lfn) gSystem->RedirectOutput(0, 0, &redH);
148  return -1;
149  }
150  }
151 
152  TH1D *h2_1 = (TH1D *)gDirectory->Get("h2_1");
153  h2_1->GetXaxis()->SetTitle("#tau[ps]");
154  h2_1->SetMarkerStyle(21);
155  h2_1->Draw();
156  c2->Update();
157  TLine *line = new TLine(0, 0, 0, c2->GetUymax());
158  line->Draw();
159 
160  // Have the number of entries on the first histogram (to cross check when running
161  // with entry lists)
162  TPaveStats *psdmd = (TPaveStats *)hdmd->GetListOfFunctions()->FindObject("stats");
163  psdmd->SetOptStat(1110);
164  c1->Modified();
165 
166  if (lfn) gSystem->RedirectOutput(0, 0, &redH);
167 
168  return 0;
169 };
170 
171 // This is the function invoked during the processing of the trees.
172 auto doH1 = [](TTreeReader &reader) {
173 
174  // Histograms
175  auto hdmd = new TH1F("hdmd", "Dm_d", 40, 0.13, 0.17);
176  auto h2 = new TH2F("h2", "ptD0 vs Dm_d", 30, 0.135, 0.165, 30, -3, 6);
177 
178  TTreeReaderValue<Float_t> fPtds_d(reader, "ptds_d");
179  TTreeReaderValue<Float_t> fEtads_d(reader, "etads_d");
180  TTreeReaderValue<Float_t> fDm_d(reader, "dm_d");
181  TTreeReaderValue<Int_t> fIk(reader, "ik");
182  TTreeReaderValue<Int_t> fIpi(reader, "ipi");
183  TTreeReaderValue<Int_t> fIpis(reader, "ipis");
184  TTreeReaderValue<Float_t> fPtd0_d(reader, "ptd0_d");
185  TTreeReaderValue<Float_t> fMd0_d(reader, "md0_d");
186  TTreeReaderValue<Float_t> fRpd0_t(reader, "rpd0_t");
187  TTreeReaderArray<Int_t> fNhitrp(reader, "nhitrp");
188  TTreeReaderArray<Float_t> fRstart(reader, "rstart");
189  TTreeReaderArray<Float_t> fRend(reader, "rend");
190  TTreeReaderArray<Float_t> fNlhk(reader, "nlhk");
191  TTreeReaderArray<Float_t> fNlhpi(reader, "nlhpi");
192  TTreeReaderValue<Int_t> fNjets(reader, "njets");
193 
194  while (reader.Next()) {
195 
196  // Return as soon as a bad entry is detected
197  if (TMath::Abs(*fMd0_d - 1.8646) >= 0.04) continue;
198  if (*fPtds_d <= 2.5) continue;
199  if (TMath::Abs(*fEtads_d) >= 1.5) continue;
200  (*fIk)--; // original fIk used f77 convention starting at 1
201  (*fIpi)--;
202 
203  if (fNhitrp.At(*fIk) * fNhitrp.At(*fIpi) <= 1) continue;
204 
205  if (fRend.At(*fIk) - fRstart.At(*fIk) <= 22) continue;
206  if (fRend.At(*fIpi) - fRstart.At(*fIpi) <= 22) continue;
207  if (fNlhk.At(*fIk) <= 0.1) continue;
208  if (fNlhpi.At(*fIpi) <= 0.1) continue;
209  (*fIpis)--;
210  if (fNlhpi.At(*fIpis) <= 0.1) continue;
211  if (*fNjets < 1) continue;
212 
213  // Fill the histograms
214  hdmd->Fill(*fDm_d);
215  h2->Fill(*fDm_d, *fRpd0_t / 0.029979 * 1.8646 / *fPtd0_d);
216  }
217 
218  // Return a list
219  auto l = new TList;
220  l->Add(hdmd);
221  l->Add(h2);
222  l->SetOwner(kFALSE);
223 
224  return l;
225 };
226 
227 // This is the function invoked during the processing of the trees to create a TEntryList
228 auto doH1fillList = [](TTreeReader &reader) {
229 
230  // Entry list
231  auto elist = new TEntryList("elist", "H1 selection from Cut");
232 
233  TTreeReaderValue<Float_t> fPtds_d(reader, "ptds_d");
234  TTreeReaderValue<Float_t> fEtads_d(reader, "etads_d");
235  TTreeReaderValue<Int_t> fIk(reader, "ik");
236  TTreeReaderValue<Int_t> fIpi(reader, "ipi");
237  TTreeReaderValue<Int_t> fIpis(reader, "ipis");
238  TTreeReaderValue<Float_t> fMd0_d(reader, "md0_d");
239  TTreeReaderArray<Int_t> fNhitrp(reader, "nhitrp");
240  TTreeReaderArray<Float_t> fRstart(reader, "rstart");
241  TTreeReaderArray<Float_t> fRend(reader, "rend");
242  TTreeReaderArray<Float_t> fNlhk(reader, "nlhk");
243  TTreeReaderArray<Float_t> fNlhpi(reader, "nlhpi");
244  TTreeReaderValue<Int_t> fNjets(reader, "njets");
245 
246  while (reader.Next()) {
247 
248  // Return as soon as a bad entry is detected
249  if (TMath::Abs(*fMd0_d - 1.8646) >= 0.04) continue;
250  if (*fPtds_d <= 2.5) continue;
251  if (TMath::Abs(*fEtads_d) >= 1.5) continue;
252  (*fIk)--; // original fIk used f77 convention starting at 1
253  (*fIpi)--;
254 
255  if (fNhitrp.At(*fIk) * fNhitrp.At(*fIpi) <= 1) continue;
256 
257  if (fRend.At(*fIk) - fRstart.At(*fIk) <= 22) continue;
258  if (fRend.At(*fIpi) - fRstart.At(*fIpi) <= 22) continue;
259  if (fNlhk.At(*fIk) <= 0.1) continue;
260  if (fNlhpi.At(*fIpi) <= 0.1) continue;
261  (*fIpis)--;
262  if (fNlhpi.At(*fIpis) <= 0.1) continue;
263  if (*fNjets < 1) continue;
264 
265  // Fill the entry list
266  elist->Enter(reader.GetCurrentEntry(), reader.GetTree());
267  }
268 
269  return elist;
270 };
271 
272 // This is the function invoked during the processing of the trees using a TEntryList
273 auto doH1useList = [](TTreeReader &reader) {
274 
275  // Histograms
276  auto hdmd = new TH1F("hdmd", "Dm_d", 40, 0.13, 0.17);
277  auto h2 = new TH2F("h2", "ptD0 vs Dm_d", 30, 0.135, 0.165, 30, -3, 6);
278 
279  TTreeReaderValue<Float_t> fDm_d(reader, "dm_d");
280  TTreeReaderValue<Float_t> fPtd0_d(reader, "ptd0_d");
281  TTreeReaderValue<Float_t> fRpd0_t(reader, "rpd0_t");
282 
283  while (reader.Next()) {
284  // Fill the histograms
285  hdmd->Fill(*fDm_d);
286  h2->Fill(*fDm_d, *fRpd0_t / 0.029979 * 1.8646 / *fPtd0_d);
287  }
288 
289  // Return a list
290  auto l = new TList;
291  l->Add(hdmd);
292  l->Add(h2);
293  l->SetOwner(kFALSE);
294 
295  return l;
296 };
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title Offset is a correction factor with respect to the "s...
Definition: TAttAxis.cxx:262
double par[1]
Definition: unuranDistr.cxx:38
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition: TH1.cxx:3126
virtual void SetParameters(const Double_t *params)
Definition: TF1.h:588
void doFit(int n, const char *fitter)
TTreeReader is a simple, robust and fast interface to read values from a TTree, TChain or TNtuple...
Definition: TTreeReader.h:43
TLine * line
return c1
Definition: legend1.C:41
R__EXTERN TStyle * gStyle
Definition: TStyle.h:402
THist< 1, float, THistStatContent, THistStatUncertainty > TH1F
Definition: THist.hxx:311
Double_t fdm2(Double_t *xx, Double_t *par)
#define gROOT
Definition: TROOT.h:375
virtual Double_t GetMean(Int_t axis=1) const
For axis = 1,2 or 3 returns the mean value of the histogram along X,Y or Z axis.
Definition: TH1.cxx:6763
tomato 1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:551
int Int_t
Definition: RtypesCore.h:41
virtual void Draw(Option_t *option="")
Default Draw method for all objects.
Definition: TObject.cxx:202
The histogram statistics painter class.
Definition: TPaveStats.h:18
Short_t Abs(Short_t d)
Definition: TMathBase.h:108
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Definition: TMath.h:628
virtual TObject * FindObject(const char *name) const
Find an object in this list using its name.
Definition: TList.cxx:501
Double_t x[n]
Definition: legend1.C:17
virtual void SetGrid(Int_t valuex=1, Int_t valuey=1)
Definition: TPad.h:323
const Double_t sigma
A doubly linked list.
Definition: TList.h:43
virtual const char * GetParName(Int_t ipar) const
Definition: TF1.h:483
virtual void SetBottomMargin(Float_t bottommargin)
Set Pad bottom margin in fraction of the pad height.
Definition: TAttPad.cxx:100
R__EXTERN TSystem * gSystem
Definition: TSystem.h:539
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2851
tomato 2-D histogram with a float per channel (see TH1 documentation)}
Definition: TH2.h:249
void SetOptFit(Int_t fit=1)
The type of information about fit parameters printed in the histogram statistics box can be selected ...
Definition: TStyle.cxx:1219
const Double_t dxbin
A simple line.
Definition: TLine.h:23
TLine * l
Definition: textangle.C:4
virtual void SetMarkerStyle(Style_t mstyle=1)
Set the marker style.
Definition: TAttMarker.h:40
tomato 1-D histogram with a double per channel (see TH1 documentation)}
Definition: TH1.h:594
const Bool_t kFALSE
Definition: RtypesCore.h:92
virtual Int_t RedirectOutput(const char *name, const char *mode="a", RedirectHandle_t *h=0)
Redirect standard output (stdout, stderr) to the specified file.
Definition: TSystem.cxx:1684
The Canvas class.
Definition: TCanvas.h:31
Double_t Exp(Double_t x)
Definition: TMath.h:622
return c2
Definition: legend2.C:14
double Double_t
Definition: RtypesCore.h:55
virtual Double_t GetEntries() const
Return the current number of entries.
Definition: TH1.cxx:4054
Double_t fdm5(Double_t *xx, Double_t *par)
Double_t GetUymax() const
Definition: TPad.h:225
double f2(const double *x)
1-Dim function class
Definition: TF1.h:150
void SetOptStat(Int_t stat=1)
The type of information printed in the histogram statistics box can be selected via the parameter mod...
Definition: TStyle.cxx:1267
virtual Double_t * GetParameters() const
Definition: TF1.h:474
#define gDirectory
Definition: TDirectory.h:211
virtual void Update()
Update canvas pad buffers.
Definition: TCanvas.cxx:2208
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition: TNamed.cxx:155
TList * GetListOfFunctions() const
Definition: TH1.h:224
THist< 2, float, THistStatContent, THistStatUncertainty > TH2F
Definition: THist.hxx:317
A List of entry numbers in a TTree or TChain.
Definition: TEntryList.h:25
virtual TFitResultPtr Fit(const char *formula, Option_t *option="", Option_t *goption="", Double_t xmin=0, Double_t xmax=0)
Fit histogram with function fname.
Definition: TH1.cxx:3564
void Modified(Bool_t flag=1)
Definition: TPad.h:410
TAxis * GetXaxis()
Definition: TH1.h:300
void SetOptStat(Int_t stat=1)
Set the stat option.
Definition: TPaveStats.cxx:303