Study of the CMS Monte Carlo and Run open data for the dimuon system  1.0
A C++ software for the measurement of the Z resonance, the angle of the negative muon in the Collins–Soper frame of the dimuon system and the forward-backward asymmetry of Drell–Yan events in pp collisions at 8 TeV
graphicalUtilities.h
Go to the documentation of this file.
1 
7 #ifndef GRAPHICALUTILITIES_H
8 #define GRAPHICALUTILITIES_H
9 
10 #include <filesystem>
11 #include <iostream>
12 #include <string>
13 #include <TH2D.h>
14 #include <TPad.h>
15 #include <TVirtualPad.h>
16 #include <TH1D.h>
17 #include <ROOT/RDataFrame.hxx>
18 #include <ROOT/RResultPtr.hxx>
19 #include <ROOT/RVec.hxx>
20 #include <Math/Vector4Dfwd.h>
21 #include <Math/Vector4D.h>
22 #include <TCanvas.h>
23 #include <TLatex.h>
24 #include <TStyle.h>
25 #include <TROOT.h>
26 #include <TLegend.h>
27 
28 void saveHistogram(TCanvas *c, string namehist, string type){
29 
41  namespace fs = std::filesystem;
42 
43  // define path where to save plot
44  type = "/" + type;
45  string subpath = "images" + type;
46  string filepath_pdf = subpath + "/" + namehist + ".pdf";
47  string filepath_png = subpath + "/" + namehist + ".png";
48 
49  //if the folder "images" doesn't exist, it will create it
50  try {
51  if (!fs::is_directory("images") || !fs::exists("images"))
52  throw("images");
53  } catch (const char *path) {
54  std::cerr << "The folder " << path << " does not exist." << std::endl;
55  std::cerr << "Creating folder..." << std::endl;
56  fs::create_directory(path);
57  std::cout << "The folder " << path << " successfully created!" << std::endl;
58  }
59 
60  //if the subfolder "type" doesn't exist, it will create it
61  try {
62  c->SaveAs(filepath_pdf.c_str());
63  c->SaveAs(filepath_png.c_str());
64  if (!fs::is_directory(subpath.c_str()) || !fs::exists(subpath.c_str()))
65  throw(subpath.c_str());
66  }catch(const char *newpath){
67  std::cerr << "The folder " << newpath << " does not exist." << std::endl;
68  std::cerr << "Creating folder..." << std::endl;
69  fs::create_directory(newpath);
70  std::cout << "The folder " << newpath << " successfully created" << std::endl;
71  c->SaveAs(filepath_pdf.c_str());
72  c->SaveAs(filepath_png.c_str());
73  }
74 }
75 
76 
77 void cosHisto(ROOT::RDF::RInterface<ROOT::Detail::RDF::RJittedFilter, void> df_MC, ROOT::RDF::RInterface<ROOT::Detail::RDF::RJittedFilter, void> df_datas, string filename, string rapiditylim, float x1, float y1, float x2, float y2, string canvasname){
78 
95  constexpr int nbins = 40; //widthbins = 0.05
96 
97  //creating new canvas
98  auto c = new TCanvas(canvasname.c_str(), "", 1000, 800);
99 
100  //creating pads
101  auto pad1 = new TPad("pad1","pad1",0., 0.33, 1., 1.00);
102  auto pad2 = new TPad("pad2","pad2",0., 0.08, 1., 0.30);
103  pad1->SetBottomMargin(0.05);
104  pad1->SetBorderMode(0);
105  pad1->SetLogy();
106  pad2->SetTopMargin(0.00001);
107  pad2->SetBottomMargin(0.25);
108  pad2->SetBorderMode(0);
109  pad1->Draw();
110  pad2->Draw();
111 
112  //fill pad1
113  pad1->cd();
114 
115  //creating a histogram
116  auto hist_MC = df_MC.Histo1D({"hist_MC", "", nbins, -1, 1}, "costheta");
117  auto hist_datas = df_datas.Histo1D({"hist_datas", "", nbins, -1, 1}, "costheta");
118  hist_datas->Scale(1./float((df_datas.Count()).GetValue()));
119  hist_MC->Scale(1./float((df_MC.Count()).GetValue()));
120 
121  //creating report
122  auto report_MC = df_MC.Report();
123  report_MC->Print();
124  auto report_datas = df_datas.Report();
125  report_datas->Print();
126 
127  //setting histogram properties
128  hist_MC->GetYaxis()->SetTitle("N_{events,norm}");
129  hist_MC->GetYaxis()->SetTitleSize(0.06);
130  hist_MC->GetYaxis()->SetTitleOffset(0.7);
131  hist_MC->GetYaxis()->SetLabelSize(0.05);
132  hist_MC->SetStats(0);
133  hist_MC->SetFillColor(kOrange-3);
134  hist_MC->SetLineColor(kOrange-3);
135  hist_MC->SetLineStyle(0);
136  hist_datas->SetMarkerStyle(8);
137  hist_datas->SetLineColor(1);
138  hist_datas->SetStats(0);
139 
140  //draw the histogram
141  hist_MC->DrawClone("HF3SAME");
142  hist_datas->DrawClone("PESAME");
143 
144  //writing the rapidity limits on the histogram
145  TLatex label;
146  label.SetTextSize(0.04);
147  label.SetNDC(true);
148  label.DrawLatex(0.15, 0.96, "CMS");
149  label.DrawLatex(0.15, 0.925, rapiditylim.c_str());
150  label.DrawLatex(0.6, 0.94, "#sqrt{s} = 8 TeV, L_{int} = 18.8 fb^{-1}");
151 
152  //writing legend
153  auto legend = new TLegend(x1, y1, x2, y2);
154  legend->AddEntry("hist_MC","MC:Z->#mu#mu","f");
155  legend->AddEntry("hist_datas","Datas","lep");
156  legend->SetBorderSize(0);
157  legend->SetFillColor(0);
158  legend->Draw();
159 
160  //fill pad2
161  pad2->cd();
162 
163  //creating divide histo Data/MC
164  hist_datas->Divide(hist_MC.GetPtr());
165 
166  //properties
167  pad2->SetGridy();
168  hist_datas->GetXaxis()->SetTitle("cos(#theta*)");
169  hist_datas->GetYaxis()->SetTitle("Data/MC");
170  hist_datas->GetXaxis()->SetLabelSize(0.12);
171  hist_datas->GetYaxis()->SetLabelSize(0.08);
172  hist_datas->GetYaxis()->SetLabelOffset(0.012);
173  hist_datas->GetYaxis()->SetTitleSize(0.13);
174  hist_datas->GetYaxis()->SetTitleOffset(0.27);
175  hist_datas->GetXaxis()->SetTitleSize(0.12);
176  hist_datas->GetXaxis()->SetTitleOffset(0.93);
177  hist_datas->SetLineColor(1);
178  hist_datas->SetStats(0);
179  hist_datas->DrawClone("PESAME");
180 
181  //save histogram in ../images/costheta/
182  saveHistogram(c, filename, "costheta");
183 
184 }
185 
186 
187 void dmMassHisto(ROOT::RDF::RInterface<ROOT::Detail::RDF::RJittedFilter, void> df_MC, ROOT::RDF::RInterface<ROOT::Detail::RDF::RJittedFilter, void> df_datas, string filename, string rapiditylim, string canvasname){
188 
201  constexpr int nbins = 100;
202 
203  //creating new canvas
204  auto c = new TCanvas(canvasname.c_str(), "", 1000, 800);
205 
206  //creating pads
207  auto pad1 = new TPad("pad1","pad1",0., 0.33, 1., 1.00);
208  auto pad2 = new TPad("pad2","pad2",0., 0.08, 1., 0.30);
209  pad1->SetBottomMargin(0.05);
210  pad1->SetBorderMode(0);
211  pad1->SetLogy();
212  pad2->SetTopMargin(0.00001);
213  pad2->SetBottomMargin(0.25);
214  pad2->SetBorderMode(0);
215  pad1->Draw();
216  pad2->Draw();
217 
218  //fill pad1
219  pad1->cd();
220 
221  //creating a histogram
222  auto hist_MC = df_MC.Histo1D({"hist_MC", "", nbins, 70, 110}, "dimuon_mass");
223  auto hist_datas = df_datas.Histo1D({"hist_datas", "", nbins, 70, 110}, "dimuon_mass");
224  hist_datas->Scale(1./float((df_datas.Count()).GetValue()));
225  hist_MC->Scale(1./float((df_MC.Count()).GetValue()));
226 
227  //creating report
228  auto report_MC = df_MC.Report();
229  report_MC->Print();
230  auto report_datas = df_datas.Report();
231  report_datas->Print();
232 
233  //setting histogram properties
234  hist_MC->GetYaxis()->SetTitle("N_{Events,norm}");
235  hist_MC->GetYaxis()->SetTitleSize(0.06);
236  hist_MC->GetYaxis()->SetTitleOffset(0.7);
237  hist_MC->GetYaxis()->SetLabelSize(0.05);
238  hist_MC->SetStats(0);
239  hist_MC->SetFillColor(kOrange-3);
240  hist_MC->SetLineColor(kOrange-3);
241  hist_MC->SetLineStyle(0);
242  hist_datas->SetMarkerStyle(8);
243  hist_datas->SetStats(0);
244  hist_datas->SetLineColor(1);
245 
246  //draw the histogram
247  hist_MC->DrawClone("HF3SAME");
248  hist_datas->DrawClone("PESAME");
249 
250  //writing the rapidity limits on the histogram
251  TLatex label;
252  label.SetTextSize(0.04);
253  label.SetTextAlign(12);
254  label.SetNDC(true);
255  label.SetTextAlign(11);
256  label.DrawLatex(0.45, 0.8, "Z");
257  label.DrawLatex(0.15, 0.84, "CMS");
258  label.DrawLatex(0.15, 0.8, rapiditylim.c_str());
259  label.SetTextAlign(31);
260  label.DrawLatex(0.90, 0.92, "#sqrt{s} = 8 TeV, L_{int} = 18.8 fb^{-1}");
261 
262  //writing legend
263  auto legend = new TLegend(0.6,0.80,0.85,0.86);
264  legend->AddEntry("hist_MC","MC:Z->#mu#mu","f");
265  legend->AddEntry("hist_datas","Datas","lep");
266  legend->SetBorderSize(0);
267  legend->SetFillColor(0);
268  legend->Draw();
269 
270  //fill pad2
271  pad2->cd();
272 
273  //creating divide histo Data/MC
274  hist_datas->Divide(hist_MC.GetPtr());
275 
276  //setting histogram properties
277  pad2->SetGridy();
278  hist_datas->SetLineColor(1);
279  hist_datas->GetXaxis()->SetTitle("m_{#mu#mu} (GeV)");
280  hist_datas->GetYaxis()->SetTitle("Data/MC");
281  hist_datas->GetXaxis()->SetLabelSize(0.12);
282  hist_datas->GetYaxis()->SetLabelSize(0.08);
283  hist_datas->GetYaxis()->SetLabelOffset(0.012);
284  hist_datas->GetYaxis()->SetTitleSize(0.13);
285  hist_datas->GetYaxis()->SetTitleOffset(0.27);
286  hist_datas->GetXaxis()->SetTitleSize(0.12);
287  hist_datas->GetXaxis()->SetTitleOffset(0.93);
288  hist_datas->SetStats(0);
289  hist_datas->DrawClone("PESAME");
290 
291  //save histogram in ../images/dimuonspectrum/
292  saveHistogram(c, filename, "dimuonspectrumZ");
293 
294 }
295 
296 
297 auto operationHist( ROOT::RDF::RResultPtr<::TH2D> & histNf, ROOT::RDF::RResultPtr<::TH2D> & histDf, ROOT::RDF::RResultPtr<::TH2D> & histNb, ROOT::RDF::RResultPtr<::TH2D> & histDb){
298 
310  //making operation on histogram
311  histNf->Add(histNb.GetPtr(), -1.0);
312  histDf->Add(histDb.GetPtr(), +1.0);
313  histNf->Divide(histDf.GetPtr());
314  histNf->Scale(0.375);
315 
316  return histNf;
317 }
318 
319 
320 void afbHist(ROOT::RDF::RInterface<ROOT::Detail::RDF::RJittedFilter, void> df_MC, ROOT::RDF::RInterface<ROOT::Detail::RDF::RJittedFilter, void> df_datas, int filenumber, string rapiditylim){
321 
333  //creating two filtered datframes, one with costheta=>0, one with costheta<0
334  auto df_cm_MC = df_MC.Filter("costheta<0", "backward");
335  auto df_cp_MC = df_MC.Filter("costheta>0", "forward");
336  auto df_cm_datas = df_datas.Filter("costheta<0", "backward");
337  auto df_cp_datas = df_datas.Filter("costheta>0", "forward");
338 
339  //creating report
340  auto report_cm_MC = df_cm_MC.Report();
341  report_cm_MC->Print();
342  auto report_cm_datas = df_cm_datas.Report();
343  report_cm_datas->Print();
344  auto report_cp_MC = df_cp_MC.Report();
345  report_cp_MC->Print();
346  auto report_cp_datas = df_cp_datas.Report();
347  report_cp_datas->Print();
348 
349  //creating four histogram 2D of the dimuonmass and rapidity weighted with wn
350  auto histDf_MC = df_cp_MC.Histo2D({"cp,Df", "", 5, 60,120,5,-2.4,2.4}, "dimuon_mass", "y", "wd");
351  auto histNf_MC = df_cp_MC.Histo2D({"cp,Nf", "", 5, 60,120,5,-2.4,2.4}, "dimuon_mass", "y", "wn");
352  auto histDb_MC = df_cm_MC.Histo2D({"cm,Db", "", 5, 60,120,5,-2.4,2.4}, "dimuon_mass", "y", "wd");
353  auto histNb_MC = df_cm_MC.Histo2D({"cm,Nb", "", 5, 60,120,5,-2.4,2.4}, "dimuon_mass", "y", "wn");
354 
355  auto histDf_datas = df_cp_datas.Histo2D({"cp,Df", "", 5, 60,120,5,-2.4,2.4}, "dimuon_mass", "y", "wd");
356  auto histNf_datas = df_cp_datas.Histo2D({"cp,Nf", "", 5, 60,120,5,-2.4,2.4}, "dimuon_mass", "y", "wn");
357  auto histDb_datas = df_cm_datas.Histo2D({"cm,Db", "", 5, 60,120,5,-2.4,2.4}, "dimuon_mass", "y", "wd");
358  auto histNb_datas = df_cm_datas.Histo2D({"cm,Nb", "", 5, 60,120,5,-2.4,2.4}, "dimuon_mass", "y", "wn");
359 
360  //normalizing histograms
361  histDf_MC->Scale(1./float((df_cp_MC.Count()).GetValue()));
362  histNf_MC->Scale(1./float((df_cp_MC.Count()).GetValue()));
363  histDb_MC->Scale(1./float((df_cm_MC.Count()).GetValue()));
364  histNb_MC->Scale(1./float((df_cm_MC.Count()).GetValue()));
365 
366  histDf_datas->Scale(1./float((df_cp_datas.Count()).GetValue()));
367  histNf_datas->Scale(1./float((df_cp_datas.Count()).GetValue()));
368  histDb_datas->Scale(1./float((df_cm_datas.Count()).GetValue()));
369  histNb_datas->Scale(1./float((df_cm_datas.Count()).GetValue()));
370 
371  //operation on histogram
372  auto hist_MC = operationHist(histNf_MC, histDf_MC, histNb_MC, histDb_MC);
373  auto hist_datas = operationHist(histNf_datas, histDf_datas, histNb_datas, histDb_datas);
374 
375  //drawing two subpads
376  auto pad1 = new TPad("pad1","pad1", 0.06, 0.33, 1., 1.00);
377  auto pad2 = new TPad("pad2","pad2", 0.06, 0.08, 1., 0.30);
378  pad1->SetBorderMode(0);
379  pad2->SetTopMargin(0.00001);
380  pad2->SetBorderMode(0);
381  pad2->SetLeftMargin(0.17);
382  pad1->SetLeftMargin(0.17);
383  pad2->SetBottomMargin(0.2);
384  pad1->Draw();
385  pad2->Draw();
386 
387  //fill pad1 with the projections
388  pad1->cd();
389 
390  //projection of the final histogram
391  auto h_MC = hist_MC->ProjectionX("MC");
392  auto h_datas = hist_datas->ProjectionX("datas");
393 
394  //graphical set
395  h_MC->SetStats(0);
396  h_datas->SetStats(0);
397 
398  h_MC->SetLineColor(2);
399  h_MC->GetXaxis()->SetTitle("m_{#mu#mu} (GeV)");
400  h_MC->GetYaxis()->SetTitle("Afb");
401 
402  h_MC->GetXaxis()->CenterTitle(true);
403  h_MC->GetYaxis()->CenterTitle(true);
404 
405  h_MC->GetXaxis()->SetTitleSize(0.13);
406  h_MC->GetXaxis()->SetNdivisions(505);
407  h_MC->GetXaxis()->SetLabelFont(42);
408  h_MC->GetXaxis()->SetLabelSize(0.11);
409  h_MC->GetXaxis()->SetLabelOffset(-0.055);
410  h_MC->GetXaxis()->SetTitleOffset(0.34);
411  h_MC->GetXaxis()->SetTickLength(0.05);
412 
413  h_MC->GetYaxis()->SetRangeUser(-1.3, 1.3);
414  h_MC->GetYaxis()->SetLabelSize(0.09);
415  h_MC->GetYaxis()->SetTitleSize(0.13);
416  h_MC->GetYaxis()->SetLabelOffset(0.01);
417  h_MC->GetYaxis()->SetTitleOffset(0.6);
418  h_MC->GetYaxis()->SetTickLength(0.05);
419 
420  h_datas->SetMarkerColor(1);
421  h_datas->SetLineColor(1);
422  h_datas->SetMarkerStyle(8);
423 
424  //drawing the final histogram
425  h_MC->DrawClone("HIST SAME");
426  h_datas->DrawClone("P SAME");
427 
428  //legend only on the first pad
429  if (filenumber == 1){
430 
431  auto legend = new TLegend(0.30, 0.30, 0.50, 0.45);
432  legend->AddEntry(h_MC, "MC", "l");
433  legend->AddEntry(h_datas, "Datas", "lep");
434  legend->SetBorderSize(0);
435  legend->SetFillColor(0);
436  legend->SetTextSize(0.13);
437  legend->Draw();
438 
439  }
440 
441  //label on the histogram
442  TLatex label;
443  label.SetTextSize(0.11);
444  label.SetTextAlign(12);
445  label.SetNDC(true);
446  label.DrawLatex(0.3, 0.2, rapiditylim.c_str());
447 
448  //residual pad
449  pad2->cd();
450  pad2->SetGridy();
451 
452  //creating histo Data-MC
453  h_datas->Add(h_MC, -1.0);
454 
455  //setting histogram properties
456  h_datas->SetStats(0);
457 
458  h_datas->GetXaxis()->SetTitle("m_{#mu#mu} (GeV)");
459  h_datas->GetXaxis()->SetLabelSize(0.10);
460  h_datas->GetXaxis()->SetNdivisions(505);
461  h_datas->GetXaxis()->SetTitleSize(0.11);
462  h_datas->GetXaxis()->SetTitleOffset(0.86);
463 
464  h_datas->GetYaxis()->SetTitle("Data-MC");
465  h_datas->GetYaxis()->SetLabelSize(0.08);
466  h_datas->GetYaxis()->SetNdivisions(505);
467  h_datas->GetYaxis()->SetTitleOffset(0.6);
468  h_datas->GetYaxis()->SetTitleSize(0.12);
469 
470  h_datas->DrawClone("PESAME");
471 }
472 #endif /* GRAPHICALUTILITIES_H */
auto operationHist(ROOT::RDF::RResultPtr<::TH2D > &histNf, ROOT::RDF::RResultPtr<::TH2D > &histDf, ROOT::RDF::RResultPtr<::TH2D > &histNb, ROOT::RDF::RResultPtr<::TH2D > &histDb)
Definition: graphicalUtilities.h:297
void afbHist(ROOT::RDF::RInterface< ROOT::Detail::RDF::RJittedFilter, void > df_MC, ROOT::RDF::RInterface< ROOT::Detail::RDF::RJittedFilter, void > df_datas, int filenumber, string rapiditylim)
Definition: graphicalUtilities.h:320
void dmMassHisto(ROOT::RDF::RInterface< ROOT::Detail::RDF::RJittedFilter, void > df_MC, ROOT::RDF::RInterface< ROOT::Detail::RDF::RJittedFilter, void > df_datas, string filename, string rapiditylim, string canvasname)
Definition: graphicalUtilities.h:187
void cosHisto(ROOT::RDF::RInterface< ROOT::Detail::RDF::RJittedFilter, void > df_MC, ROOT::RDF::RInterface< ROOT::Detail::RDF::RJittedFilter, void > df_datas, string filename, string rapiditylim, float x1, float y1, float x2, float y2, string canvasname)
Definition: graphicalUtilities.h:77
void saveHistogram(TCanvas *c, string namehist, string type)
Definition: graphicalUtilities.h:28
type
Definition: openfiles.py:68