Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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/// \authors 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
19using namespace ROOT;
20using namespace ROOT::VecOps;
22const auto z_mass = 91.2;
23
24// Reconstruct two Z candidates from four leptons of the same kind
25RVec<RVec<size_t>> reco_zz_to_4l(RVecF pt, RVecF eta, RVecF phi, RVecF mass, RVecI charge)
26{
27 RVec<RVec<size_t>> idx(2);
28 idx[0].reserve(2); idx[1].reserve(2);
29
30 // Find first lepton pair with invariant mass closest to Z mass
31 auto idx_cmb = Combinations(pt, 2);
32 auto best_mass = -1;
33 size_t best_i1 = 0; size_t best_i2 = 0;
34 for (size_t i = 0; i < idx_cmb[0].size(); i++) {
35 const auto i1 = idx_cmb[0][i];
36 const auto i2 = idx_cmb[1][i];
37 if (charge[i1] != charge[i2]) {
38 ROOT::Math::PtEtaPhiMVector p1(pt[i1], eta[i1], phi[i1], mass[i1]);
39 ROOT::Math::PtEtaPhiMVector p2(pt[i2], eta[i2], phi[i2], mass[i2]);
40 const auto this_mass = (p1 + p2).M();
41 if (std::abs(z_mass - this_mass) < std::abs(z_mass - best_mass)) {
42 best_mass = this_mass;
43 best_i1 = i1;
44 best_i2 = i2;
45 }
46 }
47 }
48 idx[0].emplace_back(best_i1);
49 idx[0].emplace_back(best_i2);
50
51 // Reconstruct second Z from remaining lepton pair
52 for (size_t i = 0; i < 4; i++) {
53 if (i != best_i1 && i != best_i2) {
54 idx[1].emplace_back(i);
55 }
56 }
57
58 // Return indices of the pairs building two Z bosons
59 return idx;
60}
61
62// Compute Z masses from four leptons of the same kind and sort ascending in distance to Z mass
63RVecF compute_z_masses_4l(const RVec<RVec<size_t>> &idx, RVecF pt, RVecF eta, RVecF phi, RVecF mass)
64{
65 RVecF z_masses(2);
66 for (size_t i = 0; i < 2; i++) {
67 const auto i1 = idx[i][0]; const auto i2 = idx[i][1];
68 ROOT::Math::PtEtaPhiMVector p1(pt[i1], eta[i1], phi[i1], mass[i1]);
69 ROOT::Math::PtEtaPhiMVector p2(pt[i2], eta[i2], phi[i2], mass[i2]);
70 z_masses[i] = (p1 + p2).M();
71 }
72 if (std::abs(z_masses[0] - z_mass) < std::abs(z_masses[1] - z_mass)) {
73 return z_masses;
74 } else {
75 return Reverse(z_masses);
76 }
77}
78
79// Compute mass of Higgs from four leptons of the same kind
80float compute_higgs_mass_4l(const RVec<RVec<size_t>> &idx, RVecF pt, RVecF eta, RVecF phi, RVecF mass)
81{
82 const auto i1 = idx[0][0]; const auto i2 = idx[0][1];
83 const auto i3 = idx[1][0]; const auto i4 = idx[1][1];
84 ROOT::Math::PtEtaPhiMVector p1(pt[i1], eta[i1], phi[i1], mass[i1]);
85 ROOT::Math::PtEtaPhiMVector p2(pt[i2], eta[i2], phi[i2], mass[i2]);
86 ROOT::Math::PtEtaPhiMVector p3(pt[i3], eta[i3], phi[i3], mass[i3]);
87 ROOT::Math::PtEtaPhiMVector p4(pt[i4], eta[i4], phi[i4], mass[i4]);
88 return (p1 + p2 + p3 + p4).M();
89}
90
91// Compute mass of two Z candidates from two electrons and two muons and sort ascending in distance to Z mass
92RVecF compute_z_masses_2el2mu(RVecF el_pt, RVecF el_eta, RVecF el_phi, RVecF el_mass, RVecF mu_pt,
93 RVecF mu_eta, RVecF mu_phi, RVecF mu_mass)
94{
95 ROOT::Math::PtEtaPhiMVector p1(mu_pt[0], mu_eta[0], mu_phi[0], mu_mass[0]);
96 ROOT::Math::PtEtaPhiMVector p2(mu_pt[1], mu_eta[1], mu_phi[1], mu_mass[1]);
97 ROOT::Math::PtEtaPhiMVector p3(el_pt[0], el_eta[0], el_phi[0], el_mass[0]);
98 ROOT::Math::PtEtaPhiMVector p4(el_pt[1], el_eta[1], el_phi[1], el_mass[1]);
99 auto mu_z = (p1 + p2).M();
100 auto el_z = (p3 + p4).M();
101 RVecF z_masses(2);
102 if (std::abs(mu_z - z_mass) < std::abs(el_z - z_mass)) {
103 z_masses[0] = mu_z;
104 z_masses[1] = el_z;
105 } else {
106 z_masses[0] = el_z;
107 z_masses[1] = mu_z;
108 }
109 return z_masses;
110}
111
112// Compute Higgs mass from two electrons and two muons
113float compute_higgs_mass_2el2mu(RVecF el_pt, RVecF el_eta, RVecF el_phi, RVecF el_mass, RVecF mu_pt, RVecF mu_eta,
114 RVecF mu_phi, RVecF mu_mass)
115{
116 ROOT::Math::PtEtaPhiMVector p1(mu_pt[0], mu_eta[0], mu_phi[0], mu_mass[0]);
117 ROOT::Math::PtEtaPhiMVector p2(mu_pt[1], mu_eta[1], mu_phi[1], mu_mass[1]);
118 ROOT::Math::PtEtaPhiMVector p3(el_pt[0], el_eta[0], el_phi[0], el_mass[0]);
119 ROOT::Math::PtEtaPhiMVector p4(el_pt[1], el_eta[1], el_phi[1], el_mass[1]);
120 return (p1 + p2 + p3 + p4).M();
121}
122
123bool filter_z_dr(const RVec<RVec<size_t>> &idx, RVecF eta, RVecF phi)
124{
125 for (size_t i = 0; i < 2; i++) {
126 const auto i1 = idx[i][0];
127 const auto i2 = idx[i][1];
128 const auto dr = DeltaR(eta[i1], eta[i2], phi[i1], phi[i2]);
129 if (dr < 0.02) {
130 return false;
131 }
132 }
133 return true;
134};
135
136bool pt_cuts(RVecF mu_pt, RVecF el_pt)
137{
138 auto mu_pt_sorted = Reverse(Sort(mu_pt));
139 if (mu_pt_sorted[0] > 20 && mu_pt_sorted[1] > 10) {
140 return true;
141 }
142 auto el_pt_sorted = Reverse(Sort(el_pt));
143 if (el_pt_sorted[0] > 20 && el_pt_sorted[1] > 10) {
144 return true;
145 }
146 return false;
147}
148
149bool dr_cuts(RVecF mu_eta, RVecF mu_phi, RVecF el_eta, RVecF el_phi)
150{
151 auto mu_dr = DeltaR(mu_eta[0], mu_eta[1], mu_phi[0], mu_phi[1]);
152 auto el_dr = DeltaR(el_eta[0], el_eta[1], el_phi[0], el_phi[1]);
153 if (mu_dr < 0.02 || el_dr < 0.02) {
154 return false;
155 }
156 return true;
157}
Class describing a generic LorentzVector in the 4D space-time, using the specified coordinate system ...
A "std::vector"-like collection of values implementing handy operation to analyse them.
Definition RVec.hxx:1492
TPaveText * pt
RVec< T > Reverse(const RVec< T > &v)
Return copy of reversed vector.
Definition RVec.hxx:2444
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:2569
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:112
RInterface<::ROOT::Detail::RDF::RNodeBase, void > RNode
This file contains a specialised ROOT message handler to test for diagnostic in unit tests.
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:431