Logo ROOT   6.16/01
Reference Guide
df103_NanoAODHiggsAnalysis.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_dataframe
3/// \notebook -draw
4/// This tutorial is a simplified but yet complex example of an analysis reconstructing
5/// the Higgs boson decaying to two Z bosons from events with four leptons. The data
6/// and simulated events are taken from CERN OpenData representing a subset of the data
7/// recorded in 2012 with the CMS detector at the LHC. The tutorials follows the Higgs
8/// to four leptons analysis published on CERN Open Data portal
9/// ([10.7483/OPENDATA.CMS.JKB8.RR42](http://opendata.cern.ch/record/5500)).
10/// The resulting plots show the invariant mass of the selected four lepton systems
11/// in different decay modes (four muons, four electrons and two of each kind)
12/// and in a combined plot indicating the decay of the Higgs boson with a mass
13/// of about 125 GeV.
14///
15/// The following steps are performed for each sample with data and simulated events
16/// in order to reconstruct the Higgs boson from the selected muons and electrons:
17/// 1. Select interesting events with multiple cuts on event properties, e.g.,
18/// number of leptons, kinematics of the letpons and quality of the tracks.
19/// 2. Reconstruct two Z bosons of which only one on the mass shell from the selected events and apply additional cuts
20/// on the reconstructed objects.
21/// 3. Reconstruct the Higgs boson from the remaining Z boson candidates and calculate
22/// its invariant mass.
23///
24/// \macro_image
25/// \macro_code
26/// \macro_output
27///
28/// \date October 2018
29/// \author Stefan Wunsch (KIT, CERN)
30
31#include "ROOT/RDataFrame.hxx"
32#include "ROOT/RVec.hxx"
34#include "TCanvas.h"
35#include "TH1D.h"
36#include "TLatex.h"
37#include "TLegend.h"
38#include "TLorentzVector.h"
39#include "TStyle.h"
40
41using namespace ROOT::VecOps;
43using rvec_f = const RVec<float> &;
44using rvec_i = const RVec<int> &;
45
46const auto z_mass = 91.2;
47
48// Select interesting events with four muons
49RNode selection_4mu(RNode df)
50{
51 auto df_ge4m = df.Filter("nMuon>=4", "At least four muons");
52 auto df_iso = df_ge4m.Filter("All(abs(Muon_pfRelIso04_all)<0.40)", "Require good isolation");
53 auto df_kin = df_iso.Filter("All(Muon_pt>5) && All(abs(Muon_eta)<2.4)", "Good muon kinematics");
54 auto df_ip3d = df_kin.Define("Muon_ip3d", "sqrt(Muon_dxy*Muon_dxy + Muon_dz*Muon_dz)");
55 auto df_sip3d = df_ip3d.Define("Muon_sip3d", "Muon_ip3d/sqrt(Muon_dxyErr*Muon_dxyErr + Muon_dzErr*Muon_dzErr)");
56 auto df_pv = df_sip3d.Filter("All(Muon_sip3d<4) && All(abs(Muon_dxy)<0.5) && All(abs(Muon_dz)<1.0)",
57 "Track close to primary vertex with small uncertainty");
58 auto df_2p2n = df_pv.Filter("nMuon==4 && Sum(Muon_charge==1)==2 && Sum(Muon_charge==-1)==2",
59 "Two positive and two negative muons");
60 return df_2p2n;
61}
62
63// Select interesting events with four electrons
64RNode selection_4el(RNode df)
65{
66 auto df_ge4el = df.Filter("nElectron>=4", "At least our electrons");
67 auto df_iso = df_ge4el.Filter("All(abs(Electron_pfRelIso03_all)<0.40)", "Require good isolation");
68 auto df_kin = df_iso.Filter("All(Electron_pt>7) && All(abs(Electron_eta)<2.5)", "Good Electron kinematics");
69 auto df_ip3d = df_kin.Define("Electron_ip3d", "sqrt(Electron_dxy*Electron_dxy + Electron_dz*Electron_dz)");
70 auto df_sip3d = df_ip3d.Define("Electron_sip3d",
71 "Electron_ip3d/sqrt(Electron_dxyErr*Electron_dxyErr + Electron_dzErr*Electron_dzErr)");
72 auto df_pv = df_sip3d.Filter("All(Electron_sip3d<4) && All(abs(Electron_dxy)<0.5) && "
73 "All(abs(Electron_dz)<1.0)",
74 "Track close to primary vertex with small uncertainty");
75 auto df_2p2n = df_pv.Filter("nElectron==4 && Sum(Electron_charge==1)==2 && Sum(Electron_charge==-1)==2",
76 "Two positive and two negative electrons");
77 return df_2p2n;
78}
79
80// Select interesting events with two electrons and two muons
81RNode selection_2el2mu(RNode df)
82{
83 auto df_ge2el2mu = df.Filter("nElectron>=2 && nMuon>=2", "At least two electrons and two muons");
84 auto df_eta = df_ge2el2mu.Filter("All(abs(Electron_eta)<2.5) && All(abs(Muon_eta)<2.4)", "Eta cuts");
85 auto pt_cuts = [](rvec_f mu_pt, rvec_f el_pt) {
86 auto mu_pt_sorted = Reverse(Sort(mu_pt));
87 if (mu_pt_sorted[0] > 20 && mu_pt_sorted[1] > 10) {
88 return true;
89 }
90 auto el_pt_sorted = Reverse(Sort(el_pt));
91 if (el_pt_sorted[0] > 20 && el_pt_sorted[1] > 10) {
92 return true;
93 }
94 return false;
95 };
96 auto df_pt = df_eta.Filter(pt_cuts, {"Muon_pt", "Electron_pt"}, "Pt cuts");
97 auto dr_cuts = [](rvec_f mu_eta, rvec_f mu_phi, rvec_f el_eta, rvec_f el_phi) {
98 auto mu_dr = sqrt(pow(mu_eta[0] - mu_eta[1], 2) + pow(mu_phi[0] - mu_phi[1], 2));
99 auto el_dr = sqrt(pow(el_eta[0] - el_eta[1], 2) + pow(el_phi[0] - el_phi[1], 2));
100 if (mu_dr < 0.02 || el_dr < 0.02) {
101 return false;
102 }
103 return true;
104 };
105 auto df_dr = df_pt.Filter(dr_cuts, {"Muon_eta", "Muon_phi", "Electron_eta", "Electron_phi"}, "Dr cuts");
106 auto df_iso = df_dr.Filter("All(abs(Electron_pfRelIso03_all)<0.40) && All(abs(Muon_pfRelIso04_all)<0.40)",
107 "Require good isolation");
108 auto df_el_ip3d = df_iso.Define("Electron_ip3d_el", "sqrt(Electron_dxy*Electron_dxy + Electron_dz*Electron_dz)");
109 auto df_el_sip3d = df_el_ip3d.Define("Electron_sip3d_el",
110 "Electron_ip3d_el/sqrt(Electron_dxyErr*Electron_dxyErr + "
111 "Electron_dzErr*Electron_dzErr)");
112 auto df_el_track = df_el_sip3d.Filter("All(Electron_sip3d_el<4) && All(abs(Electron_dxy)<0.5) && All(abs(Electron_dz)<1.0)",
113 "Electron track close to primary vertex with small uncertainty");
114 auto df_mu_ip3d = df_el_track.Define("Muon_ip3d_mu", "sqrt(Muon_dxy*Muon_dxy + Muon_dz*Muon_dz)");
115 auto df_mu_sip3d = df_mu_ip3d.Define("Muon_sip3d_mu",
116 "Muon_ip3d_mu/sqrt(Muon_dxyErr*Muon_dxyErr + Muon_dzErr*Muon_dzErr)");
117 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)",
118 "Muon track close to primary vertex with small uncertainty");
119 auto df_2p2n = df_mu_track.Filter("Sum(Electron_charge)==0 && Sum(Muon_charge)==0",
120 "Two opposite charged electron and muon pairs");
121 return df_2p2n;
122}
123
124// Reconstruct two Z candidates from four leptons of the same kind
125RVec<RVec<size_t>> reco_zz_to_4l(rvec_f pt, rvec_f eta, rvec_f phi, rvec_f mass, rvec_i charge)
126{
127 RVec<RVec<size_t>> idx(2);
128 idx[0].reserve(2); idx[1].reserve(2);
129
130 // Find first lepton pair with invariant mass closest to Z mass
131 auto idx_cmb = Combinations(pt, 2);
132 auto best_mass = -1;
133 size_t best_i1 = 0; size_t best_i2 = 0;
134 for (size_t i = 0; i < idx_cmb[0].size(); i++) {
135 const auto i1 = idx_cmb[0][i];
136 const auto i2 = idx_cmb[1][i];
137 if (charge[i1] != charge[i2]) {
139 p1.SetPtEtaPhiM(pt[i1], eta[i1], phi[i1], mass[i1]);
140 p2.SetPtEtaPhiM(pt[i2], eta[i2], phi[i2], mass[i2]);
141 const auto this_mass = (p1 + p2).M();
142 if (std::abs(z_mass - this_mass) < std::abs(z_mass - best_mass)) {
143 best_mass = this_mass;
144 best_i1 = i1;
145 best_i2 = i2;
146 }
147 }
148 }
149 idx[0].emplace_back(best_i1);
150 idx[0].emplace_back(best_i2);
151
152 // Reconstruct second Z from remaining lepton pair
153 for (size_t i = 0; i < 4; i++) {
154 if (i != best_i1 && i != best_i2) {
155 idx[1].emplace_back(i);
156 }
157 }
158
159 // Return indices of the pairs building two Z bosons
160 return idx;
161}
162
163// Compute Z masses from four leptons of the same kind and sort ascending in distance to Z mass
164RVec<float> compute_z_masses_4l(const RVec<RVec<size_t>> &idx, rvec_f pt, rvec_f eta, rvec_f phi, rvec_f mass)
165{
166 RVec<float> z_masses(2);
167 for (size_t i = 0; i < 2; i++) {
169 const auto i1 = idx[i][0]; const auto i2 = idx[i][1];
170 p1.SetPtEtaPhiM(pt[i1], eta[i1], phi[i1], mass[i1]);
171 p2.SetPtEtaPhiM(pt[i2], eta[i2], phi[i2], mass[i2]);
172 z_masses[i] = (p1 + p2).M();
173 }
174 if (std::abs(z_masses[0] - z_mass) < std::abs(z_masses[1] - z_mass)) {
175 return z_masses;
176 } else {
177 return Reverse(z_masses);
178 }
179}
180
181// Compute mass of Higgs from four leptons of the same kind
182float compute_higgs_mass_4l(const RVec<RVec<size_t>> &idx, rvec_f pt, rvec_f eta, rvec_f phi, rvec_f mass)
183{
184 TLorentzVector p1, p2, p3, p4;
185 const auto i1 = idx[0][0]; const auto i2 = idx[0][1];
186 const auto i3 = idx[1][0]; const auto i4 = idx[1][1];
187 p1.SetPtEtaPhiM(pt[i1], eta[i1], phi[i1], mass[i1]);
188 p2.SetPtEtaPhiM(pt[i2], eta[i2], phi[i2], mass[i2]);
189 p3.SetPtEtaPhiM(pt[i3], eta[i3], phi[i3], mass[i3]);
190 p4.SetPtEtaPhiM(pt[i4], eta[i4], phi[i4], mass[i4]);
191 return (p1 + p2 + p3 + p4).M();
192}
193
194// Apply selection on reconstructed Z candidates
195RNode filter_z_candidates(RNode df)
196{
197 auto df_z1_cut = df.Filter("Z_mass[0] > 40 && Z_mass[0] < 120", "Mass of first Z candidate in [40, 120]");
198 auto df_z2_cut = df_z1_cut.Filter("Z_mass[1] > 12 && Z_mass[1] < 120", "Mass of second Z candidate in [12, 120]");
199 return df_z2_cut;
200}
201
202// Reconstruct Higgs from four muons
203RNode reco_higgs_to_4mu(RNode df)
204{
205 // Filter interesting events
206 auto df_base = selection_4mu(df);
207
208 // Reconstruct Z systems
209 auto df_z_idx =
210 df_base.Define("Z_idx", reco_zz_to_4l, {"Muon_pt", "Muon_eta", "Muon_phi", "Muon_mass", "Muon_charge"});
211
212 // Cut on distance between muons building Z systems
213 auto filter_z_dr = [](const RVec<RVec<size_t>> &idx, rvec_f eta, rvec_f phi) {
214 for (size_t i = 0; i < 2; i++) {
215 const auto i1 = idx[i][0];
216 const auto i2 = idx[i][1];
217 const auto dr = sqrt(pow(eta[i1] - eta[i2], 2) + pow(phi[i1] - phi[i2], 2));
218 if (dr < 0.02) {
219 return false;
220 }
221 }
222 return true;
223 };
224 auto df_z_dr =
225 df_z_idx.Filter(filter_z_dr, {"Z_idx", "Muon_eta", "Muon_phi"}, "Delta R separation of muons building Z system");
226
227 // Compute masses of Z systems
228 auto df_z_mass =
229 df_z_dr.Define("Z_mass", compute_z_masses_4l, {"Z_idx", "Muon_pt", "Muon_eta", "Muon_phi", "Muon_mass"});
230
231 // Cut on mass of Z candidates
232 auto df_z_cut = filter_z_candidates(df_z_mass);
233
234 // Reconstruct H mass
235 auto df_h_mass =
236 df_z_cut.Define("H_mass", compute_higgs_mass_4l, {"Z_idx", "Muon_pt", "Muon_eta", "Muon_phi", "Muon_mass"});
237
238 return df_h_mass;
239}
240
241// Reconstruct Higgs from four electrons
242RNode reco_higgs_to_4el(RNode df)
243{
244 // Filter interesting events
245 auto df_base = selection_4el(df);
246
247 // Reconstruct Z systems
248 auto df_z_idx = df_base.Define("Z_idx", reco_zz_to_4l,
249 {"Electron_pt", "Electron_eta", "Electron_phi", "Electron_mass", "Electron_charge"});
250
251 // Cut on distance between Electrons building Z systems
252 auto filter_z_dr = [](const RVec<RVec<size_t>> &idx, rvec_f eta, rvec_f phi) {
253 for (size_t i = 0; i < 2; i++) {
254 const auto i1 = idx[i][0];
255 const auto i2 = idx[i][1];
256 const auto dr = sqrt(pow(eta[i1] - eta[i2], 2) + pow(phi[i1] - phi[i2], 2));
257 if (dr < 0.02) {
258 return false;
259 }
260 }
261 return true;
262 };
263 auto df_z_dr = df_z_idx.Filter(filter_z_dr, {"Z_idx", "Electron_eta", "Electron_phi"},
264 "Delta R separation of Electrons building Z system");
265
266 // Compute masses of Z systems
267 auto df_z_mass = df_z_dr.Define("Z_mass", compute_z_masses_4l,
268 {"Z_idx", "Electron_pt", "Electron_eta", "Electron_phi", "Electron_mass"});
269
270 // Cut on mass of Z candidates
271 auto df_z_cut = filter_z_candidates(df_z_mass);
272
273 // Reconstruct H mass
274 auto df_h_mass = df_z_cut.Define("H_mass", compute_higgs_mass_4l,
275 {"Z_idx", "Electron_pt", "Electron_eta", "Electron_phi", "Electron_mass"});
276
277 return df_h_mass;
278}
279
280// Compute mass of two Z candidates from two electrons and two muons and sort ascending in distance to Z mass
281RVec<float> compute_z_masses_2el2mu(rvec_f el_pt, rvec_f el_eta, rvec_f el_phi, rvec_f el_mass, rvec_f mu_pt,
282 rvec_f mu_eta, rvec_f mu_phi, rvec_f mu_mass)
283{
284 TLorentzVector p1, p2, p3, p4;
285 p1.SetPtEtaPhiM(mu_pt[0], mu_eta[0], mu_phi[0], mu_mass[0]);
286 p2.SetPtEtaPhiM(mu_pt[1], mu_eta[1], mu_phi[1], mu_mass[1]);
287 p3.SetPtEtaPhiM(el_pt[0], el_eta[0], el_phi[0], el_mass[0]);
288 p4.SetPtEtaPhiM(el_pt[1], el_eta[1], el_phi[1], el_mass[1]);
289 auto mu_z = (p1 + p2).M();
290 auto el_z = (p3 + p4).M();
291 RVec<float> z_masses(2);
292 if (std::abs(mu_z - z_mass) < std::abs(el_z - z_mass)) {
293 z_masses[0] = mu_z;
294 z_masses[1] = el_z;
295 } else {
296 z_masses[0] = el_z;
297 z_masses[1] = mu_z;
298 }
299 return z_masses;
300}
301
302// Compute Higgs mass from two electrons and two muons
303float compute_higgs_mass_2el2mu(rvec_f el_pt, rvec_f el_eta, rvec_f el_phi, rvec_f el_mass, rvec_f mu_pt, rvec_f mu_eta,
304 rvec_f mu_phi, rvec_f mu_mass)
305{
306 TLorentzVector p1, p2, p3, p4;
307 p1.SetPtEtaPhiM(mu_pt[0], mu_eta[0], mu_phi[0], mu_mass[0]);
308 p2.SetPtEtaPhiM(mu_pt[1], mu_eta[1], mu_phi[1], mu_mass[1]);
309 p3.SetPtEtaPhiM(el_pt[0], el_eta[0], el_phi[0], el_mass[0]);
310 p4.SetPtEtaPhiM(el_pt[1], el_eta[1], el_phi[1], el_mass[1]);
311 return (p1 + p2 + p3 + p4).M();
312}
313
314// Reconstruct Higgs from two electrons and two muons
315RNode reco_higgs_to_2el2mu(RNode df)
316{
317 // Filter interesting events
318 auto df_base = selection_2el2mu(df);
319
320 // Compute masses of Z systems
321 auto df_z_mass =
322 df_base.Define("Z_mass", compute_z_masses_2el2mu, {"Electron_pt", "Electron_eta", "Electron_phi", "Electron_mass",
323 "Muon_pt", "Muon_eta", "Muon_phi", "Muon_mass"});
324
325 // Cut on mass of Z candidates
326 auto df_z_cut = filter_z_candidates(df_z_mass);
327
328 // Reconstruct H mass
329 auto df_h_mass = df_z_cut.Define(
330 "H_mass", compute_higgs_mass_2el2mu,
331 {"Electron_pt", "Electron_eta", "Electron_phi", "Electron_mass", "Muon_pt", "Muon_eta", "Muon_phi", "Muon_mass"});
332
333 return df_h_mass;
334}
335
336// Plot invariant mass for signal and background processes from simulated events
337// overlay the measured data.
338template <typename T>
339void plot(T sig, T bkg, T data, const std::string &x_label, const std::string &filename)
340{
341 // Canvas and general style options
342 gStyle->SetOptStat(0);
343 gStyle->SetTextFont(42);
344 auto c = new TCanvas("c", "", 800, 700);
345 c->SetLeftMargin(0.15);
346
347 // Get signal and background histograms and stack them to show Higgs signal
348 // on top of the background process
349 auto h_sig = *sig;
350 auto h_bkg = *bkg;
351 auto h_cmb = *(TH1D*)(sig->Clone());
352 h_cmb.Add(&h_bkg);
353 h_cmb.SetTitle("");
354 h_cmb.GetXaxis()->SetTitle(x_label.c_str());
355 h_cmb.GetXaxis()->SetTitleSize(0.04);
356 h_cmb.GetYaxis()->SetTitle("N_{Events}");
357 h_cmb.GetYaxis()->SetTitleSize(0.04);
358 h_cmb.SetLineColor(kRed);
359 h_cmb.SetLineWidth(2);
360 h_cmb.SetMaximum(18);
361
362 h_bkg.SetLineWidth(2);
363 h_bkg.SetFillStyle(1001);
364 h_bkg.SetLineColor(kBlack);
365 h_bkg.SetFillColor(kAzure - 9);
366
367 // Get histogram of data points
368 auto h_data = *data;
369 h_data.SetLineWidth(1);
370 h_data.SetMarkerStyle(20);
371 h_data.SetMarkerSize(1.0);
372 h_data.SetMarkerColor(kBlack);
373 h_data.SetLineColor(kBlack);
374
375 // Draw histograms
376 h_cmb.DrawClone("HIST");
377 h_bkg.DrawClone("HIST SAME");
378 h_data.DrawClone("PE1 SAME");
379
380 // Add legend
381 TLegend legend(0.62, 0.70, 0.82, 0.88);
382 legend.SetFillColor(0);
383 legend.SetBorderSize(0);
384 legend.SetTextSize(0.03);
385 legend.AddEntry(&h_data, "Data", "PE1");
386 legend.AddEntry(&h_bkg, "ZZ", "f");
387 legend.AddEntry(&h_cmb, "m_{H} = 125 GeV", "f");
388 legend.Draw();
389
390 // Add header
391 TLatex cms_label;
392 cms_label.SetTextSize(0.04);
393 cms_label.DrawLatexNDC(0.16, 0.92, "#bf{CMS Open Data}");
394 TLatex header;
395 header.SetTextSize(0.03);
396 header.DrawLatexNDC(0.63, 0.92, "#sqrt{s} = 8 TeV, L_{int} = 11.6 fb^{-1}");
397
398 // Save plot
399 c->SaveAs(filename.c_str());
400}
401
402void df103_NanoAODHiggsAnalysis()
403{
404 // Enable multi-threading
406
407 // Create dataframes for signal, background and data samples
408
409 // Signal: Higgs -> 4 leptons
410 ROOT::RDataFrame df_sig_4l("Events",
411 "root://eospublic.cern.ch//eos/root-eos/cms_opendata_2012_nanoaod/SMHiggsToZZTo4L.root");
412
413 // Background: ZZ -> 4 leptons
414 // Note that additional background processes from the original paper with minor contribution were left out for this
415 // tutorial.
416 ROOT::RDataFrame df_bkg_4mu("Events",
417 "root://eospublic.cern.ch//eos/root-eos/cms_opendata_2012_nanoaod/ZZTo4mu.root");
418 ROOT::RDataFrame df_bkg_4el("Events",
419 "root://eospublic.cern.ch//eos/root-eos/cms_opendata_2012_nanoaod/ZZTo4e.root");
420 ROOT::RDataFrame df_bkg_2el2mu("Events",
421 "root://eospublic.cern.ch//eos/root-eos/cms_opendata_2012_nanoaod/ZZTo2e2mu.root");
422
423 // CMS data taken in 2012 (11.6 fb^-1 integrated luminosity)
424 ROOT::RDataFrame df_data_doublemu(
425 "Events", {"root://eospublic.cern.ch//eos/root-eos/cms_opendata_2012_nanoaod/Run2012B_DoubleMuParked.root",
426 "root://eospublic.cern.ch//eos/root-eos/cms_opendata_2012_nanoaod/Run2012C_DoubleMuParked.root"});
427 ROOT::RDataFrame df_data_doubleel(
428 "Events", {"root://eospublic.cern.ch//eos/root-eos/cms_opendata_2012_nanoaod/Run2012B_DoubleElectron.root",
429 "root://eospublic.cern.ch//eos/root-eos/cms_opendata_2012_nanoaod/Run2012C_DoubleElectron.root"});
430
431 // Reconstruct Higgs to 4 muons
432 auto df_sig_4mu_reco = reco_higgs_to_4mu(df_sig_4l);
433 const auto luminosity = 11580.0; // Integrated luminosity of the data samples
434 const auto xsec_SMHiggsToZZTo4L = 0.0065; // H->4l: Standard Model cross-section
435 const auto nevt_SMHiggsToZZTo4L = 299973.0; // H->4l: Number of simulated events
436 const auto nbins = 36; // Number of bins for the invariant mass spectrum
437 auto df_h_sig_4mu = df_sig_4mu_reco
438 .Define("weight", [&]() { return luminosity * xsec_SMHiggsToZZTo4L / nevt_SMHiggsToZZTo4L; }, {})
439 .Histo1D({"h_sig_4mu", "", nbins, 70, 180}, "H_mass", "weight");
440
441 const auto scale_ZZTo4l = 1.386; // ZZ->4mu: Scale factor for ZZ to four leptons
442 const auto xsec_ZZTo4mu = 0.077; // ZZ->4mu: Standard Model cross-section
443 const auto nevt_ZZTo4mu = 1499064.0; // ZZ->4mu: Number of simulated events
444 auto df_bkg_4mu_reco = reco_higgs_to_4mu(df_bkg_4mu);
445 auto df_h_bkg_4mu = df_bkg_4mu_reco
446 .Define("weight", [&]() { return luminosity * xsec_ZZTo4mu * scale_ZZTo4l / nevt_ZZTo4mu; }, {})
447 .Histo1D({"h_bkg_4mu", "", nbins, 70, 180}, "H_mass", "weight");
448
449 auto df_data_4mu_reco = reco_higgs_to_4mu(df_data_doublemu);
450 auto df_h_data_4mu = df_data_4mu_reco
451 .Define("weight", []() { return 1.0; }, {})
452 .Histo1D({"h_data_4mu", "", nbins, 70, 180}, "H_mass", "weight");
453
454 // Reconstruct Higgs to 4 electrons
455 auto df_sig_4el_reco = reco_higgs_to_4el(df_sig_4l);
456 auto df_h_sig_4el = df_sig_4el_reco
457 .Define("weight", [&]() { return luminosity * xsec_SMHiggsToZZTo4L / nevt_SMHiggsToZZTo4L; }, {})
458 .Histo1D({"h_sig_4el", "", nbins, 70, 180}, "H_mass", "weight");
459
460 const auto xsec_ZZTo4el = xsec_ZZTo4mu; // ZZ->4el: Standard Model cross-section
461 const auto nevt_ZZTo4el = 1499093.0; // ZZ->4el: Number of simulated events
462 auto df_bkg_4el_reco = reco_higgs_to_4el(df_bkg_4el);
463 auto df_h_bkg_4el = df_bkg_4el_reco
464 .Define("weight", [&]() { return luminosity * xsec_ZZTo4el * scale_ZZTo4l / nevt_ZZTo4el; }, {})
465 .Histo1D({"h_bkg_4el", "", nbins, 70, 180}, "H_mass", "weight");
466
467 auto df_data_4el_reco = reco_higgs_to_4el(df_data_doubleel);
468 auto df_h_data_4el = df_data_4el_reco.Define("weight", []() { return 1.0; }, {})
469 .Histo1D({"h_data_4el", "", nbins, 70, 180}, "H_mass", "weight");
470
471 // Reconstruct Higgs to 2 electrons and 2 muons
472 auto df_sig_2el2mu_reco = reco_higgs_to_2el2mu(df_sig_4l);
473 auto df_h_sig_2el2mu = df_sig_2el2mu_reco
474 .Define("weight", [&]() { return luminosity * xsec_SMHiggsToZZTo4L / nevt_SMHiggsToZZTo4L; }, {})
475 .Histo1D({"h_sig_2el2mu", "", nbins, 70, 180}, "H_mass", "weight");
476
477 const auto xsec_ZZTo2el2mu = 0.18; // ZZ->2el2mu: Standard Model cross-section
478 const auto nevt_ZZTo2el2mu = 1497445.0; // ZZ->2el2mu: Number of simulated events
479 auto df_bkg_2el2mu_reco = reco_higgs_to_2el2mu(df_bkg_2el2mu);
480 auto df_h_bkg_2el2mu = df_bkg_2el2mu_reco
481 .Define("weight", [&]() { return luminosity * xsec_ZZTo2el2mu * scale_ZZTo4l / nevt_ZZTo2el2mu; }, {})
482 .Histo1D({"h_bkg_2el2mu", "", nbins, 70, 180}, "H_mass", "weight");
483
484 auto df_data_2el2mu_reco = reco_higgs_to_2el2mu(df_data_doublemu);
485 auto df_h_data_2el2mu = df_data_2el2mu_reco.Define("weight", []() { return 1.0; }, {})
486 .Histo1D({"h_data_2el2mu_doublemu", "", nbins, 70, 180}, "H_mass", "weight");
487
488 // Produce histograms for different channels and make plots
489 plot(df_h_sig_4mu, df_h_bkg_4mu, df_h_data_4mu, "m_{4#mu} (GeV)", "higgs_4mu.pdf");
490 plot(df_h_sig_4el, df_h_bkg_4el, df_h_data_4el, "m_{4e} (GeV)", "higgs_4el.pdf");
491 plot(df_h_sig_2el2mu, df_h_bkg_2el2mu, df_h_data_2el2mu, "m_{2e2#mu} (GeV)", "higgs_2el2mu.pdf");
492
493 // Combine channels for final plot
494 auto h_data_4l = df_h_data_4mu.GetPtr();
495 h_data_4l->Add(df_h_data_4el.GetPtr());
496 h_data_4l->Add(df_h_data_2el2mu.GetPtr());
497 auto h_sig_4l = df_h_sig_4mu.GetPtr();
498 h_sig_4l->Add(df_h_sig_4el.GetPtr());
499 h_sig_4l->Add(df_h_sig_2el2mu.GetPtr());
500 auto h_bkg_4l = df_h_bkg_4mu.GetPtr();
501 h_bkg_4l->Add(df_h_bkg_4el.GetPtr());
502 h_bkg_4l->Add(df_h_bkg_2el2mu.GetPtr());
503 plot(h_sig_4l, h_bkg_4l, h_data_4l, "m_{4l} (GeV)", "higgs_4l.pdf");
504}
505
506int main()
507{
508 df103_NanoAODHiggsAnalysis();
509}
#define c(i)
Definition: RSha256.hxx:101
static double p3(double t, double a, double b, double c, double d)
static double p1(double t, double a, double b)
static double p2(double t, double a, double b, double c)
@ kRed
Definition: Rtypes.h:63
@ kBlack
Definition: Rtypes.h:62
@ kAzure
Definition: Rtypes.h:64
double pow(double, double)
double sqrt(double)
R__EXTERN TStyle * gStyle
Definition: TStyle.h:406
ROOT's RDataFrame offers a high level interface for analyses of data stored in TTrees,...
Definition: RDataFrame.hxx:41
A "std::vector"-like collection of values implementing handy operation to analyse them.
Definition: RVec.hxx:221
virtual void SetTextFont(Font_t tfont=62)
Set the text font.
Definition: TAttText.h:45
virtual void SetTextSize(Float_t tsize=1)
Set the text size.
Definition: TAttText.h:46
The Canvas class.
Definition: TCanvas.h:31
1-D histogram with a double per channel (see TH1 documentation)}
Definition: TH1.h:614
To draw Mathematical Formula.
Definition: TLatex.h:18
TLatex * DrawLatexNDC(Double_t x, Double_t y, const char *text)
Draw this TLatex with new coordinates in NDC.
Definition: TLatex.cxx:1927
This class displays a legend box (TPaveText) containing several legend entries.
Definition: TLegend.h:23
TLorentzVector is a general four-vector class, which can be used either for the description of positi...
void SetPtEtaPhiM(Double_t pt, Double_t eta, Double_t phi, Double_t m)
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:1444
TPaveText * pt
int main(int argc, char **argv)
double T(double x)
Definition: ChebyshevPol.h:34
RInterface<::ROOT::Detail::RDF::RNodeBase, void > RNode
RVec< T > Reverse(const RVec< T > &v)
Return copy of reversed vector.
Definition: RVec.hxx:940
RVec< RVec< typename RVec< T1 >::size_type > > Combinations(const RVec< T1 > &v1, const RVec< T2 > &v2)
Return the indices that represent all combinations of the elements of two RVecs.
Definition: RVec.hxx:1008
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:576
void Sort(Index n, const Element *a, Index *index, Bool_t down=kTRUE)
Definition: TMathBase.h:362