CbmRoot
CbmAnaDielectronTaskDrawAll.cxx
Go to the documentation of this file.
1 
8 
9 #include "CbmDrawHist.h"
10 #include "CbmHistManager.h"
11 #include "CbmUtils.h"
12 
13 #include <iomanip>
14 #include <iostream>
15 #include <string>
16 
17 #include "TCanvas.h"
18 #include "TClass.h"
19 #include "TEllipse.h"
20 #include "TF1.h"
21 #include "TFile.h"
22 #include "TH1.h"
23 #include "TH1D.h"
24 #include "TH2D.h"
25 #include "TKey.h"
26 #include "TLatex.h"
27 #include "TMath.h"
28 #include "TStyle.h"
29 #include "TSystem.h"
30 #include "TText.h"
31 #include <TLegend.h>
32 
33 using namespace std;
34 using namespace Cbm;
35 
37  const string& fileNameInmed,
38  const string& fileNameQgp,
39  const string& fileNameOmega,
40  const string& fileNamePhi,
41  const string& fileNameOmegaDalitz,
42  const string& outputDir,
43  Bool_t useMvd) {
44  fOutputDir = outputDir;
45  fUseMvd = useMvd;
46  fDrawQgp = (fileNameQgp != "");
47 
48  //SetDefaultDrawStyle();
49  vector<string> fileNames = {fileNameInmed,
50  fileNameQgp,
51  fileNameOmega,
52  fileNamePhi,
53  fileNameOmegaDalitz};
54 
55  fHM.resize(fNofSignals);
56  for (int i = 0; i < fNofSignals; i++) {
57  fHM[i] = new CbmHistManager();
58  if (!fDrawQgp && i == kQgp) continue;
59  TFile* file = new TFile(fileNames[i].c_str());
60  fHM[i]->ReadFromFile(file);
61  Int_t nofEvents = (int) H1(i, "fh_event_number")->GetEntries();
62  fHM[i]->ScaleByPattern(".*", 1. / nofEvents);
63  cout << "nofEvents = " << nofEvents << endl;
64  }
65 
66  // index: AnalysisSteps
67  fh_mean_bg_minv.resize(CbmLmvmHist::fNofAnaSteps);
68  fh_mean_eta_minv.resize(CbmLmvmHist::fNofAnaSteps);
69  fh_mean_pi0_minv.resize(CbmLmvmHist::fNofAnaSteps);
70  fh_sum_s_minv.resize(CbmLmvmHist::fNofAnaSteps);
71  fh_mean_eta_minv_pt.resize(CbmLmvmHist::fNofAnaSteps);
72  fh_mean_pi0_minv_pt.resize(CbmLmvmHist::fNofAnaSteps);
73  fh_mean_sbg_vs_minv.resize(CbmLmvmHist::fNofAnaSteps);
74 
75  FillMeanHist();
76  FillSumSignalsHist();
77  CalcCutEffRange(0.0, 0.2);
78  CalcCutEffRange(0.2, 0.6);
79  CalcCutEffRange(0.6, 1.2);
80  SBgRangeAll();
81  DrawSBgSignals();
82  DrawMinvAll();
83  DrawMinvPtAll();
84  DrawSBgVsMinv();
85  SaveHist();
86  SaveCanvasToImage();
87 }
88 
89 
90 TH1D* CbmAnaDielectronTaskDrawAll::H1(int signalType, const string& name) {
91  return (TH1D*) fHM[signalType]->H1(name);
92 }
93 
94 TH2D* CbmAnaDielectronTaskDrawAll::H2(int signalType, const string& name) {
95  return (TH2D*) fHM[signalType]->H1(name);
96 }
97 
99  TH1D* sInmed =
100  (TH1D*) H1(kInmed, "fh_signal_minv_" + CbmLmvmHist::fAnaSteps[step])
101  ->Clone();
102  TH1D* sQgp =
103  (fDrawQgp)
104  ? (TH1D*) H1(kQgp, "fh_signal_minv_" + CbmLmvmHist::fAnaSteps[step])
105  ->Clone()
106  : nullptr;
107  TH1D* sOmega =
108  (TH1D*) H1(kOmega, "fh_signal_minv_" + CbmLmvmHist::fAnaSteps[step])
109  ->Clone();
110  TH1D* sPhi =
111  (TH1D*) H1(kPhi, "fh_signal_minv_" + CbmLmvmHist::fAnaSteps[step])->Clone();
112  TH1D* sEta = fh_mean_eta_minv[step];
113  TH1D* sPi0 = fh_mean_pi0_minv[step];
114  TH1D* sOmegaDalitz =
115  (TH1D*) H1(kOmegaD, "fh_signal_minv_" + CbmLmvmHist::fAnaSteps[step])
116  ->Clone();
117 
118  TH1D* coctail = (TH1D*) sInmed->Clone();
119  if (fDrawQgp) coctail->Add(sQgp);
120  coctail->Add(sOmega);
121  coctail->Add(sPhi);
122  coctail->Add(sEta);
123  coctail->Add(sPi0);
124  coctail->Add(sOmegaDalitz);
125 
126  return coctail;
127 }
128 
130  fHM[0]->CreateCanvas("minv_all_mc", "minv_all_mc", 800, 800);
131  DrawMinv(kMc);
132 
133  fHM[0]->CreateCanvas("minv_all_acc", "minv_all_acc", 800, 800);
134  DrawMinv(kAcc);
135 
136  fHM[0]->CreateCanvas("minv_all_ptcut", "minv_all_ptcut", 800, 800);
137  DrawMinv(kPtCut);
138 
139  fHM[0]->CreateCanvas("minv_all_ttcut", "minv_all_ttcut", 800, 800);
140  DrawMinv(kTtCut);
141 
142  fHM[0]->CreateCanvas("minv_all_elid", "minv_all_elid", 800, 800);
143  DrawMinv(kElId);
144 }
145 
147  TH1D* sInmed =
148  (TH1D*) H1(kInmed, "fh_signal_minv_" + CbmLmvmHist::fAnaSteps[step])
149  ->Clone();
150  TH1D* sQgp =
151  (fDrawQgp)
152  ? (TH1D*) H1(kQgp, "fh_signal_minv_" + CbmLmvmHist::fAnaSteps[step])
153  ->Clone()
154  : nullptr;
155  TH1D* sOmega =
156  (TH1D*) H1(kOmega, "fh_signal_minv_" + CbmLmvmHist::fAnaSteps[step])
157  ->Clone();
158  TH1D* sPhi =
159  (TH1D*) H1(kPhi, "fh_signal_minv_" + CbmLmvmHist::fAnaSteps[step])->Clone();
160  TH1D* bg = (TH1D*) fh_mean_bg_minv[step]->Clone();
161  TH1D* sEta = (TH1D*) fh_mean_eta_minv[step]->Clone();
162  TH1D* sPi0 = (TH1D*) fh_mean_pi0_minv[step]->Clone();
163  TH1D* sOmegaDalitz =
164  (TH1D*) H1(kOmegaD, "fh_signal_minv_" + CbmLmvmHist::fAnaSteps[step])
165  ->Clone();
166 
167  TH1D* coctail = GetCoctailMinv(step);
168 
169  TH1D* sbg = (TH1D*) bg->Clone();
170  sbg->Add(sInmed);
171  if (fDrawQgp) sbg->Add(sQgp);
172  sbg->Add(sOmega);
173  sbg->Add(sPhi);
174  sbg->Add(sEta);
175  sbg->Add(sPi0);
176  sbg->Add(sOmegaDalitz);
177 
178 
179  int nRebin = 20;
180  sbg->Rebin(nRebin);
181  coctail->Rebin(nRebin);
182  bg->Rebin(nRebin);
183  sPi0->Rebin(nRebin);
184  sEta->Rebin(nRebin);
185  sOmegaDalitz->Rebin(nRebin);
186  sOmega->Rebin(nRebin);
187  sInmed->Rebin(nRebin);
188  if (fDrawQgp) sQgp->Rebin(nRebin);
189  sPhi->Rebin(nRebin);
190 
191  /*sbg->Scale(1./nRebin);
192  coctail->Scale(1./nRebin);
193  bg->Scale(1./nRebin);
194  sPi0->Scale(1./nRebin);
195  sEta->Scale(1./nRebin);
196  sOmegaDalitz->Scale(1./nRebin);
197  sOmega->Scale(1./nRebin);
198  sInmed->Scale(1./nRebin);
199  sQgp->Scale(1./nRebin);
200  sPhi->Scale(1./nRebin);*/
201 
202  double binWidth = sbg->GetBinWidth(1);
203  sbg->Scale(1. / binWidth);
204  coctail->Scale(1. / binWidth);
205  bg->Scale(1. / binWidth);
206  sPi0->Scale(1. / binWidth);
207  sEta->Scale(1. / binWidth);
208  sOmegaDalitz->Scale(1. / binWidth);
209  sOmega->Scale(1. / binWidth);
210  sInmed->Scale(1. / binWidth);
211  if (fDrawQgp) sQgp->Scale(1. / binWidth);
212  sPhi->Scale(1. / binWidth);
213 
214 
215  sbg->SetMinimum(5e-8);
216  sbg->SetMaximum(2e-2);
217  sbg->GetXaxis()->SetRangeUser(0, 2.);
218  bg->GetXaxis()->SetRangeUser(0, 2.);
219  coctail->GetXaxis()->SetRangeUser(0, 2.);
220  sPi0->GetXaxis()->SetRangeUser(0, 2.);
221  sEta->GetXaxis()->SetRangeUser(0, 2.);
222  sOmegaDalitz->GetXaxis()->SetRangeUser(0, 2.);
223  sOmega->GetXaxis()->SetRangeUser(0, 2.);
224  sInmed->GetXaxis()->SetRangeUser(0, 2.);
225  if (fDrawQgp) sQgp->GetXaxis()->SetRangeUser(0, 2.);
226  sPhi->GetXaxis()->SetRangeUser(0, 2.);
227 
228  /*
229  if (step == kMc) {
230  DrawH1({coctail, sPi0, sEta, sOmegaDalitz, sOmega, sInmed, sQgp, sPhi},
231  {"", "", "", "", "", "", "", ""}, kLinear, kLog, false, 0.8, 0.8, 0.99, 0.99);
232  } else {
233  DrawH1({sbg, bg, coctail, sPi0, sEta, sOmegaDalitz, sOmega, sInmed, sQgp, sPhi},
234  {"", "", "", "", "", "", "", "", "", ""}, kLinear, kLog, false, 0.8, 0.8, 0.99, 0.99);
235  }
236  */
237 
238  if (step == kMc) {
239  if (fDrawQgp) {
240  DrawH1({coctail, sPi0, sEta, sOmegaDalitz, sOmega, sInmed, sQgp, sPhi},
241  {"", "", "", "", "", "", "", ""},
242  kLinear,
243  kLog,
244  false,
245  0.8,
246  0.8,
247  0.99,
248  0.99,
249  "HIST L");
250  } else {
251  DrawH1({coctail, sPi0, sEta, sOmegaDalitz, sOmega, sInmed, sPhi},
252  {"", "", "", "", "", "", ""},
253  kLinear,
254  kLog,
255  false,
256  0.8,
257  0.8,
258  0.99,
259  0.99,
260  "HIST L");
261  }
262  } else {
263  if (fDrawQgp) {
264  DrawH1({sbg,
265  bg,
266  coctail,
267  sPi0,
268  sEta,
269  sOmegaDalitz,
270  sOmega,
271  sInmed,
272  sQgp,
273  sPhi},
274  {"", "", "", "", "", "", "", "", "", ""},
275  kLinear,
276  kLog,
277  false,
278  0.8,
279  0.8,
280  0.99,
281  0.99,
282  "HIST L");
283  } else {
284  DrawH1({sbg, bg, coctail, sPi0, sEta, sOmegaDalitz, sOmega, sInmed, sPhi},
285  {"", "", "", "", "", "", "", "", ""},
286  kLinear,
287  kLog,
288  false,
289  0.8,
290  0.8,
291  0.99,
292  0.99,
293  "HIST L");
294  }
295  }
296 
297  string yTitle = "dN/dM_{ee} [GeV/c^{2}]^{-1}";
298  coctail->GetYaxis()->SetTitle(yTitle.c_str());
299  sbg->GetYaxis()->SetTitle(yTitle.c_str());
300  coctail->GetYaxis()->SetLabelSize(0.05);
301  sbg->GetYaxis()->SetLabelSize(0.05);
302 
303  sInmed->SetFillColor(kMagenta - 3);
304  sInmed->SetLineColor(kMagenta - 2);
305  sInmed->SetLineStyle(0);
306  sInmed->SetLineWidth(3);
307  sInmed->SetFillStyle(3344);
308 
309  if (fDrawQgp) {
310  sQgp->SetFillColor(kOrange - 2);
311  sQgp->SetLineColor(kOrange - 3);
312  sQgp->SetLineStyle(0);
313  sQgp->SetLineWidth(3);
314  sQgp->SetFillStyle(3444);
315  }
316 
317 
318  sOmega->SetFillColor(kOrange + 7);
319  sOmega->SetLineColor(kOrange + 4);
320  sOmega->SetLineStyle(0);
321  sOmega->SetLineWidth(2);
322 
323  sPhi->SetFillColor(kAzure + 2);
324  sPhi->SetLineColor(kAzure + 3);
325  sPhi->SetLineStyle(0);
326  sPhi->SetLineWidth(2);
327  sPhi->SetFillStyle(3112);
328  gStyle->SetHatchesLineWidth(1);
329  gStyle->SetHatchesSpacing(1.);
330 
331 
332  bg->SetFillColor(kGray);
333  bg->SetLineColor(kBlack);
334  bg->SetLineStyle(0);
335  bg->SetLineWidth(1);
336 
337 
338  sEta->SetFillColor(kRed - 4);
339  sEta->SetLineColor(kRed + 2);
340  sEta->SetLineStyle(0);
341  sEta->SetLineWidth(2);
342 
343 
344  sPi0->SetFillColor(kGreen - 3);
345  sPi0->SetLineColor(kGreen + 3);
346  sPi0->SetLineStyle(0);
347  sPi0->SetLineWidth(2);
348 
349  sOmegaDalitz->SetFillColor(kCyan + 2);
350  sOmegaDalitz->SetLineColor(kCyan + 4);
351  sOmegaDalitz->SetLineStyle(0);
352  sOmegaDalitz->SetLineWidth(2);
353 
354  sbg->SetFillColor(kBlack);
355  sbg->SetLineColor(kBlack);
356  sbg->SetLineStyle(0);
357  sbg->SetLineWidth(1);
358 
359  coctail->SetLineColor(kRed + 2);
360  coctail->SetFillStyle(0);
361  coctail->SetLineWidth(3);
362 
363  if (step != kMc) {
364  TLegend* legend = new TLegend(0.7, 0.6, 0.99, 0.99);
365  legend->SetFillColor(kWhite);
366  legend->SetTextSize(0.04);
367  legend->AddEntry(sOmega, "#omega #rightarrow e^{+}e^{-}", "f");
368  legend->AddEntry(sOmegaDalitz, "#omega #rightarrow #pi^{0}e^{+}e^{-}", "f");
369  legend->AddEntry(sPhi, "#phi #rightarrow e^{+}e^{-}", "f");
370  legend->AddEntry(sInmed, "in-medium #rho", "f");
371  if (fDrawQgp) legend->AddEntry(sQgp, "QGP radiation", "f");
372  legend->AddEntry(sEta, "#eta #rightarrow #gammae^{+}e^{-}", "f");
373  legend->AddEntry(sPi0, "#pi^{0} #rightarrow #gammae^{+}e^{-}", "f");
374  legend->AddEntry(coctail, "Coctail", "f");
375  legend->AddEntry(bg, "Background", "f");
376  legend->AddEntry(sbg, "Coctail+BG", "f");
377  legend->Draw();
378  } else {
379  TLegend* legend = new TLegend(0.7, 0.7, 0.99, 0.99);
380  legend->SetFillColor(kWhite);
381  legend->SetTextSize(0.04);
382  legend->AddEntry(sOmega, "#omega #rightarrow e^{+}e^{-}", "f");
383  legend->AddEntry(sOmegaDalitz, "#omega #rightarrow #pi^{0}e^{+}e^{-}", "f");
384  legend->AddEntry(sPhi, "#phi #rightarrow e^{+}e^{-}", "f");
385  legend->AddEntry(sInmed, "in-medium #rho", "f");
386  if (fDrawQgp) legend->AddEntry(sQgp, "QGP radiation", "f");
387  legend->AddEntry(sEta, "#eta #rightarrow #gammae^{+}e^{-}", "f");
388  legend->AddEntry(sPi0, "#pi^{0} #rightarrow #gammae^{+}e^{-}", "f");
389  legend->AddEntry(coctail, "Coctail", "f");
390  legend->Draw();
391  }
392 
393 
394  gPad->SetLogy(true);
395 }
396 
398  TH1D* bg = (TH1D*) fh_mean_bg_minv[kTtCut]->Clone();
399  TH1D* coctail = GetCoctailMinv(kTtCut);
400  fh_mean_sbg_vs_minv[kTtCut] =
401  new TH1D(("fh_sbg_vs_minv_" + CbmLmvmHist::fAnaSteps[kTtCut]).c_str(),
402  ("fh_sbg_vs_minv_" + CbmLmvmHist::fAnaSteps[kTtCut]
403  + ";M_{ee} [GeV/c^{2}];Cocktail/Background")
404  .c_str(),
405  bg->GetNbinsX(),
406  bg->GetXaxis()->GetXmin(),
407  bg->GetXaxis()->GetXmax());
408  fh_mean_sbg_vs_minv[kTtCut]->Divide(coctail, bg, 1., 1., "B");
409  fh_mean_sbg_vs_minv[kTtCut]->Rebin(20);
410  fh_mean_sbg_vs_minv[kTtCut]->Scale(1. / 20.);
411  fh_mean_sbg_vs_minv[kTtCut]->GetXaxis()->SetRangeUser(0, 2.);
412 
413  bg = (TH1D*) fh_mean_bg_minv[kPtCut]->Clone();
414  coctail = GetCoctailMinv(kPtCut);
415  fh_mean_sbg_vs_minv[kPtCut] =
416  new TH1D(("fh_sbg_vs_minv_" + CbmLmvmHist::fAnaSteps[kPtCut]).c_str(),
417  ("fh_sbg_vs_minv_" + CbmLmvmHist::fAnaSteps[kPtCut]
418  + ";M_{ee} [GeV/c^{2}];Cocktail/Background")
419  .c_str(),
420  bg->GetNbinsX(),
421  bg->GetXaxis()->GetXmin(),
422  bg->GetXaxis()->GetXmax());
423  fh_mean_sbg_vs_minv[kPtCut]->Divide(coctail, bg, 1., 1., "B");
424  fh_mean_sbg_vs_minv[kPtCut]->Rebin(20);
425  fh_mean_sbg_vs_minv[kPtCut]->Scale(1. / 20.);
426  fh_mean_sbg_vs_minv[kPtCut]->GetXaxis()->SetRangeUser(0, 2.);
427 
428  fHM[0]->CreateCanvas("lmvm_sbg_vs_minv", "lmvm_sbg_vs_minv", 800, 800);
429  DrawH1({fh_mean_sbg_vs_minv[kTtCut], fh_mean_sbg_vs_minv[kPtCut]},
430  {"Without Pt cut", "With Pt cut"},
431  kLinear,
432  kLog,
433  true,
434  0.6,
435  0.85,
436  0.99,
437  0.99);
438  gPad->SetLogy(true);
439 }
440 
442  fHM[0]->CreateCanvas("minv_pt_ptcut", "minv_pt_ptcut", 800, 800);
443  DrawMinvPt(kPtCut);
444 }
445 
447  TH2D* sInmed =
448  (TH2D*) H2(kInmed, "fh_signal_minv_pt_" + CbmLmvmHist::fAnaSteps[step])
449  ->Clone();
450  TH2D* sQgp =
451  (fDrawQgp)
452  ? (TH2D*) H2(kQgp, "fh_signal_minv_pt_" + CbmLmvmHist::fAnaSteps[step])
453  ->Clone()
454  : nullptr;
455  TH2D* sOmega =
456  (TH2D*) H2(kOmega, "fh_signal_minv_pt_" + CbmLmvmHist::fAnaSteps[step])
457  ->Clone();
458  TH2D* sOmegaDalitz =
459  (TH2D*) H2(kOmegaD, "fh_signal_minv_pt_" + CbmLmvmHist::fAnaSteps[step])
460  ->Clone();
461  TH2D* sPhi =
462  (TH2D*) H2(kPhi, "fh_signal_minv_pt_" + CbmLmvmHist::fAnaSteps[step])
463  ->Clone();
464  TH2D* sEta = fh_mean_eta_minv_pt[step];
465  TH2D* sPi0 = fh_mean_pi0_minv_pt[step];
466 
467  TH2D* coctail = (TH2D*) sInmed->Clone();
468  if (fDrawQgp) coctail->Add(sQgp);
469  coctail->Add(sOmega);
470  coctail->Add(sPhi);
471  coctail->Add(sOmegaDalitz);
472  coctail->Add(sEta);
473  coctail->Add(sPi0);
474 
475  DrawH2(coctail);
476 }
477 
479  for (int step = 0; step < CbmLmvmHist::fNofAnaSteps; step++) {
480  for (int iS = 0; iS < fNofSignals; iS++) {
481  if (!fDrawQgp && iS == kQgp) continue;
482  if (iS == 0) {
483  fh_mean_bg_minv[step] =
484  (TH1D*) H1(iS, "fh_bg_minv_" + CbmLmvmHist::fAnaSteps[step])->Clone();
485  fh_mean_eta_minv[step] =
486  (TH1D*) H1(iS, "fh_eta_minv_" + CbmLmvmHist::fAnaSteps[step])
487  ->Clone();
488  fh_mean_pi0_minv[step] =
489  (TH1D*) H1(iS, "fh_pi0_minv_" + CbmLmvmHist::fAnaSteps[step])
490  ->Clone();
491  fh_mean_eta_minv_pt[step] =
492  (TH2D*) H2(iS, "fh_eta_minv_pt_" + CbmLmvmHist::fAnaSteps[step])
493  ->Clone();
494  fh_mean_pi0_minv_pt[step] =
495  (TH2D*) H2(iS, "fh_pi0_minv_pt_" + CbmLmvmHist::fAnaSteps[step])
496  ->Clone();
497  } else {
498  fh_mean_bg_minv[step]->Add(
499  (TH1D*) H1(iS, "fh_bg_minv_" + CbmLmvmHist::fAnaSteps[step])
500  ->Clone());
501  fh_mean_eta_minv[step]->Add(
502  (TH1D*) H1(iS, "fh_eta_minv_" + CbmLmvmHist::fAnaSteps[step])
503  ->Clone());
504  fh_mean_pi0_minv[step]->Add(
505  (TH1D*) H1(iS, "fh_pi0_minv_" + CbmLmvmHist::fAnaSteps[step])
506  ->Clone());
507  fh_mean_eta_minv_pt[step]->Add(
508  (TH2D*) H2(iS, "fh_eta_minv_pt_" + CbmLmvmHist::fAnaSteps[step])
509  ->Clone());
510  fh_mean_pi0_minv_pt[step]->Add(
511  (TH2D*) H2(iS, "fh_pi0_minv_pt_" + CbmLmvmHist::fAnaSteps[step])
512  ->Clone());
513  }
514  }
515  //TODO: nofSignals without QGP
516  fh_mean_bg_minv[step]->Scale(1. / (double) fNofSignals);
517  fh_mean_eta_minv[step]->Scale(1. / (double) fNofSignals);
518  fh_mean_pi0_minv[step]->Scale(1. / (double) fNofSignals);
519  fh_mean_eta_minv_pt[step]->Scale(1. / (double) fNofSignals);
520  fh_mean_pi0_minv_pt[step]->Scale(1. / (double) fNofSignals);
521  }
522 }
523 
525  if (fOutputDir != "") {
526  gSystem->mkdir(fOutputDir.c_str(), true);
527  TFile* f = TFile::Open(string(fOutputDir + "/draw_all_hist.root").c_str(),
528  "RECREATE");
529  for (int i = 0; i < CbmLmvmHist::fNofAnaSteps; i++) {
530  fh_mean_bg_minv[i]->Write();
531  fh_mean_eta_minv[i]->Write();
532  fh_mean_pi0_minv[i]->Write();
533  }
534  fh_mean_sbg_vs_minv[kTtCut]->Write();
535  fh_mean_sbg_vs_minv[kPtCut]->Write();
536  f->Close();
537  }
538 }
539 
541  for (int step = 0; step < CbmLmvmHist::fNofAnaSteps; step++) {
542  fh_sum_s_minv[step] =
543  (TH1D*) H1(kInmed, "fh_signal_minv_" + CbmLmvmHist::fAnaSteps[step])
544  ->Clone();
545  if (fDrawQgp)
546  fh_sum_s_minv[step]->Add(
547  (TH1D*) H1(kQgp, "fh_signal_minv_" + CbmLmvmHist::fAnaSteps[step])
548  ->Clone());
549  fh_sum_s_minv[step]->Add(
550  (TH1D*) H1(kOmega, "fh_signal_minv_" + CbmLmvmHist::fAnaSteps[step])
551  ->Clone());
552  fh_sum_s_minv[step]->Add(
553  (TH1D*) H1(kPhi, "fh_signal_minv_" + CbmLmvmHist::fAnaSteps[step])
554  ->Clone());
555  fh_sum_s_minv[step]->Add(
556  (TH1D*) H1(kOmegaD, "fh_signal_minv_" + CbmLmvmHist::fAnaSteps[step])
557  ->Clone());
558  fh_sum_s_minv[step]->Add((TH1D*) fh_mean_eta_minv[step]->Clone());
559  fh_sum_s_minv[step]->Add((TH1D*) fh_mean_pi0_minv[step]->Clone());
560  }
561 }
562 
564  Double_t maxMinv) {
565  stringstream ss1;
566  ss1 << minMinv << "_" << maxMinv;
567  TH1D* grS = new TH1D(("grS_" + ss1.str()).c_str(),
568  ";Analysis step;Efficiency [%]",
570  0,
572  TH1D* grB = new TH1D(("grB_" + ss1.str()).c_str(),
573  ";Analysis step;Efficiency [%]",
575  0,
577  int x = 1;
578  for (int step = kElId; step < CbmLmvmHist::fNofAnaSteps; step++) {
579  if (!fUseMvd && (step == kMvd1Cut || step == kMvd2Cut)) { continue; }
580  Int_t x1 = fh_sum_s_minv[step]->FindBin(minMinv);
581  Int_t x2 = fh_sum_s_minv[step]->FindBin(maxMinv);
582 
583  double yS = 100. * fh_sum_s_minv[step]->Integral(x1, x2)
584  / fh_sum_s_minv[kElId]->Integral(x1, x2);
585  double yB = 100. * fh_mean_bg_minv[step]->Integral(x1, x2)
586  / fh_mean_bg_minv[kElId]->Integral(x1, x2);
587 
588  grB->GetXaxis()->SetBinLabel(x, CbmLmvmHist::fAnaStepsLatex[step].c_str());
589  grB->SetBinContent(x, yB);
590  grS->SetBinContent(x, yS);
591  x++;
592  }
593 
594  grB->GetXaxis()->SetLabelSize(0.06);
595  grB->GetXaxis()->SetRange(1, x - 1);
596  grS->GetXaxis()->SetRange(1, x - 1);
597 
598  stringstream ss;
599  ss << "lmvm_cut_eff_" << minMinv << "_" << maxMinv;
600  fHM[0]->CreateCanvas(ss.str().c_str(), ss.str().c_str(), 700, 700);
601  DrawH1(
602  {grB, grS}, {"BG", "Signal"}, kLinear, kLinear, true, 0.75, 0.85, 1.0, 1.0);
603  grS->SetLineWidth(4);
604  grB->SetLineWidth(4);
605  grB->SetMinimum(1);
606  grB->SetMaximum(105);
607 
608  stringstream ss2;
609  ss2 << minMinv << "<M [GeV/c^2]<" << maxMinv;
610  TText* t = new TText(0.5, 110, ss2.str().c_str());
611  t->Draw();
612 }
613 
614 
615 TH1D* CbmAnaDielectronTaskDrawAll::SBgRange(Double_t min, Double_t max) {
616  stringstream ss;
617  ss << "lmvm_s_bg_region_" << min << "_" << max;
618  TH1D* h_s_bg = new TH1D(ss.str().c_str(),
619  string(ss.str() + ";Analysis steps;S/BG").c_str(),
621  0,
623  h_s_bg->GetXaxis()->SetLabelSize(0.06);
624  int x = 1;
625  for (int step = kElId; step < CbmLmvmHist::fNofAnaSteps; step++) {
626  if (!fUseMvd && (step == kMvd1Cut || step == kMvd2Cut)) { continue; }
627  Int_t bin1 = fh_sum_s_minv[step]->FindBin(min);
628  Int_t bin2 = fh_sum_s_minv[step]->FindBin(max);
629  double y = fh_sum_s_minv[step]->Integral(bin1, bin2)
630  / fh_mean_bg_minv[step]->Integral(bin1, bin2);
631 
632  h_s_bg->GetXaxis()->SetBinLabel(x,
633  CbmLmvmHist::fAnaStepsLatex[step].c_str());
634  h_s_bg->SetBinContent(x, y);
635  // replace "." with "_"
636  string str = ss.str();
637  for (string::iterator it = str.begin(); it < str.end(); it++) {
638  if (*it == '.') *it = '_';
639  }
640  x++;
641  }
642  h_s_bg->GetXaxis()->SetRange(1, x - 1);
643  return h_s_bg;
644 }
645 
647  TH1D* h_00_02 = SBgRange(0.0, 0.2);
648  TH1D* h_02_06 = SBgRange(0.2, 0.6);
649  TH1D* h_06_12 = SBgRange(0.6, 1.2);
650 
651  fHM[0]->CreateCanvas("lmvm_s_bg_ranges", "lmvm_s_bg_ranges", 700, 700);
652  DrawH1(
653  {h_00_02, h_02_06, h_06_12},
654  {"0.0<M [GeV/c^{2}]<0.2", "0.2<M [GeV/c^{2}]<0.6", "0.6<M [GeV/c^{2}]<1.2"},
655  kLinear,
656  kLog,
657  true,
658  0.25,
659  0.8,
660  0.75,
661  0.99);
662 
663  h_00_02->SetMinimum(1e-3);
664  h_00_02->SetMaximum(3);
665  h_00_02->SetLineWidth(4);
666  h_02_06->SetLineWidth(4);
667  h_06_12->SetLineWidth(4);
668 
669  TH1D* h_05_06 = SBgRange(0.5, 0.6);
670  fHM[0]->CreateCanvas(
671  "lmvm_s_bg_ranges_05_06", "lmvm_s_bg_ranges_05_06", 700, 700);
672  DrawH1(h_05_06, kLinear, kLinear);
673  h_05_06->SetMinimum(1e-3);
674  h_05_06->SetMaximum(2e-2);
675  h_05_06->SetLineWidth(4);
676 }
677 
679  // Double_t y[CbmLmvmHist::fNofAnaSteps];
680  TCanvas* cFit =
681  fHM[0]->CreateCanvas("lmvm_signal_fit", "lmvm_signal_fit", 600, 600);
682  TCanvas* cDashboard =
683  fHM[0]->CreateCanvas("lmvm_dashboard", "lmvm_dashboard", 1000, 900);
684  int iDash = 2;
685  TLatex* latex = new TLatex();
686  latex->SetTextSize(0.03);
687  latex->DrawLatex(0.05, 0.95, "signal");
688  latex->DrawLatex(0.2, 0.95, "step");
689  latex->DrawLatex(0.4, 0.95, "eff, %");
690  latex->DrawLatex(0.55, 0.95, "S/BG");
691  latex->DrawLatex(0.7, 0.95, "mean");
692  latex->DrawLatex(0.85, 0.95, "sigma");
693  TString str;
694  for (int iF = 0; iF < fNofSignals - 1; iF++) {
695  if (!fDrawQgp && iF == kQgp) continue;
696  string signalName = CbmLmvmHist::fSignalNames[iF];
697  cout << "Signal: " << signalName << endl;
698  stringstream ss;
699  ss << "lmvm_s_bg_cuts_" << signalName;
700 
701  TH1D* h_s_bg = new TH1D(ss.str().c_str(),
702  string(ss.str() + ";Analysis steps;S/BG").c_str(),
704  0,
706  h_s_bg->GetXaxis()->SetLabelSize(0.06);
707  h_s_bg->SetLineWidth(4);
708  int x = 1;
709  iDash++; // empty string after each signal
710  for (int step = 0; step < CbmLmvmHist::fNofAnaSteps; step++) {
711  if (step < kElId) continue;
712  if (!fUseMvd && (step == kMvd1Cut || step == kMvd2Cut)) { continue; }
713 
714  TH1D* s = (TH1D*) H1(iF, "fh_signal_minv_" + CbmLmvmHist::fAnaSteps[step])
715  ->Clone();
716  TH1D* bg = (TH1D*) fh_mean_bg_minv[step]->Clone();
717  cFit->cd();
718  if (iF == kPhi) {
719  if (s->GetEntries() > 0) s->Fit("gaus", "Q", "", 0.95, 1.05);
720  } else if (iF == kOmega) {
721  if (s->GetEntries() > 0) s->Fit("gaus", "Q", "", 0.69, 0.81);
722  } else {
723  if (s->GetEntries() > 0) s->Fit("gaus", "Q");
724  }
725 
726  TF1* func = s->GetFunction("gaus");
727  Double_t mean = (func != NULL) ? func->GetParameter("Mean") : 0.;
728  Double_t sigma = (func != NULL) ? func->GetParameter("Sigma") : 0.;
729  Int_t minInd = s->FindBin(mean - 2. * sigma);
730  Int_t maxInd = s->FindBin(mean + 2. * sigma);
731 
732  Double_t sumSignal = 0.;
733  Double_t sumBg = 0.;
734  for (Int_t i = minInd + 1; i <= maxInd - 1; i++) {
735  sumSignal += s->GetBinContent(i);
736  sumBg += bg->GetBinContent(i);
737  }
738  Double_t sbg = sumSignal / sumBg;
739  double eff =
740  100.
741  * H1(iF, "fh_signal_pty_" + CbmLmvmHist::fAnaSteps[step])->GetEntries()
742  / H1(iF, "fh_signal_pty_" + CbmLmvmHist::fAnaSteps[kMc])->GetEntries();
743 
744  bool isOmegaOrPhi = (iF == kPhi || iF == kOmega);
745  cDashboard->cd();
746  latex->DrawLatex(0.05, 1.0 - iDash * 0.033, signalName.c_str());
747  latex->DrawLatex(
748  0.2, 1.0 - iDash * 0.033, CbmLmvmHist::fAnaSteps[step].c_str());
749  str.Form("%.2f", eff);
750  latex->DrawLatex(0.4, 1.0 - iDash * 0.033, str.Data());
751  str.Form("%.3f", sumSignal / sumBg);
752  latex->DrawLatex(
753  0.55, 1.0 - iDash * 0.033, (isOmegaOrPhi) ? str.Data() : "-");
754  str.Form("%.1f", 1000. * mean);
755  latex->DrawLatex(
756  0.7, 1.0 - iDash * 0.033, (isOmegaOrPhi) ? str.Data() : "-");
757  str.Form("%.1f", 1000. * sigma);
758  latex->DrawLatex(
759  0.85, 1.0 - iDash * 0.033, (isOmegaOrPhi) ? str.Data() : "-");
760 
761  h_s_bg->GetXaxis()->SetBinLabel(
762  x, CbmLmvmHist::fAnaStepsLatex[step].c_str());
763  if (sbg < 1000.) h_s_bg->SetBinContent(x, sbg);
764  x++;
765  iDash++;
766  }
767  h_s_bg->GetXaxis()->SetRange(1, x - 1);
768  fHM[0]->CreateCanvas(ss.str().c_str(), ss.str().c_str(), 800, 800);
769  DrawH1(h_s_bg);
770  h_s_bg->SetLineWidth(4);
771 
772  cDashboard->Draw();
773  }
774 }
775 
777  cout << "Images output dir:" << fOutputDir << endl;
778  fHM[0]->SaveCanvasToImage(fOutputDir, "eps;png");
779 }
780 
CbmAnaDielectronTaskDrawAll.h
kMc
@ kMc
Definition: CbmLmvmHist.h:16
CbmAnaDielectronTaskDrawAll::H2
TH2D * H2(int signalType, const std::string &name)
Definition: CbmAnaDielectronTaskDrawAll.cxx:94
CbmAnaDielectronTaskDrawAll::DrawMinv
void DrawMinv(CbmLmvmAnalysisSteps step)
Draw invariant mass spectra for all signal types for specified analysis step.
Definition: CbmAnaDielectronTaskDrawAll.cxx:146
f
float f
Definition: L1/vectors/P4_F32vec4.h:24
kElId
@ kElId
Definition: CbmLmvmHist.h:20
CbmAnaDielectronTaskDrawAll::DrawMinvPtAll
void DrawMinvPtAll()
Draw invariant mass vs Pt histograms.
Definition: CbmAnaDielectronTaskDrawAll.cxx:441
kTtCut
@ kTtCut
Definition: CbmLmvmHist.h:26
CbmLmvmHist::fNofAnaSteps
static const int fNofAnaSteps
Definition: CbmLmvmHist.h:49
CbmAnaDielectronTaskDrawAll::SaveCanvasToImage
void SaveCanvasToImage()
Save all created canvases to images.
Definition: CbmAnaDielectronTaskDrawAll.cxx:776
CbmAnaDielectronTaskDrawAll::DrawMinvAll
void DrawMinvAll()
Draw invariant mass histograms.
Definition: CbmAnaDielectronTaskDrawAll.cxx:129
kPtCut
@ kPtCut
Definition: CbmLmvmHist.h:27
i
int i
Definition: L1/vectors/P4_F32vec4.h:25
CbmAnaDielectronTaskDrawAll::DrawSBgSignals
void DrawSBgSignals()
Draw S/BG vs plots for different signals.
Definition: CbmAnaDielectronTaskDrawAll.cxx:678
kOmega
@ kOmega
Definition: CbmAnaDielectronTaskDrawAll.h:25
CbmAnaDielectronTaskDrawAll::SaveHist
void SaveHist()
Save histograms for the study report.
Definition: CbmAnaDielectronTaskDrawAll.cxx:524
CbmLmvmHist::fSignalNames
static const std::vector< std::string > fSignalNames
Definition: CbmLmvmHist.h:54
kOmegaD
@ kOmegaD
Definition: CbmAnaDielectronTaskDrawAll.h:25
CbmAnaDielectronTaskDrawAll::FillSumSignalsHist
void FillSumSignalsHist()
Fill sum signals.
Definition: CbmAnaDielectronTaskDrawAll.cxx:540
CbmDrawHist.h
Helper functions for drawing 1D and 2D histograms and graphs.
CbmHistManager.h
Histogram manager.
CbmAnaDielectronTaskDrawAll::DrawMinvPt
void DrawMinvPt(CbmLmvmAnalysisSteps step)
Draw invariant mass spectra vs Pt for all signal types for specified analysis step.
Definition: CbmAnaDielectronTaskDrawAll.cxx:446
CbmAnaDielectronTaskDrawAll::SBgRangeAll
void SBgRangeAll()
Draw S/BG vs plots for different mass ranges.
Definition: CbmAnaDielectronTaskDrawAll.cxx:646
min
friend F32vec4 min(const F32vec4 &a, const F32vec4 &b)
Definition: L1/vectors/P4_F32vec4.h:33
CbmLmvmHist::fAnaSteps
static const std::vector< std::string > fAnaSteps
Definition: CbmLmvmHist.h:50
kQgp
@ kQgp
Definition: CbmAnaDielectronTaskDrawAll.h:25
DrawH1
void DrawH1(TH1 *hist, HistScale logx, HistScale logy, const string &drawOpt, Int_t color, Int_t lineWidth, Int_t lineStyle, Int_t markerSize, Int_t markerStyle)
Definition: CbmDrawHist.cxx:49
CbmLmvmHist::fAnaStepsLatex
static const std::vector< std::string > fAnaStepsLatex
Definition: CbmLmvmHist.h:51
CbmHistManager
Histogram manager.
Definition: CbmHistManager.h:41
kLinear
@ kLinear
Definition: CbmDrawHist.h:79
ClassImp
ClassImp(CbmAnaDielectronTaskDrawAll)
CbmAnaDielectronTaskDrawAll::DrawHistosFromFile
void DrawHistosFromFile(const std::string &fileNameInmed, const std::string &fileNameQgp, const std::string &fileNameOmega, const std::string &fileNamePhi, const std::string &fileNameOmegaDalitz, const std::string &outputDir="", Bool_t useMvd=false)
Implement functionality of drawing histograms in the macro from the specified files,...
Definition: CbmAnaDielectronTaskDrawAll.cxx:36
kInmed
@ kInmed
Definition: CbmAnaDielectronTaskDrawAll.h:25
CbmUtils.h
CbmAnaDielectronTaskDrawAll::FillMeanHist
void FillMeanHist()
It creates a mean histogram from 4 files.
Definition: CbmAnaDielectronTaskDrawAll.cxx:478
CbmAnaDielectronTaskDrawAll::CalcCutEffRange
void CalcCutEffRange(Double_t minMinv, Double_t maxMinv)
Calculate cut efficiency in specified invariant mass region.
Definition: CbmAnaDielectronTaskDrawAll.cxx:563
kMvd1Cut
@ kMvd1Cut
Definition: CbmLmvmHist.h:22
CbmAnaDielectronTaskDrawAll::SBgRange
TH1D * SBgRange(Double_t minMinv, Double_t maxMinv)
Create S/BG vs cuts for specified invariant mass range.
Definition: CbmAnaDielectronTaskDrawAll.cxx:615
CbmAnaDielectronTaskDrawAll::H1
TH1D * H1(int signalType, const std::string &name)
Definition: CbmAnaDielectronTaskDrawAll.cxx:90
kMvd2Cut
@ kMvd2Cut
Definition: CbmLmvmHist.h:23
x
Double_t x
Definition: CbmMvdSensorDigiToHitTask.cxx:68
y
Double_t y
Definition: CbmMvdSensorDigiToHitTask.cxx:68
DrawH2
void DrawH2(TH2 *hist, HistScale logx, HistScale logy, HistScale logz, const string &drawOpt)
Definition: CbmDrawHist.cxx:84
CbmAnaDielectronTaskDrawAll
Definition: CbmAnaDielectronTaskDrawAll.h:27
max
friend F32vec4 max(const F32vec4 &a, const F32vec4 &b)
Definition: L1/vectors/P4_F32vec4.h:36
CbmAnaDielectronTaskDrawAll::GetCoctailMinv
TH1D * GetCoctailMinv(CbmLmvmAnalysisSteps step)
Create and return cotail vs. minv.
Definition: CbmAnaDielectronTaskDrawAll.cxx:98
CbmAnaDielectronTaskDrawAll::DrawSBgVsMinv
void DrawSBgVsMinv()
Draw S/Bg vs minv.
Definition: CbmAnaDielectronTaskDrawAll.cxx:397
kAcc
@ kAcc
Definition: CbmLmvmHist.h:17
kLog
@ kLog
Definition: CbmDrawHist.h:78
kPhi
@ kPhi
Definition: CbmAnaDielectronTaskDrawAll.h:25
CbmLmvmAnalysisSteps
CbmLmvmAnalysisSteps
Definition: CbmLmvmHist.h:15
Cbm
Definition: CbmGeometryUtils.cxx:31