Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
hist103_NanoAODHiggsAnalysis.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_histv7
3/// \notebook -draw
4/// An example of complex analysis with RDataFrame and RHist: reconstructing the Higgs boson.
5///
6/// This tutorial is based on df103_NanoAODHiggsAnalysis.C.
7/// It is a simplified but yet complex example of an analysis reconstructing
8/// the Higgs boson decaying to two Z bosons from events with four leptons. The data
9/// and simulated events are taken from CERN OpenData representing a subset of the data
10/// recorded in 2012 with the CMS detector at the LHC. The tutorials follows the Higgs
11/// to four leptons analysis published on CERN Open Data portal
12/// ([10.7483/OPENDATA.CMS.JKB8.RR42](http://opendata.cern.ch/record/5500)).
13/// The resulting plots show the invariant mass of the selected four lepton systems
14/// in different decay modes (four muons, four electrons and two of each kind)
15/// and in a combined plot indicating the decay of the Higgs boson with a mass
16/// of about 125 GeV.
17///
18/// The following steps are performed for each sample with data and simulated events
19/// in order to reconstruct the Higgs boson from the selected muons and electrons:
20/// 1. Select interesting events with multiple cuts on event properties, e.g.,
21/// number of leptons, kinematics of the leptons and quality of the tracks.
22/// 2. Reconstruct two Z bosons of which only one on the mass shell from the selected events and apply additional cuts
23/// on the reconstructed objects.
24/// 3. Reconstruct the Higgs boson from the remaining Z boson candidates and calculate
25/// its invariant mass.
26///
27/// The tutorial has the fast mode enabled by default, which reads the data from already skimmed
28/// datasets with a total size of only 51MB. If the fast mode is disabled, the tutorial runs over
29/// the full dataset with a size of 12GB.
30///
31/// \macro_image
32/// \macro_code
33/// \macro_output
34///
35/// \date October 2018, February 2026
36/// \author Stefan Wunsch (KIT, CERN), the ROOT Team
37
38#include <ROOT/RDataFrame.hxx>
39#include <ROOT/RHist.hxx>
40#include <ROOT/RVec.hxx>
42#include <TCanvas.h>
43#include <TH1D.h>
44#include <TLatex.h>
45#include <TLegend.h>
48#include <TStyle.h>
49#include <string>
50
51using namespace ROOT::VecOps;
53using cRVecF = const ROOT::RVecF &;
54
55const auto z_mass = 91.2;
56
57// Select interesting events with four muons
58RNode selection_4mu(RNode df)
59{
60 auto df_ge4m = df.Filter("nMuon>=4", "At least four muons");
61 auto df_iso = df_ge4m.Filter("All(abs(Muon_pfRelIso04_all)<0.40)", "Require good isolation");
62 auto df_kin = df_iso.Filter("All(Muon_pt>5) && All(abs(Muon_eta)<2.4)", "Good muon kinematics");
63 auto df_ip3d = df_kin.Define("Muon_ip3d", "sqrt(Muon_dxy*Muon_dxy + Muon_dz*Muon_dz)");
64 auto df_sip3d = df_ip3d.Define("Muon_sip3d", "Muon_ip3d/sqrt(Muon_dxyErr*Muon_dxyErr + Muon_dzErr*Muon_dzErr)");
65 auto df_pv = df_sip3d.Filter("All(Muon_sip3d<4) && All(abs(Muon_dxy)<0.5) && All(abs(Muon_dz)<1.0)",
66 "Track close to primary vertex with small uncertainty");
67 auto df_2p2n = df_pv.Filter("nMuon==4 && Sum(Muon_charge==1)==2 && Sum(Muon_charge==-1)==2",
68 "Two positive and two negative muons");
69 return df_2p2n;
70}
71
72// Select interesting events with four electrons
73RNode selection_4el(RNode df)
74{
75 auto df_ge4el = df.Filter("nElectron>=4", "At least our electrons");
76 auto df_iso = df_ge4el.Filter("All(abs(Electron_pfRelIso03_all)<0.40)", "Require good isolation");
77 auto df_kin = df_iso.Filter("All(Electron_pt>7) && All(abs(Electron_eta)<2.5)", "Good Electron kinematics");
78 auto df_ip3d = df_kin.Define("Electron_ip3d", "sqrt(Electron_dxy*Electron_dxy + Electron_dz*Electron_dz)");
79 auto df_sip3d = df_ip3d.Define(
80 "Electron_sip3d", "Electron_ip3d/sqrt(Electron_dxyErr*Electron_dxyErr + Electron_dzErr*Electron_dzErr)");
81 auto df_pv = df_sip3d.Filter("All(Electron_sip3d<4) && All(abs(Electron_dxy)<0.5) && "
82 "All(abs(Electron_dz)<1.0)",
83 "Track close to primary vertex with small uncertainty");
84 auto df_2p2n = df_pv.Filter("nElectron==4 && Sum(Electron_charge==1)==2 && Sum(Electron_charge==-1)==2",
85 "Two positive and two negative electrons");
86 return df_2p2n;
87}
88
89// Select interesting events with two electrons and two muons
90RNode selection_2el2mu(RNode df)
91{
92 auto df_ge2el2mu = df.Filter("nElectron>=2 && nMuon>=2", "At least two electrons and two muons");
93 auto df_eta = df_ge2el2mu.Filter("All(abs(Electron_eta)<2.5) && All(abs(Muon_eta)<2.4)", "Eta cuts");
94 auto pt_cuts = [](cRVecF mu_pt, cRVecF el_pt) {
96 if (mu_pt_sorted[0] > 20 && mu_pt_sorted[1] > 10) {
97 return true;
98 }
100 if (el_pt_sorted[0] > 20 && el_pt_sorted[1] > 10) {
101 return true;
102 }
103 return false;
104 };
105 auto df_pt = df_eta.Filter(pt_cuts, {"Muon_pt", "Electron_pt"}, "Pt cuts");
107 auto mu_dr = DeltaR(mu_eta[0], mu_eta[1], mu_phi[0], mu_phi[1]);
108 auto el_dr = DeltaR(el_eta[0], el_eta[1], el_phi[0], el_phi[1]);
109 if (mu_dr < 0.02 || el_dr < 0.02) {
110 return false;
111 }
112 return true;
113 };
114 auto df_dr = df_pt.Filter(dr_cuts, {"Muon_eta", "Muon_phi", "Electron_eta", "Electron_phi"}, "Dr cuts");
115 auto df_iso = df_dr.Filter("All(abs(Electron_pfRelIso03_all)<0.40) && All(abs(Muon_pfRelIso04_all)<0.40)",
116 "Require good isolation");
117 auto df_el_ip3d = df_iso.Define("Electron_ip3d_el", "sqrt(Electron_dxy*Electron_dxy + Electron_dz*Electron_dz)");
118 auto df_el_sip3d = df_el_ip3d.Define("Electron_sip3d_el", "Electron_ip3d_el/sqrt(Electron_dxyErr*Electron_dxyErr + "
119 "Electron_dzErr*Electron_dzErr)");
120 auto df_el_track =
121 df_el_sip3d.Filter("All(Electron_sip3d_el<4) && All(abs(Electron_dxy)<0.5) && All(abs(Electron_dz)<1.0)",
122 "Electron track close to primary vertex with small uncertainty");
123 auto df_mu_ip3d = df_el_track.Define("Muon_ip3d_mu", "sqrt(Muon_dxy*Muon_dxy + Muon_dz*Muon_dz)");
124 auto df_mu_sip3d =
125 df_mu_ip3d.Define("Muon_sip3d_mu", "Muon_ip3d_mu/sqrt(Muon_dxyErr*Muon_dxyErr + Muon_dzErr*Muon_dzErr)");
126 auto df_mu_track = df_mu_sip3d.Filter("All(Muon_sip3d_mu<4) && All(abs(Muon_dxy)<0.5) && All(abs(Muon_dz)<1.0)",
127 "Muon track close to primary vertex with small uncertainty");
128 auto df_2p2n = df_mu_track.Filter("Sum(Electron_charge)==0 && Sum(Muon_charge)==0",
129 "Two opposite charged electron and muon pairs");
130 return df_2p2n;
131}
132
133// Reconstruct two Z candidates from four leptons of the same kind
135{
136 RVec<RVec<size_t>> idx(2);
137 idx[0].reserve(2);
138 idx[1].reserve(2);
139
140 // Find first lepton pair with invariant mass closest to Z mass
141 auto idx_cmb = Combinations(pt, 2);
142 auto best_mass = -1;
143 size_t best_i1 = 0;
144 size_t best_i2 = 0;
145 for (size_t i = 0; i < idx_cmb[0].size(); i++) {
146 const auto i1 = idx_cmb[0][i];
147 const auto i2 = idx_cmb[1][i];
148 if (charge[i1] != charge[i2]) {
149 ROOT::Math::PtEtaPhiMVector p1(pt[i1], eta[i1], phi[i1], mass[i1]);
150 ROOT::Math::PtEtaPhiMVector p2(pt[i2], eta[i2], phi[i2], mass[i2]);
151 const auto this_mass = (p1 + p2).M();
152 if (std::abs(z_mass - this_mass) < std::abs(z_mass - best_mass)) {
154 best_i1 = i1;
155 best_i2 = i2;
156 }
157 }
158 }
159 idx[0].emplace_back(best_i1);
160 idx[0].emplace_back(best_i2);
161
162 // Reconstruct second Z from remaining lepton pair
163 for (size_t i = 0; i < 4; i++) {
164 if (i != best_i1 && i != best_i2) {
165 idx[1].emplace_back(i);
166 }
167 }
168
169 // Return indices of the pairs building two Z bosons
170 return idx;
171}
172
173// Compute Z masses from four leptons of the same kind and sort ascending in distance to Z mass
175{
177 for (size_t i = 0; i < 2; i++) {
178 const auto i1 = idx[i][0];
179 const auto i2 = idx[i][1];
180 ROOT::Math::PtEtaPhiMVector p1(pt[i1], eta[i1], phi[i1], mass[i1]);
181 ROOT::Math::PtEtaPhiMVector p2(pt[i2], eta[i2], phi[i2], mass[i2]);
182 z_masses[i] = (p1 + p2).M();
183 }
184 if (std::abs(z_masses[0] - z_mass) < std::abs(z_masses[1] - z_mass)) {
185 return z_masses;
186 } else {
187 return Reverse(z_masses);
188 }
189}
190
191// Compute mass of Higgs from four leptons of the same kind
192float compute_higgs_mass_4l(const RVec<RVec<size_t>> &idx, cRVecF pt, cRVecF eta, cRVecF phi, cRVecF mass)
193{
194 const auto i1 = idx[0][0];
195 const auto i2 = idx[0][1];
196 const auto i3 = idx[1][0];
197 const auto i4 = idx[1][1];
198 ROOT::Math::PtEtaPhiMVector p1(pt[i1], eta[i1], phi[i1], mass[i1]);
199 ROOT::Math::PtEtaPhiMVector p2(pt[i2], eta[i2], phi[i2], mass[i2]);
200 ROOT::Math::PtEtaPhiMVector p3(pt[i3], eta[i3], phi[i3], mass[i3]);
201 ROOT::Math::PtEtaPhiMVector p4(pt[i4], eta[i4], phi[i4], mass[i4]);
202 return (p1 + p2 + p3 + p4).M();
203}
204
205// Apply selection on reconstructed Z candidates
207{
208 auto df_z1_cut = df.Filter("Z_mass[0] > 40 && Z_mass[0] < 120", "Mass of first Z candidate in [40, 120]");
209 auto df_z2_cut = df_z1_cut.Filter("Z_mass[1] > 12 && Z_mass[1] < 120", "Mass of second Z candidate in [12, 120]");
210 return df_z2_cut;
211}
212
213// Reconstruct Higgs from four muons
214RNode reco_higgs_to_4mu(RNode df)
215{
216 // Filter interesting events
217 auto df_base = selection_4mu(df);
218
219 // Reconstruct Z systems
220 auto df_z_idx =
221 df_base.Define("Z_idx", reco_zz_to_4l, {"Muon_pt", "Muon_eta", "Muon_phi", "Muon_mass", "Muon_charge"});
222
223 // Cut on distance between muons building Z systems
224 auto filter_z_dr = [](const RVec<RVec<size_t>> &idx, cRVecF eta, cRVecF phi) {
225 for (size_t i = 0; i < 2; i++) {
226 const auto i1 = idx[i][0];
227 const auto i2 = idx[i][1];
228 const auto dr = DeltaR(eta[i1], eta[i2], phi[i1], phi[i2]);
229 if (dr < 0.02) {
230 return false;
231 }
232 }
233 return true;
234 };
235 auto df_z_dr =
236 df_z_idx.Filter(filter_z_dr, {"Z_idx", "Muon_eta", "Muon_phi"}, "Delta R separation of muons building Z system");
237
238 // Compute masses of Z systems
239 auto df_z_mass =
240 df_z_dr.Define("Z_mass", compute_z_masses_4l, {"Z_idx", "Muon_pt", "Muon_eta", "Muon_phi", "Muon_mass"});
241
242 // Cut on mass of Z candidates
244
245 // Reconstruct H mass
246 auto df_h_mass =
247 df_z_cut.Define("H_mass", compute_higgs_mass_4l, {"Z_idx", "Muon_pt", "Muon_eta", "Muon_phi", "Muon_mass"});
248
249 return df_h_mass;
250}
251
252// Reconstruct Higgs from four electrons
253RNode reco_higgs_to_4el(RNode df)
254{
255 // Filter interesting events
256 auto df_base = selection_4el(df);
257
258 // Reconstruct Z systems
259 auto df_z_idx = df_base.Define("Z_idx", reco_zz_to_4l,
260 {"Electron_pt", "Electron_eta", "Electron_phi", "Electron_mass", "Electron_charge"});
261
262 // Cut on distance between Electrons building Z systems
263 auto filter_z_dr = [](const RVec<RVec<size_t>> &idx, cRVecF eta, cRVecF phi) {
264 for (size_t i = 0; i < 2; i++) {
265 const auto i1 = idx[i][0];
266 const auto i2 = idx[i][1];
267 const auto dr = DeltaR(eta[i1], eta[i2], phi[i1], phi[i2]);
268 if (dr < 0.02) {
269 return false;
270 }
271 }
272 return true;
273 };
274 auto df_z_dr = df_z_idx.Filter(filter_z_dr, {"Z_idx", "Electron_eta", "Electron_phi"},
275 "Delta R separation of Electrons building Z system");
276
277 // Compute masses of Z systems
278 auto df_z_mass = df_z_dr.Define("Z_mass", compute_z_masses_4l,
279 {"Z_idx", "Electron_pt", "Electron_eta", "Electron_phi", "Electron_mass"});
280
281 // Cut on mass of Z candidates
283
284 // Reconstruct H mass
285 auto df_h_mass = df_z_cut.Define("H_mass", compute_higgs_mass_4l,
286 {"Z_idx", "Electron_pt", "Electron_eta", "Electron_phi", "Electron_mass"});
287
288 return df_h_mass;
289}
290
291// Compute mass of two Z candidates from two electrons and two muons and sort ascending in distance to Z mass
294{
299 auto mu_z = (p1 + p2).M();
300 auto el_z = (p3 + p4).M();
302 if (std::abs(mu_z - z_mass) < std::abs(el_z - z_mass)) {
303 z_masses[0] = mu_z;
304 z_masses[1] = el_z;
305 } else {
306 z_masses[0] = el_z;
307 z_masses[1] = mu_z;
308 }
309 return z_masses;
310}
311
312// Compute Higgs mass from two electrons and two muons
315{
320 return (p1 + p2 + p3 + p4).M();
321}
322
323// Reconstruct Higgs from two electrons and two muons
325{
326 // Filter interesting events
327 auto df_base = selection_2el2mu(df);
328
329 // Compute masses of Z systems
330 auto df_z_mass = df_base.Define(
331 "Z_mass", compute_z_masses_2el2mu,
332 {"Electron_pt", "Electron_eta", "Electron_phi", "Electron_mass", "Muon_pt", "Muon_eta", "Muon_phi", "Muon_mass"});
333
334 // Cut on mass of Z candidates
336
337 // Reconstruct H mass
338 auto df_h_mass = df_z_cut.Define(
340 {"Electron_pt", "Electron_eta", "Electron_phi", "Electron_mass", "Muon_pt", "Muon_eta", "Muon_phi", "Muon_mass"});
341
342 return df_h_mass;
343}
344
345// Plot invariant mass for signal and background processes from simulated events
346// overlay the measured data.
347template <typename T>
348void plot(T sig, T bkg, T data, const std::string &x_label, const std::string &filename)
349{
350 // Canvas and general style options
351 gStyle->SetTextFont(42);
352 auto c = new TCanvas("", "", 800, 700);
353 c->SetLeftMargin(0.15);
354
355 // Get signal and background histograms and stack them to show Higgs signal
356 // on top of the background process
357 auto h_bkg = ROOT::Experimental::Hist::ConvertToTH1D(*bkg).release();
358 auto h_cmb = ROOT::Experimental::Hist::ConvertToTH1D(*sig).release();
359
360 h_cmb->Add(h_bkg);
361 h_cmb->SetTitle("");
362 h_cmb->GetXaxis()->SetTitle(x_label.c_str());
363 h_cmb->GetXaxis()->SetTitleSize(0.04);
364 h_cmb->GetYaxis()->SetTitle("N_{Events}");
365 h_cmb->GetYaxis()->SetTitleSize(0.04);
366 h_cmb->SetLineColor(kRed);
367 h_cmb->SetLineWidth(2);
368 h_cmb->SetMaximum(18);
369 h_cmb->SetStats(kFALSE);
370
371 h_bkg->SetLineWidth(2);
372 h_bkg->SetFillStyle(1001);
373 h_bkg->SetLineColor(kBlack);
374 h_bkg->SetFillColor(kAzure - 9);
375
376 // Get histogram of data points
377 auto h_data = ROOT::Experimental::Hist::ConvertToTH1D(*data).release();
378 h_data->SetLineWidth(1);
379 h_data->SetMarkerStyle(20);
380 h_data->SetMarkerSize(1.0);
381 h_data->SetMarkerColor(kBlack);
382 h_data->SetLineColor(kBlack);
383
384 // Draw histograms
385 h_cmb->Draw("HIST");
386 h_bkg->Draw("HIST SAME");
387 h_data->Draw("PE1 SAME");
388
389 // Add legend
390 auto legend = new TLegend(0.62, 0.70, 0.82, 0.88);
391 legend->SetFillColor(0);
392 legend->SetBorderSize(0);
393 legend->SetTextSize(0.03);
394 legend->AddEntry(h_data, "Data", "pe");
395 legend->AddEntry(h_bkg, "ZZ", "f");
396 legend->AddEntry(h_cmb, "m_{H} = 125 GeV", "f");
397 legend->Draw();
398
399 // Add header
401 cms_label.SetTextSize(0.04);
402 cms_label.DrawLatexNDC(0.16, 0.92, "#bf{CMS Open Data}");
403 TLatex header;
404 header.SetTextSize(0.03);
405 header.DrawLatexNDC(0.63, 0.92, "#sqrt{s} = 8 TeV, L_{int} = 11.6 fb^{-1}");
406
407 // Save plot
408 c->SaveAs(filename.c_str());
409}
410
411void hist103_NanoAODHiggsAnalysis(const bool run_fast = true)
412{
413 // Enable multi-threading
415
416 // In fast mode, take samples from */cms_opendata_2012_nanoaod_skimmed/*, which has
417 // the preselections from the selection_* functions already applied.
418 std::string path = "root://eospublic.cern.ch//eos/root-eos/cms_opendata_2012_nanoaod/";
419 if (run_fast)
420 path = "root://eospublic.cern.ch//eos/root-eos/cms_opendata_2012_nanoaod_skimmed/";
421
422 // Create dataframes for signal, background and data samples
423
424 // Signal: Higgs -> 4 leptons
425 ROOT::RDataFrame df_sig_4l("Events", path + "SMHiggsToZZTo4L.root");
426
427 // Background: ZZ -> 4 leptons
428 // Note that additional background processes from the original paper with minor contribution were left out for this
429 // tutorial.
430 ROOT::RDataFrame df_bkg_4mu("Events", path + "ZZTo4mu.root");
431 ROOT::RDataFrame df_bkg_4el("Events", path + "ZZTo4e.root");
432 ROOT::RDataFrame df_bkg_2el2mu("Events", path + "ZZTo2e2mu.root");
433
434 // CMS data taken in 2012 (11.6 fb^-1 integrated luminosity)
436 {path + "Run2012B_DoubleMuParked.root", path + "Run2012C_DoubleMuParked.root"});
438 {path + "Run2012B_DoubleElectron.root", path + "Run2012C_DoubleElectron.root"});
439
440 // Reconstruct Higgs to 4 muons
442 const auto luminosity = 11580.0; // Integrated luminosity of the data samples
443 const auto xsec_SMHiggsToZZTo4L = 0.0065; // H->4l: Standard Model cross-section
444 const auto nevt_SMHiggsToZZTo4L = 299973.0; // H->4l: Number of simulated events
445 const ROOT::Experimental::RRegularAxis axis(36, {70, 180});
446 auto df_h_sig_4mu =
447 df_sig_4mu_reco.Define("weight", [&] { return luminosity * xsec_SMHiggsToZZTo4L / nevt_SMHiggsToZZTo4L; })
448 .Hist({axis}, {"H_mass"}, "weight");
449
450 const auto scale_ZZTo4l = 1.386; // ZZ->4mu: Scale factor for ZZ to four leptons
451 const auto xsec_ZZTo4mu = 0.077; // ZZ->4mu: Standard Model cross-section
452 const auto nevt_ZZTo4mu = 1499064.0; // ZZ->4mu: Number of simulated events
454 auto df_h_bkg_4mu =
455 df_bkg_4mu_reco.Define("weight", [&] { return luminosity * xsec_ZZTo4mu * scale_ZZTo4l / nevt_ZZTo4mu; })
456 .Hist({axis}, {"H_mass"}, "weight");
457
459 auto df_h_data_4mu = df_data_4mu_reco.Define("weight", [] { return 1.0; }).Hist({axis}, {"H_mass"}, "weight");
460
461 // Reconstruct Higgs to 4 electrons
463 auto df_h_sig_4el =
464 df_sig_4el_reco.Define("weight", [&] { return luminosity * xsec_SMHiggsToZZTo4L / nevt_SMHiggsToZZTo4L; })
465 .Hist({axis}, {"H_mass"}, "weight");
466
467 const auto xsec_ZZTo4el = xsec_ZZTo4mu; // ZZ->4el: Standard Model cross-section
468 const auto nevt_ZZTo4el = 1499093.0; // ZZ->4el: Number of simulated events
470 auto df_h_bkg_4el =
471 df_bkg_4el_reco.Define("weight", [&] { return luminosity * xsec_ZZTo4el * scale_ZZTo4l / nevt_ZZTo4el; })
472 .Hist({axis}, {"H_mass"}, "weight");
473
475 auto df_h_data_4el = df_data_4el_reco.Define("weight", [] { return 1.0; }).Hist({axis}, {"H_mass"}, "weight");
476
477 // Reconstruct Higgs to 2 electrons and 2 muons
479 auto df_h_sig_2el2mu =
480 df_sig_2el2mu_reco.Define("weight", [&] { return luminosity * xsec_SMHiggsToZZTo4L / nevt_SMHiggsToZZTo4L; })
481 .Hist({axis}, {"H_mass"}, "weight");
482
483 const auto xsec_ZZTo2el2mu = 0.18; // ZZ->2el2mu: Standard Model cross-section
484 const auto nevt_ZZTo2el2mu = 1497445.0; // ZZ->2el2mu: Number of simulated events
486 auto df_h_bkg_2el2mu =
487 df_bkg_2el2mu_reco.Define("weight", [&] { return luminosity * xsec_ZZTo2el2mu * scale_ZZTo4l / nevt_ZZTo2el2mu; })
488 .Hist({axis}, {"H_mass"}, "weight");
489
491 auto df_h_data_2el2mu = df_data_2el2mu_reco.Define("weight", [] { return 1.0; }).Hist({axis}, {"H_mass"}, "weight");
492
493 // RunGraphs allows to run the event loops of the separate RDataFrame graphs
494 // concurrently. This results in an improved usage of the available resources
495 // if each separate RDataFrame can not utilize all available resources, e.g.,
496 // because not enough data is available.
499
500 // Make plots
501 plot(df_h_sig_4mu, df_h_bkg_4mu, df_h_data_4mu, "m_{4#mu} (GeV)", "hist103_higgs_4mu.pdf");
502 plot(df_h_sig_4el, df_h_bkg_4el, df_h_data_4el, "m_{4e} (GeV)", "hist103_higgs_4el.pdf");
503 plot(df_h_sig_2el2mu, df_h_bkg_2el2mu, df_h_data_2el2mu, "m_{2e2#mu} (GeV)", "hist103_higgs_2el2mu.pdf");
504
505 // Combine channels for final plot
506 auto h_data_4l = df_h_data_4mu.GetPtr();
509 auto h_sig_4l = df_h_sig_4mu.GetPtr();
510 h_sig_4l->Add(*df_h_sig_4el);
512 auto h_bkg_4l = df_h_bkg_4mu.GetPtr();
513 h_bkg_4l->Add(*df_h_bkg_4el);
515 plot(h_sig_4l, h_bkg_4l, h_data_4l, "m_{4l} (GeV)", "hist103_higgs_4l.pdf");
516}
517
518int main()
519{
520 hist103_NanoAODHiggsAnalysis(/*fast=*/true);
521}
int main()
Definition Prototype.cxx:12
#define c(i)
Definition RSha256.hxx:101
constexpr Bool_t kFALSE
Definition RtypesCore.h:108
@ kRed
Definition Rtypes.h:67
@ kBlack
Definition Rtypes.h:66
@ kAzure
Definition Rtypes.h:68
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
winID h TVirtualViewer3D TVirtualGLPainter char TVirtualGLPainter plot
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char filename
R__EXTERN TStyle * gStyle
Definition TStyle.h:442
A regular axis with equidistant bins in the interval .
ROOT's RDataFrame offers a modern, high-level interface for analysis of data stored in TTree ,...
virtual void SetTextFont(Font_t tfont=62)
Set the text font.
Definition TAttText.h:48
virtual void SetTextSize(Float_t tsize=1)
Set the text size.
Definition TAttText.h:49
The Canvas class.
Definition TCanvas.h:23
To draw Mathematical Formula.
Definition TLatex.h:20
TLatex * DrawLatexNDC(Double_t x, Double_t y, const char *text)
Draw this TLatex with new coordinates in NDC.
Definition TLatex.cxx:1963
This class displays a legend box (TPaveText) containing several legend entries.
Definition TLegend.h:23
virtual void SaveAs(const char *filename="", Option_t *option="") const
Save this object in the file specified by filename.
Definition TObject.cxx:706
TPaveText * pt
Vector1::Scalar DeltaR(const Vector1 &v1, const Vector2 &v2)
Find difference in pseudorapidity (Eta) and Phi between two generic vectors The only requirements on ...
Definition VectorUtil.h:113
RVec< T > Reverse(const RVec< T > &v)
Return copy of reversed vector.
Definition RVec.hxx:2476
RVec< RVec< std::size_t > > Combinations(const std::size_t size1, const std::size_t size2)
Return the indices that represent all combinations of the elements of two RVecs.
Definition RVec.hxx:2601
std::unique_ptr< TH1D > ConvertToTH1D(const RHistEngine< double > &engine)
Convert a one-dimensional histogram to TH1D.
unsigned int RunGraphs(std::vector< RResultHandle > handles)
Run the event loops of multiple RDataFrames concurrently.
RInterface<::ROOT::Detail::RDF::RNodeBase > RNode
void EnableImplicitMT(UInt_t numthreads=0)
Enable ROOT's implicit multi-threading for all objects and methods that provide an internal paralleli...
Definition TROOT.cxx:544
void Sort(Index n, const Element *a, Index *index, Bool_t down=kTRUE)
Sort the n elements of the array a of generic templated type Element.
Definition TMathBase.h:413