Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
df103_NanoAODHiggsAnalysis_python.h
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_dataframe
3 /// Header file with functions needed to execute the Python version
4 /// of the NanoAOD Higgs tutorial. The header is declared to the
5 /// ROOT C++ interpreter prior to the start of the analysis via the
6 /// `ROOT.gInterpreter.Declare()` function.
7 ///
8 /// \date July 2019
9 /// \author Stefan Wunsch (KIT, CERN), Vincenzo Eduardo Padulano (UniMiB, CERN)
10 
11 #include "ROOT/RDataFrame.hxx"
12 #include "ROOT/RVec.hxx"
13 #include "TCanvas.h"
14 #include "TH1D.h"
15 #include "TLatex.h"
16 #include "Math/Vector4D.h"
17 #include "TStyle.h"
18 
19 using namespace ROOT::VecOps;
20 using RNode = ROOT::RDF::RNode;
21 using rvec_f = const RVec<float> &;
22 using rvec_i = const RVec<int> &;
23 const auto z_mass = 91.2;
24 
25 // Reconstruct two Z candidates from four leptons of the same kind
26 RVec<RVec<size_t>> reco_zz_to_4l(rvec_f pt, rvec_f eta, rvec_f phi, rvec_f mass, rvec_i charge)
27 {
28  RVec<RVec<size_t>> idx(2);
29  idx[0].reserve(2); idx[1].reserve(2);
30 
31  // Find first lepton pair with invariant mass closest to Z mass
32  auto idx_cmb = Combinations(pt, 2);
33  auto best_mass = -1;
34  size_t best_i1 = 0; size_t best_i2 = 0;
35  for (size_t i = 0; i < idx_cmb[0].size(); i++) {
36  const auto i1 = idx_cmb[0][i];
37  const auto i2 = idx_cmb[1][i];
38  if (charge[i1] != charge[i2]) {
39  ROOT::Math::PtEtaPhiMVector p1(pt[i1], eta[i1], phi[i1], mass[i1]);
40  ROOT::Math::PtEtaPhiMVector p2(pt[i2], eta[i2], phi[i2], mass[i2]);
41  const auto this_mass = (p1 + p2).M();
42  if (std::abs(z_mass - this_mass) < std::abs(z_mass - best_mass)) {
43  best_mass = this_mass;
44  best_i1 = i1;
45  best_i2 = i2;
46  }
47  }
48  }
49  idx[0].emplace_back(best_i1);
50  idx[0].emplace_back(best_i2);
51 
52  // Reconstruct second Z from remaining lepton pair
53  for (size_t i = 0; i < 4; i++) {
54  if (i != best_i1 && i != best_i2) {
55  idx[1].emplace_back(i);
56  }
57  }
58 
59  // Return indices of the pairs building two Z bosons
60  return idx;
61 }
62 
63 // Compute Z masses from four leptons of the same kind and sort ascending in distance to Z mass
64 RVec<float> compute_z_masses_4l(const RVec<RVec<size_t>> &idx, rvec_f pt, rvec_f eta, rvec_f phi, rvec_f mass)
65 {
66  RVec<float> z_masses(2);
67  for (size_t i = 0; i < 2; i++) {
68  const auto i1 = idx[i][0]; const auto i2 = idx[i][1];
69  ROOT::Math::PtEtaPhiMVector p1(pt[i1], eta[i1], phi[i1], mass[i1]);
70  ROOT::Math::PtEtaPhiMVector p2(pt[i2], eta[i2], phi[i2], mass[i2]);
71  z_masses[i] = (p1 + p2).M();
72  }
73  if (std::abs(z_masses[0] - z_mass) < std::abs(z_masses[1] - z_mass)) {
74  return z_masses;
75  } else {
76  return Reverse(z_masses);
77  }
78 }
79 
80 // Compute mass of Higgs from four leptons of the same kind
81 float compute_higgs_mass_4l(const RVec<RVec<size_t>> &idx, rvec_f pt, rvec_f eta, rvec_f phi, rvec_f mass)
82 {
83  const auto i1 = idx[0][0]; const auto i2 = idx[0][1];
84  const auto i3 = idx[1][0]; const auto i4 = idx[1][1];
85  ROOT::Math::PtEtaPhiMVector p1(pt[i1], eta[i1], phi[i1], mass[i1]);
86  ROOT::Math::PtEtaPhiMVector p2(pt[i2], eta[i2], phi[i2], mass[i2]);
87  ROOT::Math::PtEtaPhiMVector p3(pt[i3], eta[i3], phi[i3], mass[i3]);
88  ROOT::Math::PtEtaPhiMVector p4(pt[i4], eta[i4], phi[i4], mass[i4]);
89  return (p1 + p2 + p3 + p4).M();
90 }
91 
92 // Compute mass of two Z candidates from two electrons and two muons and sort ascending in distance to Z mass
93 RVec<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,
94  rvec_f mu_eta, rvec_f mu_phi, rvec_f mu_mass)
95 {
96  ROOT::Math::PtEtaPhiMVector p1(mu_pt[0], mu_eta[0], mu_phi[0], mu_mass[0]);
97  ROOT::Math::PtEtaPhiMVector p2(mu_pt[1], mu_eta[1], mu_phi[1], mu_mass[1]);
98  ROOT::Math::PtEtaPhiMVector p3(el_pt[0], el_eta[0], el_phi[0], el_mass[0]);
99  ROOT::Math::PtEtaPhiMVector p4(el_pt[1], el_eta[1], el_phi[1], el_mass[1]);
100  auto mu_z = (p1 + p2).M();
101  auto el_z = (p3 + p4).M();
102  RVec<float> z_masses(2);
103  if (std::abs(mu_z - z_mass) < std::abs(el_z - z_mass)) {
104  z_masses[0] = mu_z;
105  z_masses[1] = el_z;
106  } else {
107  z_masses[0] = el_z;
108  z_masses[1] = mu_z;
109  }
110  return z_masses;
111 }
112 
113 // Compute Higgs mass from two electrons and two muons
114 float 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,
115  rvec_f mu_phi, rvec_f mu_mass)
116 {
117  ROOT::Math::PtEtaPhiMVector p1(mu_pt[0], mu_eta[0], mu_phi[0], mu_mass[0]);
118  ROOT::Math::PtEtaPhiMVector p2(mu_pt[1], mu_eta[1], mu_phi[1], mu_mass[1]);
119  ROOT::Math::PtEtaPhiMVector p3(el_pt[0], el_eta[0], el_phi[0], el_mass[0]);
120  ROOT::Math::PtEtaPhiMVector p4(el_pt[1], el_eta[1], el_phi[1], el_mass[1]);
121  return (p1 + p2 + p3 + p4).M();
122 }
123 
124 bool filter_z_dr(const RVec<RVec<size_t>> &idx, rvec_f eta, rvec_f phi)
125 {
126  for (size_t i = 0; i < 2; i++) {
127  const auto i1 = idx[i][0];
128  const auto i2 = idx[i][1];
129  const auto dr = DeltaR(eta[i1], eta[i2], phi[i1], phi[i2]);
130  if (dr < 0.02) {
131  return false;
132  }
133  }
134  return true;
135 };
136 
137 bool pt_cuts(rvec_f mu_pt, rvec_f el_pt)
138 {
139  auto mu_pt_sorted = Reverse(Sort(mu_pt));
140  if (mu_pt_sorted[0] > 20 && mu_pt_sorted[1] > 10) {
141  return true;
142  }
143  auto el_pt_sorted = Reverse(Sort(el_pt));
144  if (el_pt_sorted[0] > 20 && el_pt_sorted[1] > 10) {
145  return true;
146  }
147  return false;
148 }
149 
150 bool dr_cuts(rvec_f mu_eta, rvec_f mu_phi, rvec_f el_eta, rvec_f el_phi)
151 {
152  auto mu_dr = DeltaR(mu_eta[0], mu_eta[1], mu_phi[0], mu_phi[1]);
153  auto el_dr = DeltaR(el_eta[0], el_eta[1], el_phi[0], el_phi[1]);
154  if (mu_dr < 0.02 || el_dr < 0.02) {
155  return false;
156  }
157  return true;
158 }