Rivet analyses referenceATLAS_2019_I1718132Control region measurements for leptoquark search at 13 TeVExperiment: ATLAS (LHC) Inspire ID: 1718132 Status: VALIDATED Authors:
Beam energies: (6500.0, 6500.0) GeV Run details:
Searches for scalar leptoquarks pair-produced in proton-proton collisions at $\sqrt{s} = 13$ TeV at the Large Hadron Collider are performed by the ATLAS experiment. A data set corresponding to an integrated luminosity of 36.1 fb${}^{-1}$ is used. Final states containing two electrons or two muons and two or more jets are studied, as are states with one electron or muon, missing transverse momentum and two or more jets. No statistically significant excess above the Standard Model expectation is observed. The observed and expected lower limits on the leptoquark mass at 95% confidence level extend up to 1.25 TeV for first- and second-generation leptoquarks, as postulated in the minimal Buchmuller-Ruckl-Wyler model, assuming a branching ratio into a charged lepton and a quark of 50%. In addition, measurements of particle-level fiducial and differential cross sections are presented for the $Z\to ee$, $Z\to\mu\mu$ and $t\bar{t}$ processes in several regions related to the search control regions. Predictions from a range of generators are compared with the measurements, and good agreement is seen for many of the observables. However, the predictions for the $Z\to\ell\ell$ measurements in observables sensitive to jet energies disagree with the data. Source code: ATLAS_2019_I1718132.cc 1// -*- C++ -*-
2#include "Rivet/Analysis.hh"
3#include "Rivet/Projections/FastJets.hh"
4#include "Rivet/Projections/FinalState.hh"
5#include "Rivet/Projections/PromptFinalState.hh"
6#include "Rivet/Projections/LeptonFinder.hh"
7
8namespace Rivet {
9
10
11 /// @brief leptoquark search at 13 TeV
12 ///
13 /// @note This base class contains a "mode" variable to specify lepton channel
14 class ATLAS_2019_I1718132 : public Analysis {
15 public:
16
17 /// Constructor
18 RIVET_DEFAULT_ANALYSIS_CTOR(ATLAS_2019_I1718132);
19
20 /// Book histograms and initialise projections before the run
21 void init() {
22
23 // default to widest cut, electrons and muons.
24 _mode = 3;
25 if ( getOption("LMODE") == "ELEL" ) _mode = 1;
26 if ( getOption("LMODE") == "MUMU" ) _mode = 2;
27 if ( getOption("LMODE") == "ELMU" ) _mode = 3;
28
29 // Lepton cuts
30 Cut baseline_lep_cuts = Cuts::abseta < 2.5 && Cuts::pT >= 40*GeV;
31
32 // All final state particles
33 const FinalState fs;
34
35 // Get photons to dress leptons
36 FinalState photons(Cuts::abspid == PID::PHOTON);
37
38 // Find and dress the electrons and muons
39 PromptFinalState bare_leps(Cuts::abspid == PID::ELECTRON || Cuts::abspid == PID::MUON);
40 LeptonFinder dressed_leps(bare_leps, photons, 0.1, baseline_lep_cuts);
41 declare(dressed_leps, "leptons");
42
43 //and finally the jets:
44 FastJets jets(fs, JetAlg::ANTIKT, 0.4, JetMuons::NONE);
45 declare(jets, "jets");
46
47 size_t offset = _mode;
48 book(_h["JetPt_leading"], 1 + offset, 1, 1);
49 book(_h["JetPt_subleading"], 7 + offset, 1, 1);
50 book(_h["minDeltaPhiJ0_L"], 13 + offset, 1, 1);
51 book(_h["minDeltaPhiJ1_L"], 19 + offset, 1, 1);
52 book(_h["DeltaEtaJJ"], 25 + offset, 1, 1);
53 book(_h["DeltaPhiJJ"], 31 + offset, 1, 1);
54 book(_h["DeltaPhiLL"], 37 + offset, 1, 1);
55 book(_h["DiJetMass"], 43 + offset, 1, 1);
56 book(_h["DiLepPt"], 49 + offset, 1, 1);
57 book(_h["Ht"], 55 + offset, 1, 1);
58 book(_h["St"], 61 + offset, 1, 1);
59
60 offset = _mode + 3;
61 book(_hST["JetPt_leading"], 1 + offset, 1, 1);
62 book(_hST["JetPt_subleading"], 7 + offset, 1, 1);
63 book(_hST["minDeltaPhiJ0_L"], 13 + offset, 1, 1);
64 book(_hST["minDeltaPhiJ1_L"], 19 + offset, 1, 1);
65 book(_hST["DeltaEtaJJ"], 25 + offset, 1, 1);
66 book(_hST["DeltaPhiJJ"], 31 + offset, 1, 1);
67 book(_hST["DeltaPhiLL"], 37 + offset, 1, 1);
68 book(_hST["DiJetMass"], 43 + offset, 1, 1);
69 book(_hST["DiLepPt"], 49 + offset, 1, 1);
70 book(_hST["Ht"], 55 + offset, 1, 1);
71 book(_hST["St"], 61 + offset, 1, 1);
72
73 }
74
75 /// Perform the per-event analysis
76 void analyze(const Event& event) {
77
78 // Get the selected leptons:
79 DressedLeptons leptons = apply<LeptonFinder>(event, "leptons").dressedLeptons();
80
81 // get the selected jets:
82 Jets jets = apply<JetFinder>(event, "jets").jetsByPt(Cuts::pT > 60*GeV && Cuts::absrap < 2.5);
83
84 // exclude jets which are actually electrons
85
86
87 // This would be super-sweet, but unfortunately lxplus default gcc4.8 doesn't like it :(
88 /*idiscardIfAny(jets, leptons, [](const Jet& jet, const DressedLepton& lep) {
89 return lep.abspid() == PID::ELECTRON and deltaR(jet, lep) < 0.2;
90 });*/
91
92 for (const DressedLepton& lep : leptons) {
93 idiscard(jets, [&](const Jet& jet) {
94 return lep.abspid() == PID::ELECTRON and deltaR(jet, lep) < 0.2;
95 });
96 }
97
98 // remove cases where muons are too close to a jet
99 if (_mode == 3) { // el-mu case
100 /*idiscardIfAny(leptons, jets, [](const DressedLepton& lep, const Jet& jet) {
101 return lep.abspid() == PID::MUON and deltaR(jet, lep) < 0.4;
102 });*/
103 for (const DressedLepton& lep : leptons) {
104 idiscard(jets, [&](const Jet& jet) {
105 return lep.abspid() == PID::MUON and deltaR(jet, lep) < 0.4;
106 });
107 }
108 }
109
110 // make sure we have the right number of jets
111 if (jets.size() < 2) vetoEvent;
112
113 // make sure we have right number of leptons
114 size_t requiredMuons = 0, requiredElecs = 0;
115 if (_mode == 1) requiredElecs = 2;
116 if (_mode == 2) requiredMuons = 2;
117 if (_mode == 3) requiredMuons = requiredElecs = 1;
118
119 size_t nEl = count(leptons, [](const DressedLepton& lep) { return lep.abspid() == PID::ELECTRON; });
120 if (nEl != requiredElecs) vetoEvent;
121
122 size_t nMu = count(leptons, [](const DressedLepton& lep) { return lep.abspid() == PID::MUON; });
123 if (nMu != requiredMuons) vetoEvent;
124
125 // make sure leptons are in right order (should be OK byt better safe than sorry!)
126 std::sort(leptons.begin(), leptons.end(), cmpMomByPt);
127
128 // calculate all observables and store in dict
129 const double jetpt_leading = jets[0].pT()/GeV;
130 const double jetpt_subleading = jets[1].pT()/GeV;
131 const double mll = (leptons[0].momentum() + leptons[1].momentum()).mass()/GeV;
132 const double st = (leptons[0].pT() + leptons[1].pT() + jets[0].pT() + jets[1].pT())/GeV;
133 const double ht = (jets[0].pT() + jets[1].pT())/GeV;
134 const double dijetmass = (jets[0].mom() + jets[1].mom()).mass()/GeV;
135 const double dileppt = (leptons[0].momentum() + leptons[1].momentum()).pT()/GeV;
136 const double deltaphijj = deltaPhi(jets[0], jets[1]);
137 const double deltaetajj = deltaEta(jets[0], jets[1]);
138 const double deltaphill = deltaPhi(leptons[0], leptons[1]);
139 const double mindeltaphij0_l = min(deltaPhi(jets[0], leptons[0]), deltaPhi(jets[0], leptons[1]));
140 const double mindeltaphij1_l = min(deltaPhi(jets[1], leptons[0]), deltaPhi(jets[1], leptons[1]));
141
142 // add Z-mass window cut if needed
143 bool addMllCut = _mode == 1 || _mode == 2;
144 if (addMllCut && (mll < 70. || 110. < mll)) vetoEvent;
145
146 // fill output histos
147 _h["JetPt_leading"]->fill(jetpt_leading);
148 _h["JetPt_subleading"]->fill(jetpt_subleading);
149 _h["St"]->fill(st);
150 _h["Ht"]->fill(ht);
151 _h["DiJetMass"]->fill(dijetmass);
152 _h["DiLepPt"]->fill(dileppt);
153 _h["DeltaPhiJJ"]->fill(deltaphijj);
154 _h["DeltaEtaJJ"]->fill(deltaetajj);
155 _h["DeltaPhiLL"]->fill(deltaphill);
156 _h["minDeltaPhiJ0_L"]->fill(mindeltaphij0_l);
157 _h["minDeltaPhiJ1_L"]->fill(mindeltaphij1_l);
158
159 // "extreme" ST cut
160 if (st > 600.) {
161 // fill output histos
162 _hST["JetPt_leading"]->fill(jetpt_leading);
163 _hST["JetPt_subleading"]->fill(jetpt_subleading);
164 _hST["St"]->fill(st);
165 _hST["Ht"]->fill(ht);
166 _hST["DiJetMass"]->fill(dijetmass);
167 _hST["DiLepPt"]->fill(dileppt);
168 _hST["DeltaPhiJJ"]->fill(deltaphijj);
169 _hST["DeltaEtaJJ"]->fill(deltaetajj);
170 _hST["DeltaPhiLL"]->fill(deltaphill);
171 _hST["minDeltaPhiJ0_L"]->fill(mindeltaphij0_l);
172 _hST["minDeltaPhiJ1_L"]->fill(mindeltaphij1_l);
173 }
174
175 }
176
177
178 void finalize() {
179 const double sf = crossSectionPerEvent();
180 scale(_h, sf); scale(_hST, sf);
181 }
182
183
184 protected:
185
186 size_t _mode;
187
188 private:
189
190 map<string, Histo1DPtr> _h, _hST;
191
192 };
193
194
195 RIVET_DECLARE_PLUGIN(ATLAS_2019_I1718132);
196
197}
|