rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

ATLAS_2019_I1750330

Semileptonic ttbar at 13 TeV
Experiment: ATLAS (LHC)
Inspire ID: 1750330
Status: VALIDATED
Authors:
  • Federica Fabbri
  • Francesco La Ruffa
  • Christian Gutschow
References: Beams: p+ p+
Beam energies: (6500.0, 6500.0) GeV
Run details:
  • pp -> non-allhadronic ttbar production at 13 TeV

Single- and double-differential cross-section measurements are presented for the production of top-quark pairs, in the lepton + jets channel at particle and parton level. Two topologies, resolved and boosted, are considered and the results are presented as a function of several kinematic variables characterising the top and the system and jet multiplicities. The study was performed using data from $pp$ collisions at centre-of-mass energy of 13 TeV collected in 2015 and 2016 by the ATLAS detector at the CERN Large Hadron Collider (LHC), corresponding to an integrated luminosity of 36 fb$^{-1}$. Due to the large $t\bar{t}$ cross-section at the LHC, such measurements allow a detailed study of the properties of top-quark production and decay, enabling precision tests of several Monte Carlo generators and fixed-order Standard Model predictions. Overall, there is good agreement between the theoretical predictions and the data.

Source code: ATLAS_2019_I1750330.cc
  1// -*- C++ -*
  2#include "Rivet/Analysis.hh"
  3#include "Rivet/Projections/FinalState.hh"
  4#include "Rivet/Projections/VetoedFinalState.hh"
  5#include "Rivet/Projections/InvisibleFinalState.hh"
  6#include "Rivet/Projections/PromptFinalState.hh"
  7#include "Rivet/Projections/LeptonFinder.hh"
  8#include "Rivet/Projections/MissingMomentum.hh"
  9#include "Rivet/Projections/FastJets.hh"
 10
 11namespace Rivet {
 12
 13
 14  /// @brief Semileptonic ttbar at 13 TeV
 15  class ATLAS_2019_I1750330 : public Analysis {
 16  public:
 17
 18    RIVET_DEFAULT_ANALYSIS_CTOR(ATLAS_2019_I1750330);
 19
 20    void init() {
 21
 22      _doBoosted = true, _doResolved = true;
 23      if ( getOption("TYPE") == "BOOSTED" )  _doResolved = false;
 24      else if ( getOption("TYPE") == "RESOLVED" )  _doBoosted = false;
 25
 26      Cut eta_full =  (Cuts::abseta < 5.0);
 27      Cut lep_cuts = (Cuts::abseta < 2.5) && (Cuts::pT > 27*GeV);
 28      const FinalState fs(eta_full);
 29
 30      FinalState all_photons(fs, Cuts::abspid == PID::PHOTON);
 31
 32      PromptFinalState photons(all_photons, TauDecaysAs::NONPROMPT);
 33      declare(photons, "photons");
 34
 35      PromptFinalState electrons(Cuts::abspid == PID::ELECTRON, TauDecaysAs::PROMPT);
 36      declare(electrons, "electrons");
 37
 38      LeptonFinder dressedelectrons(electrons, photons, 0.1, lep_cuts);
 39      declare(dressedelectrons, "dressedelectrons");
 40
 41      LeptonFinder ewdressedelectrons(electrons, all_photons, 0.1, eta_full);
 42      declare(ewdressedelectrons, "ewdressedelectrons");
 43
 44      PromptFinalState muons(Cuts::abspid == PID::MUON, TauDecaysAs::PROMPT);
 45      declare(muons, "muons");
 46
 47      LeptonFinder dressedmuons(muons, photons, 0.1, lep_cuts);
 48      declare(dressedmuons, "dressedmuons");
 49
 50      LeptonFinder ewdressedmuons(muons, all_photons, 0.1, eta_full);
 51      declare(ewdressedmuons, "ewdressedmuons");
 52
 53      InvisibleFinalState neutrinos(OnlyPrompt::YES, TauDecaysAs::PROMPT);
 54
 55      VetoedFinalState vfs(fs);
 56      vfs.addVetoOnThisFinalState(dressedelectrons);
 57      vfs.addVetoOnThisFinalState(dressedmuons);
 58      vfs.addVetoOnThisFinalState(neutrinos);
 59      FastJets jets(vfs, JetAlg::ANTIKT, 0.4, JetMuons::ALL, JetInvisibles::ALL);
 60      declare(jets, "boosted_jets");
 61
 62      VetoedFinalState vfs_res(fs);
 63      vfs_res.addVetoOnThisFinalState(ewdressedelectrons);
 64      vfs_res.addVetoOnThisFinalState(ewdressedmuons);
 65      vfs_res.addVetoOnThisFinalState(neutrinos);
 66      FastJets jets_res(vfs_res, JetAlg::ANTIKT, 0.4, JetMuons::ALL, JetInvisibles::ALL);
 67      declare(jets_res, "resolved_jets");
 68
 69      declare(MissingMomentum(), "MissingMomentum");
 70
 71      // Bins for 2D resolved
 72      std::vector<double> ttbar_m_2D_bins = {200,400,550,700,1000,2000};
 73      std::vector<double> top_had_pt_2D_bins = {0,60,120,200,300,1000};
 74      std::vector<double> ttbar_pt_2D_bins = {0,30,80,190,800};
 75      std::vector<double> top_had_abs_y_2D_bins = {0,0.7,1.4,2.5};
 76      std::vector<double> ttbar_abs_y_2D_bins = {0.0,0.4,0.8,1.2,2.5};
 77
 78      std::vector<double> n_jet_bins = {3.5, 4.5, 5.5, 6.5, 7.5};
 79      std::vector<double> n_jet_bins_for_ttbar_m = {3.5, 4.5, 5.5, 6.5};
 80      std::vector<double> n_extrajet_bins = {-0.5, 0.5, 1.5, 2.5, 3.5, 4.5, 5.5, 6.5};
 81
 82      //Bins for 2D boosted
 83      std::vector<double> eta_2D_bins = {0,1,2};
 84      std::vector<double> etattbar_2D_bins = {0,60,120,200,300,1000};
 85      std::vector<double> pttbar_2D_bins = {0.0, 40.0, 150.0, 1000.0};
 86      std::vector<double> mtt_2D_bins = {490.0, 1160, 3000.0};
 87      std::vector<double> eta_external_2D_bins = {0.0, 0.65, 1.3, 2.0};
 88      std::vector<double> ptt_external_mtt_2D_bins = {0.0, 40.0, 150.0, 1000.0};
 89      std::vector<double> Htt_external_2D_bins = {350.0, 780.0, 2500.0};
 90      std::vector<double> eta_external_ptt_2D_bins = {0.0, 0.65, 2.0};
 91
 92      std::vector<double> n_jet_pttop_bins = {-0.5,1.5,2.5,3.5};
 93      std::vector<double> n_jet_ptttbar_bins = {-0.5,1.5,3.5};
 94      std::vector<double> n_jet_Pout_bins = {-0.5,1.5,3.5};
 95      std::vector<double> n_jet_mtt_bins = {-0.5,0.5,1.5,2.5};
 96
 97      //Resolved histograms (digits correspond to "Table ID" from HepData)
 98      book2D("ttbar_m_top_had_pt_multi_norm", ttbar_m_2D_bins,54);
 99      book2D("ttbar_m_top_had_pt_multi", ttbar_m_2D_bins, 74);
100
101      book2D("ttbar_m_ttbar_pt_multi_norm", ttbar_m_2D_bins, 94);
102      book2D("ttbar_m_ttbar_pt_multi", ttbar_m_2D_bins, 114);
103
104      book2D("top_had_pt_absPout_multi_norm", top_had_pt_2D_bins, 134);
105      book2D("top_had_pt_absPout_multi", top_had_pt_2D_bins,154);
106
107      book2D("top_had_pt_jet_n_multi_norm", n_jet_bins,174);
108      book2D("top_had_pt_jet_n_multi", n_jet_bins, 188 );
109
110      book2D("ttbar_m_jet_n_multi_norm", n_jet_bins_for_ttbar_m,202);
111      book2D("ttbar_m_jet_n_multi", n_jet_bins_for_ttbar_m,211);
112
113      book2D("ttbar_pt_jet_n_multi_norm", n_jet_bins, 220 );
114      book2D("ttbar_pt_jet_n_multi", n_jet_bins, 234 );
115
116      book2D("absPout_jet_n_multi_norm", n_jet_bins, 248 );
117      book2D("absPout_jet_n_multi", n_jet_bins, 262 );
118
119      book2D("deltaPhi_tt_jet_n_multi_norm", n_jet_bins, 276 );
120      book2D("deltaPhi_tt_jet_n_multi", n_jet_bins, 290 );
121
122      book2D("HT_tt_jet_n_multi_norm", n_jet_bins, 304 );
123      book2D("HT_tt_jet_n_multi", n_jet_bins, 318 );
124
125      book2D("top_had_abs_y_jet_n_multi_norm", n_jet_bins, 332 );
126      book2D("top_had_abs_y_jet_n_multi", n_jet_bins, 346 );
127
128      book2D("ttbar_abs_y_jet_n_multi_norm", n_jet_bins, 360 );
129      book2D("ttbar_abs_y_jet_n_multi", n_jet_bins, 374 );
130
131      book2D("chi_tt_jet_n_multi_norm", n_jet_bins, 388 );
132      book2D("chi_tt_jet_n_multi", n_jet_bins, 402 );
133
134      book2D("top_had_abs_y_top_had_pt_multi_norm", top_had_abs_y_2D_bins,416 );
135      book2D("top_had_abs_y_top_had_pt_multi", top_had_abs_y_2D_bins, 425);
136
137      book2D("ttbar_abs_y_ttbar_pt_multi_norm", ttbar_abs_y_2D_bins, 434);
138      book2D("ttbar_abs_y_ttbar_pt_multi", ttbar_abs_y_2D_bins, 448);
139
140      book2D("ttbar_abs_y_ttbar_m_multi_norm", ttbar_abs_y_2D_bins, 462);
141      book2D("ttbar_abs_y_ttbar_m_multi", ttbar_abs_y_2D_bins, 476);
142
143      book2D("ttbar_pt_top_had_pt_multi_norm", ttbar_pt_2D_bins, 490);
144      book2D("ttbar_pt_top_had_pt_multi", ttbar_pt_2D_bins, 504);
145
146      book_hist("top_had_pt",1);
147      book_hist("top_had_abs_y_fine",5);
148      book_hist("leading_top_pt",9);
149      book_hist("subleading_top_pt",13);
150      book_hist("ttbar_m",17);
151      book_hist("ttbar_pt",21);
152      book_hist("absPout",25);
153      book_hist("deltaPhi_tt",29);
154      book_hist("HT_tt",33);
155      book_disc("extrajet_n",37);
156      book_hist("ttbar_abs_y_fine",41);
157      book_hist("abs_y_boost",45);
158      book_hist("chi_tt",49);
159
160      //Boosted histograms (digits correspond to "Table ID" from HepData)
161      book2D("boosted_rc_pttop_etatop_multi", eta_2D_bins, 922);
162      book2D("boosted_rc_pttop_etattbar_multi", eta_2D_bins, 912);
163      book2D("boosted_rc_pttop_ptttbar_multi", pttbar_2D_bins, 898);
164      book2D("boosted_rc_pttop_mttbar_multi", mtt_2D_bins, 932);
165      book2D("boosted_rc_mttbar_etattbar_multi", eta_external_2D_bins, 974);
166      book2D("boosted_rc_mttbar_ptttbar_multi", ptt_external_mtt_2D_bins, 956);
167      book2D("boosted_rc_mttbar_HT_multi", Htt_external_2D_bins, 942);
168      book2D("boosted_rc_pttop_extrajet_multi", n_jet_pttop_bins, 992);
169      book2D("boosted_rc_ptttbar_extrajet_multi", n_jet_ptttbar_bins, 1006);
170      book2D("boosted_rc_mttbar_extrajet_multi", n_jet_mtt_bins, 1020);
171
172      book2D("boosted_rc_pttop_etatop_multi_norm", eta_2D_bins, 917);
173      book2D("boosted_rc_pttop_etattbar_multi_norm", eta_2D_bins, 907);
174      book2D("boosted_rc_pttop_ptttbar_multi_norm", pttbar_2D_bins, 889);
175      book2D("boosted_rc_pttop_mttbar_multi_norm", mtt_2D_bins, 927);
176      book2D("boosted_rc_mttbar_etattbar_multi_norm", eta_external_2D_bins, 965);
177      book2D("boosted_rc_mttbar_ptttbar_multi_norm", ptt_external_mtt_2D_bins, 947);
178      book2D("boosted_rc_mttbar_HT_multi_norm", Htt_external_2D_bins, 937);
179      book2D("boosted_rc_pttop_extrajet_multi_norm", n_jet_pttop_bins, 983);
180      book2D("boosted_rc_ptttbar_extrajet_multi_norm", n_jet_ptttbar_bins, 1001);
181      book2D("boosted_rc_mttbar_extrajet_multi_norm", n_jet_mtt_bins, 1011);
182
183      book_hist("hadTop_boosted_rc_pt",840);
184      book_hist("hadTop_boosted_rc_y",844);
185      book_hist("LeadingTop_boosted_rc_pt",848);
186      book_hist("SubLeadingTop_boosted_rc_pt",852);
187      book_hist("boosted_rc_Pout_lep",872);
188      book_hist("boosted_rc_chi_tt",868);
189      book_hist("boosted_rc_HT",876);
190      book_disc("hadTop_boosted_rc_subjets",884);
191      book_disc("boosted_rc_extrajet",880);
192      book_hist("ttbar_boosted_rc_m",864);
193      book_hist("ttbar_boosted_rc_pt",856);
194      book_hist("ttbar_boosted_rc_Rapidity",860);
195
196    }
197
198
199    void analyze(const Event& event) {
200      if (_doResolved)  Resolved_selection(event);
201      if (_doBoosted)   Boosted_selection(event);
202    }
203
204
205    void Resolved_selection(const Event& event) {
206
207      // Get the selected objects, using the projections.
208      DressedLeptons electrons = apply<LeptonFinder>(event, "dressedelectrons").dressedLeptons();
209      DressedLeptons muons     = apply<LeptonFinder>(event, "dressedmuons").dressedLeptons();
210      const Jets& jets = apply<FastJets>(event, "resolved_jets").jetsByPt(Cuts::pT > 25*GeV && Cuts::abseta < 2.5);
211      FourMomentum met = apply<MissingMomentum>(event, "MissingMomentum").missingMomentum();
212
213      Jets bjets, lightjets;
214
215      // OVERLAP REMOVAL
216      idiscardIfAnyDeltaRLess(muons, jets, 0.4);
217      idiscardIfAnyDeltaRLess(electrons, jets, 0.4);
218
219      // b-tagging
220      // If there are more than 2 b-tagged jets, the extra b-tagged jets will be treat as light jets
221      for (const Jet& jet : jets) {
222        bool b_tagged = jet.bTagged(Cuts::pT > 5*GeV);
223        if ( b_tagged && bjets.size() < 2)  bjets +=jet;
224        else lightjets += jet;
225      }
226
227      bool single_electron = electrons.size() == 1 && muons.empty();
228      bool single_muon     = muons.size() == 1 && electrons.empty();
229
230      DressedLepton *lepton = NULL;
231      if (single_electron)   lepton = &electrons[0];
232      else if (single_muon)  lepton = &muons[0];
233
234      if (!single_electron && !single_muon) vetoEvent;
235      bool num_b_tagged_jets = (bjets.size() == 2);
236      if (!num_b_tagged_jets) vetoEvent;
237
238      if (lightjets.size() < 2) vetoEvent;
239
240      FourMomentum pbjet1; //Momentum of bjet1
241      FourMomentum pbjet2; //Momentum of bjet
242      if ( deltaR(bjets[0], *lepton) <= deltaR(bjets[1], *lepton) ) {
243        pbjet1 = bjets[0].momentum();
244        pbjet2 = bjets[1].momentum();
245      } else {
246        pbjet1 = bjets[1].momentum();
247        pbjet2 = bjets[0].momentum();
248      }
249
250      double bestWmass = 1000.0*TeV;
251      double mWPDG = 80.399*GeV;
252      int Wj1index = -1, Wj2index = -1;
253      for (unsigned int i = 0; i < (lightjets.size() - 1); ++i) {
254        for (unsigned int j = i + 1; j < lightjets.size(); ++j) {
255          double wmass = (lightjets[i].momentum() + lightjets[j].momentum()).mass();
256          if (fabs(wmass - mWPDG) < fabs(bestWmass - mWPDG)) {
257            bestWmass = wmass;
258            Wj1index = i;
259            Wj2index = j;
260          }
261        }
262      }
263
264      FourMomentum pjet1 = lightjets[Wj1index].momentum();
265      FourMomentum pjet2 = lightjets[Wj2index].momentum();
266
267      // compute hadronic W boson
268      FourMomentum pWhadron = pjet1 + pjet2;
269      double pz = computeneutrinoz(lepton->momentum(), met);
270      FourMomentum ppseudoneutrino( sqrt(sqr(met.px()) + sqr(met.py()) + sqr(pz)), met.px(), met.py(), pz);
271
272      //compute leptonic, hadronic, combined pseudo-top
273      FourMomentum ppseudotoplepton = lepton->momentum() + ppseudoneutrino + pbjet1;
274      FourMomentum ppseudotophadron = pbjet2 + pWhadron;
275      FourMomentum pttbar = ppseudotoplepton + ppseudotophadron;
276
277      Vector3 z_versor(0,0,1);
278      Vector3 vpseudotophadron = ppseudotophadron.vector3();
279      Vector3 vpseudotoplepton = ppseudotoplepton.vector3();
280
281      // Variables
282      double ystar = (ppseudotophadron.pt() > ppseudotoplepton.pt()) ? 0.5 * (ppseudotophadron.rap()-ppseudotoplepton.rap()) : 0.5*(ppseudotoplepton.rap()-ppseudotophadron.rap());
283      double chi_ttbar = exp(2 * fabs(ystar));
284      double deltaPhi_ttbar = deltaPhi(ppseudotoplepton,ppseudotophadron);
285      double HT_ttbar = ppseudotophadron.pt() + ppseudotoplepton.pt();
286      double Yboost = 0.5 * (ppseudotophadron.rapidity() + ppseudotoplepton.rapidity());
287      double Pout = vpseudotophadron.dot((vpseudotoplepton.cross(z_versor))/(vpseudotoplepton.cross(z_versor).mod()));
288      double absPout = fabs(Pout);
289      double Leading_top_pt = (ppseudotophadron.pt() > ppseudotoplepton.pt()) ? ppseudotophadron.pt() : ppseudotoplepton.pt();
290      double Subleading_top_pt = (ppseudotophadron.pt() > ppseudotoplepton.pt()) ? ppseudotoplepton.pt() : ppseudotophadron.pt();
291      int jet_multiplicity = jets.size();
292      int extrajet_n = jet_multiplicity - 4;
293      int new_jet_multi = TransformJetMultiplicity(jet_multiplicity);
294      int new_jet_multi_for_ttbar_m = TransformJetMultiplicity_for_ttbar_m(jet_multiplicity);
295      const string new_extrajet_multi = TransformExtrajetMultiplicity(extrajet_n);
296
297      _h_multi["top_had_pt_absPout_multi"]->fill(ppseudotophadron.pt()/GeV, absPout);
298      _h_multi["ttbar_m_top_had_pt_multi"]->fill(pttbar.mass()/GeV, ppseudotophadron.pt()/GeV);
299      _h_multi["ttbar_m_ttbar_pt_multi"]->fill(pttbar.mass()/GeV, pttbar.pt()/GeV);
300      _h_multi["ttbar_pt_top_had_pt_multi"]->fill(pttbar.pt()/GeV, ppseudotophadron.pt()/GeV);
301      _h_multi["ttbar_abs_y_ttbar_pt_multi"]->fill(pttbar.absrap(), pttbar.pt()/GeV);
302      _h_multi["ttbar_abs_y_ttbar_m_multi"]->fill(pttbar.absrap(), pttbar.mass()/GeV);
303      _h_multi["top_had_abs_y_top_had_pt_multi"]->fill(ppseudotophadron.absrap(), ppseudotophadron.pt()/GeV);
304
305      _h_multi["ttbar_pt_jet_n_multi"]->fill(new_jet_multi, pttbar.pt()/GeV);
306      _h_multi["ttbar_m_jet_n_multi"]->fill(new_jet_multi_for_ttbar_m, pttbar.mass()/GeV);
307      _h_multi["chi_tt_jet_n_multi"]->fill(new_jet_multi, chi_ttbar);
308      _h_multi["absPout_jet_n_multi"]->fill(new_jet_multi, absPout);
309      _h_multi["deltaPhi_tt_jet_n_multi"]->fill(new_jet_multi, deltaPhi_ttbar);
310      _h_multi["HT_tt_jet_n_multi"]->fill(new_jet_multi, HT_ttbar/GeV);
311      _h_multi["top_had_pt_jet_n_multi"]->fill(new_jet_multi, ppseudotophadron.pt()/GeV);
312      _h_multi["top_had_abs_y_jet_n_multi"]->fill(new_jet_multi, ppseudotophadron.absrap());
313      _h_multi["ttbar_abs_y_jet_n_multi"]->fill(new_jet_multi, pttbar.absrap());
314
315      _h_multi["top_had_pt_absPout_multi_norm"]->fill(ppseudotophadron.pt()/GeV, absPout);
316      _h_multi["ttbar_m_top_had_pt_multi_norm"]->fill(pttbar.mass()/GeV, ppseudotophadron.pt()/GeV);
317      _h_multi["ttbar_m_ttbar_pt_multi_norm"]->fill(pttbar.mass()/GeV, pttbar.pt()/GeV);
318      _h_multi["ttbar_pt_top_had_pt_multi_norm"]->fill(pttbar.pt()/GeV, ppseudotophadron.pt()/GeV);
319      _h_multi["ttbar_abs_y_ttbar_pt_multi_norm"]->fill(pttbar.absrap(), pttbar.pt()/GeV);
320      _h_multi["ttbar_abs_y_ttbar_m_multi_norm"]->fill(pttbar.absrap(), pttbar.mass()/GeV);
321      _h_multi["top_had_abs_y_top_had_pt_multi_norm"]->fill(ppseudotophadron.absrap(), ppseudotophadron.pt()/GeV);
322
323      _h_multi["ttbar_pt_jet_n_multi_norm"]->fill(new_jet_multi, pttbar.pt());
324      _h_multi["ttbar_m_jet_n_multi_norm"]->fill(new_jet_multi_for_ttbar_m, pttbar.mass());
325      _h_multi["chi_tt_jet_n_multi_norm"]->fill(new_jet_multi, chi_ttbar);
326      _h_multi["absPout_jet_n_multi_norm"]->fill(new_jet_multi, absPout);
327      _h_multi["deltaPhi_tt_jet_n_multi_norm"]->fill(new_jet_multi, deltaPhi_ttbar);
328      _h_multi["HT_tt_jet_n_multi_norm"]->fill(new_jet_multi, HT_ttbar/GeV);
329      _h_multi["top_had_pt_jet_n_multi_norm"]->fill(new_jet_multi, ppseudotophadron.pt()/GeV);
330      _h_multi["top_had_abs_y_jet_n_multi_norm"]->fill(new_jet_multi, ppseudotophadron.absrap());
331      _h_multi["ttbar_abs_y_jet_n_multi_norm"]->fill(new_jet_multi, pttbar.absrap());
332
333      _h["chi_tt"]->fill(chi_ttbar);
334      _h["deltaPhi_tt"]->fill(deltaPhi_ttbar);
335      _h["HT_tt"]->fill(HT_ttbar/GeV);
336      _h["absPout"]->fill(absPout);
337      _h["abs_y_boost"]->fill(fabs(Yboost));
338      _h["top_had_pt"]->fill(ppseudotophadron.pt()/GeV);
339      _h["top_had_abs_y_fine"]->fill(ppseudotophadron.absrap());
340      _h["ttbar_pt"]->fill(pttbar.pt()/GeV);
341      _h["ttbar_m"]->fill(pttbar.mass()/GeV);
342      _h["ttbar_abs_y_fine"]->fill(pttbar.absrap());
343      _h["leading_top_pt"]->fill(Leading_top_pt/GeV);
344      _h["subleading_top_pt"]->fill(Subleading_top_pt/GeV);
345      _d["extrajet_n"]->fill(new_extrajet_multi);
346
347      _h["chi_tt_norm"]->fill(chi_ttbar);
348      _h["deltaPhi_tt_norm"]->fill(deltaPhi_ttbar);
349      _h["HT_tt_norm"]->fill(HT_ttbar/GeV);
350      _h["absPout_norm"]->fill(absPout);
351      _h["abs_y_boost_norm"]->fill(fabs(Yboost));
352      _h["top_had_pt_norm"]->fill(ppseudotophadron.pt()/GeV);
353      _h["top_had_abs_y_fine_norm"]->fill(ppseudotophadron.absrap());
354      _h["ttbar_pt_norm"]->fill(pttbar.pt()/GeV);
355      _h["ttbar_m_norm"]->fill(pttbar.mass()/GeV);
356      _h["ttbar_abs_y_fine_norm"]->fill(pttbar.absrap());
357      _h["leading_top_pt_norm"]->fill(Leading_top_pt/GeV);
358      _h["subleading_top_pt_norm"]->fill(Subleading_top_pt/GeV);
359      _d["extrajet_n_norm"]->fill(new_extrajet_multi);
360    }
361
362    void Boosted_selection(const Event& event) {
363
364      //Projections
365      DressedLeptons electrons = apply<LeptonFinder>(event, "dressedelectrons").dressedLeptons();
366      DressedLeptons muons = apply<LeptonFinder>(event, "dressedmuons").dressedLeptons();
367      const Jets& jets = apply<FastJets>(event, "boosted_jets").jetsByPt(Cuts::pT > 25*GeV && Cuts::abseta <= 2.5);
368      const FourMomentum& met = apply<MissingMomentum>(event, "MissingMomentum").missingMomentum();
369
370      if (jets.size() < 2)  vetoEvent;
371      PseudoJets smallRjets;
372      for (const Jet& jet : jets) {
373        smallRjets += jet.pseudojet();
374        bool b_tagged = jet.bTagged(Cuts::pT > 5*GeV);
375        smallRjets[smallRjets.size()-1].set_user_index(b_tagged); // cheeky, but works
376      }
377
378      idiscardIfAnyDeltaRLess(muons, jets, 0.4);
379      idiscardIfAnyDeltaRLess(electrons, jets, 0.4);
380
381      fastjet::Filter trimmer(fastjet::JetDefinition(fastjet::antikt_algorithm, 1.0), fastjet::SelectorPtFractionMin(0.05));
382      fastjet::ClusterSequence antikt_cs(smallRjets, fastjet::JetDefinition(fastjet::antikt_algorithm, 1.0));
383      PseudoJets reclustered_jets = antikt_cs.inclusive_jets();
384
385      // trim the jets
386      Jets TrimmedJets;
387      for (const PseudoJet& pjet : reclustered_jets) {
388        PseudoJet ptrim = trimmer(pjet);
389        if (ptrim.perp() < 350*GeV)  continue;
390        if (fabs(ptrim.eta()) > 2.0) continue;
391        bool bTagged  = false;
392        Particles constituents;
393        for (const PseudoJet& c : ptrim.constituents()) {
394          // we only care about the number of subjets, so
395          // fine to treat as Particles with dummy PID
396          constituents += Particle(0, momentum(c));
397          bTagged |= c.user_index();
398        }
399        ptrim.set_user_index(bTagged);
400        TrimmedJets += Jet(ptrim, constituents);
401      }
402      Cut trim_selection = Cuts::abseta < 2.0 && Cuts::pT > 200*GeV && Cuts::massIn(120*GeV, 220*GeV);
403      iselect(isortByPt(TrimmedJets), trim_selection);
404      if (TrimmedJets.empty())  vetoEvent;
405
406
407      // SINGLE LEPTON
408      bool single_electron=(electrons.size() == 1) && (muons.empty());
409      bool single_muon=(muons.size() == 1) && (electrons.empty());
410
411      DressedLepton *lepton = NULL;
412      if (single_electron)   lepton = &electrons[0];
413      else if (single_muon)  lepton = &muons[0];
414      if (!single_electron && !single_muon) vetoEvent;
415
416      //MET
417      if (met.pT() < 20*GeV)  vetoEvent;
418
419      //MET+MWT
420      double transmass = TransMass(lepton->pt(), lepton->phi(), met.pt(), met.phi());
421      if ((met.pT() + transmass) < 60*GeV)   vetoEvent;
422
423      size_t subjets = 0;
424      bool btag_hadside=false;
425      bool hasHadTopCandidate = false;
426      FourMomentum HadTopCandidate;
427      for (const Jet& rc_jet : TrimmedJets) {
428        FourMomentum rc_jet_mom = rc_jet.mom();
429        if (rc_jet_mom.pt() < 350*GeV)   continue;
430        double dPhi_lepJet = fabs(deltaPhi(rc_jet_mom.phi(), lepton->phi()));
431        if (dPhi_lepJet < 1.)  continue;
432        if (rc_jet.pseudojet().user_index()) {
433          btag_hadside=true;
434        }
435        HadTopCandidate = momentum(rc_jet);
436        subjets = rc_jet.constituents().size();
437        hasHadTopCandidate = true;
438        break;
439      }
440      if (!hasHadTopCandidate)  vetoEvent;
441
442      Jets LepTopCandidates = discard(jets, [&](const Jet& j) {
443          return deltaR(j, HadTopCandidate) < 1.5 || deltaR(j, *lepton) > 2.0;
444        });
445      if (LepTopCandidates.empty()) vetoEvent;
446
447      FourMomentum ltop;
448      bool btag_lepside=false;
449      for (const Jet& jet : LepTopCandidates) {
450        if (jet.bTagged(Cuts::pT > 5*GeV)) {
451          btag_lepside = true;
452          ltop = jet.mom();
453          break;
454        }
455      }
456      if (!btag_hadside && !btag_lepside)  vetoEvent;
457      if (!btag_lepside)   ltop = LepTopCandidates[0].momentum();
458      double pz = computeneutrinoz(lepton->momentum(), met);
459      FourMomentum neutrino( sqrt(sqr(met.px()) + sqr(met.py()) + sqr(pz)), met.px(), met.py(), pz);
460      FourMomentum LeptonicTop = lepton->momentum() + neutrino + ltop;
461      FourMomentum HadronicTop = HadTopCandidate;
462      FourMomentum pttbar = HadronicTop + LeptonicTop;
463
464      Vector3 z_versor(0,0,1);
465      Vector3 vpseudotophadron = HadronicTop.vector3();
466      Vector3 vpseudotoplepton = LeptonicTop.vector3();
467      // Variables
468      double ystar = (HadronicTop.pt() > LeptonicTop.pt()) ? 0.5 * (HadronicTop.rap()-LeptonicTop.rap()) : 0.5*(LeptonicTop.rap()-HadronicTop.rap());
469      double chi_ttbar = exp(2 * fabs(ystar));
470      double pt_leading = (HadronicTop.pt() > LeptonicTop.pt()) ? HadronicTop.pt() : LeptonicTop.pt();
471      double pt_subleading = (HadronicTop.pt() > LeptonicTop.pt()) ? LeptonicTop.pt() : HadronicTop.pt();
472      double HT_ttbar = HadronicTop.pt() + LeptonicTop.pt();
473      double absPout_lep = fabs(vpseudotoplepton.dot((vpseudotophadron.cross(z_versor))/(vpseudotophadron.cross(z_versor).mod())));
474      size_t extrajet = smallRjets.size() - subjets -1;
475
476      const string new_subjets_multi = TransformExtrajetMultiplicity_boosted(subjets);
477      const string new_extrajet_multi = TransformExtrajetMultiplicity_boosted(extrajet);
478      size_t new_extrajet_multi_pttop = TransformJetMultiplicity_pttop(extrajet);
479      size_t new_extrajet_multi_ptttbar = TransformJetMultiplicity_ptttbar(extrajet);
480      size_t new_extrajet_multi_mttbar = TransformJetMultiplicity_mttbar(extrajet);
481
482      _h_multi["boosted_rc_pttop_etatop_multi"]->fill(HadronicTop.absrap(), HadronicTop.pt()/GeV);
483      _h_multi["boosted_rc_pttop_etattbar_multi"]->fill(pttbar.absrap(), HadronicTop.pt()/GeV);
484      _h_multi["boosted_rc_pttop_ptttbar_multi"]->fill(pttbar.pt()/GeV, HadronicTop.pt()/GeV);
485      _h_multi["boosted_rc_pttop_mttbar_multi"]->fill(pttbar.mass()/GeV, HadronicTop.pt()/GeV);
486      _h_multi["boosted_rc_mttbar_etattbar_multi"]->fill(pttbar.absrap(), pttbar.mass()/GeV);
487      _h_multi["boosted_rc_mttbar_ptttbar_multi"]->fill(pttbar.pt()/GeV, pttbar.mass()/GeV);
488      _h_multi["boosted_rc_mttbar_HT_multi"]->fill(HT_ttbar, pttbar.mass()/GeV);
489
490      _h_multi["boosted_rc_pttop_extrajet_multi"]->fill(new_extrajet_multi_pttop, HadronicTop.pt()/GeV);
491      _h_multi["boosted_rc_ptttbar_extrajet_multi"]->fill(new_extrajet_multi_ptttbar, pttbar.pt()/GeV);
492      _h_multi["boosted_rc_mttbar_extrajet_multi"]->fill(new_extrajet_multi_mttbar, pttbar.mass()/GeV);
493
494      _h_multi["boosted_rc_pttop_etatop_multi_norm"]->fill(HadronicTop.absrap(), HadronicTop.pt()/GeV);
495      _h_multi["boosted_rc_pttop_etattbar_multi_norm"]->fill(pttbar.absrap(), HadronicTop.pt()/GeV);
496      _h_multi["boosted_rc_pttop_ptttbar_multi_norm"]->fill(pttbar.pt()/GeV, HadronicTop.pt()/GeV);
497      _h_multi["boosted_rc_pttop_mttbar_multi_norm"]->fill(pttbar.mass()/GeV, HadronicTop.pt()/GeV);
498      _h_multi["boosted_rc_mttbar_etattbar_multi_norm"]->fill(pttbar.absrap(), pttbar.mass()/GeV);
499      _h_multi["boosted_rc_mttbar_ptttbar_multi_norm"]->fill(pttbar.pt()/GeV, pttbar.mass()/GeV);
500      _h_multi["boosted_rc_mttbar_HT_multi_norm"]->fill(HT_ttbar/GeV, pttbar.mass()/GeV);
501      _h_multi["boosted_rc_pttop_extrajet_multi_norm"]->fill(new_extrajet_multi_pttop, HadronicTop.pt()/GeV);
502      _h_multi["boosted_rc_ptttbar_extrajet_multi_norm"]->fill(new_extrajet_multi_ptttbar, pttbar.pt()/GeV);
503      _h_multi["boosted_rc_mttbar_extrajet_multi_norm"]->fill(new_extrajet_multi_mttbar, pttbar.mass()/GeV);
504
505      _h["hadTop_boosted_rc_pt"]->fill(HadronicTop.pt()/GeV);
506      _h["hadTop_boosted_rc_y"]->fill(HadronicTop.absrap());
507      _h["LeadingTop_boosted_rc_pt"]->fill(pt_leading/GeV);
508      _h["SubLeadingTop_boosted_rc_pt"]->fill(pt_subleading/GeV);
509      _h["boosted_rc_Pout_lep"]->fill(absPout_lep);
510      _h["boosted_rc_chi_tt"]->fill(chi_ttbar);
511      _h["boosted_rc_HT"]->fill(HT_ttbar/GeV);
512      _d["hadTop_boosted_rc_subjets"]->fill(new_subjets_multi);
513      _d["boosted_rc_extrajet"]->fill(new_extrajet_multi);
514      _h["ttbar_boosted_rc_m"]->fill(pttbar.mass()/GeV);
515      _h["ttbar_boosted_rc_pt"]->fill(pttbar.pt()/GeV);
516      _h["ttbar_boosted_rc_Rapidity"]->fill(pttbar.absrapidity());
517
518      _h["hadTop_boosted_rc_pt_norm"]->fill(HadronicTop.pt()/GeV);
519      _h["hadTop_boosted_rc_y_norm"]->fill(HadronicTop.absrap());
520      _h["LeadingTop_boosted_rc_pt_norm"]->fill(pt_leading/GeV);
521      _h["SubLeadingTop_boosted_rc_pt_norm"]->fill(pt_subleading/GeV);
522      _h["boosted_rc_Pout_lep_norm"]->fill(absPout_lep);
523      _h["boosted_rc_chi_tt_norm"]->fill(chi_ttbar);
524      _h["boosted_rc_HT_norm"]->fill(HT_ttbar/GeV);
525      _d["hadTop_boosted_rc_subjets_norm"]->fill(new_subjets_multi);
526      _d["boosted_rc_extrajet_norm"]->fill(new_extrajet_multi);
527      _h["ttbar_boosted_rc_m_norm"]->fill(pttbar.mass()/GeV);
528      _h["ttbar_boosted_rc_pt_norm"]->fill(pttbar.pt()/GeV);
529      _h["ttbar_boosted_rc_Rapidity_norm"]->fill(pttbar.absrap());
530    }
531
532
533    void finalize() {
534      // Normalize to cross-section
535      const double sf = crossSection()/picobarn / sumOfWeights();
536      for (auto& hit : _h) {
537        if (hit.first.find("_norm") != string::npos)  normalize(hit.second, 1.0, false);
538        else                                          scale(hit.second, sf);
539      }
540      for (auto& hit : _d) {
541        scale(hit.second, sf);
542        if (hit.first.find("_norm") != string::npos)  normalize(hit.second, 1.0, false);
543      }
544      for (auto& hit : _h_multi) {
545        if (hit.first.find("_norm") != string::npos) {
546          normalizeGroup(hit.second, 1.0, false);
547        }
548        else {
549          scale(hit.second, sf);
550        }
551      }
552      divByGroupWidth(_h_multi);
553    }
554
555
556  private:
557
558    bool _doBoosted, _doResolved;
559
560
561    double TransMass(double ptLep, double phiLep, double met, double phiMet) {
562      return std::sqrt(2.0*ptLep*met*( 1 - std::cos( phiLep-phiMet ) ) );
563    }
564
565
566    double computeneutrinoz(const FourMomentum& lepton, const FourMomentum& met) const {
567      //computing z component of neutrino momentum given lepton and met
568      double pzneutrino;
569      double m_W = 80.399; // in GeV, given in the paper
570      double k = (( sqr( m_W ) - sqr( lepton.mass() ) ) / 2 ) + (lepton.px() * met.px() + lepton.py() * met.py());
571      double a = sqr ( lepton.E() )- sqr ( lepton.pz() );
572      double b = -2*k*lepton.pz();
573      double c = sqr( lepton.E() ) * sqr( met.pT() ) - sqr( k );
574      double discriminant = sqr(b) - 4 * a * c;
575      double quad[2] = { (- b - sqrt(discriminant)) / (2 * a), (- b + sqrt(discriminant)) / (2 * a) }; //two possible quadratic solns
576      if (discriminant < 0)  pzneutrino = - b / (2 * a); //if the discriminant is negative
577      else { //if the discriminant is greater than or equal to zero, take the soln with smallest absolute value
578        double absquad[2];
579        for (int n=0; n<2; ++n)  absquad[n] = fabs(quad[n]);
580        if (absquad[0] < absquad[1])  pzneutrino = quad[0];
581        else                          pzneutrino = quad[1];
582      }
583      return pzneutrino;
584    }
585
586
587
588
589    void book2D(const string& name, const std::vector<double>& doubleDiff_bins, size_t table){
590      book(_h_multi[name], doubleDiff_bins);
591      for (auto& b : _h_multi[name]->bins()) {
592        book(b, table+b.index(), 1, 1);
593      }
594    }
595
596
597    void book_hist(const string& name, size_t table) {
598      // HepData entry has dummy "Table of Contents",
599      // so need to offset everything by one unit
600      book(_h[name], table+3, 1, 1);
601      book(_h[name+"_norm"], table+1, 1, 1);
602    }
603
604    void book_disc(const string& name, size_t table) {
605      // HepData entry has dummy "Table of Contents",
606      // so need to offset everything by one unit
607      book(_d[name], table+3, 1, 1);
608      book(_d[name+"_norm"], table+1, 1, 1);
609    }
610
611    size_t TransformJetMultiplicity(size_t jet_n) const { return jet_n > 7 ? 7 : jet_n; }
612
613    string TransformExtrajetMultiplicity(size_t jet_n) const {
614      if (jet_n == 0)       return "0.0"s;
615      else if (jet_n == 1)  return "1.0"s;
616      else if (jet_n == 2)  return "2.0"s;
617      else if (jet_n == 3)  return "3.0"s;
618      else if (jet_n == 4)  return "4.0"s;
619      else if (jet_n == 5)  return "5.0"s;
620      else                  return "$\\geq$6.0"s;
621    }
622
623    string TransformExtrajetMultiplicity_boosted(size_t jet_n) const {
624      if (jet_n == 0)       return "0.0"s;
625      else if (jet_n == 1)  return "1.0"s;
626      else if (jet_n == 2)  return "2.0"s;
627      else if (jet_n == 3)  return "3.0"s;
628      else                  return "$\\geq$4.0"s;
629    }
630
631    size_t TransformJetMultiplicity_for_ttbar_m(size_t jet_n) const { return jet_n > 6 ? 6 : jet_n; }
632
633    size_t TransformJetMultiplicity_pttop(size_t jet_n) const {
634      if (jet_n < 2)  return 0;
635      if (jet_n == 2)               return 2;
636      if (jet_n > 2)                return 3;
637      return jet_n;
638    }
639
640    size_t TransformJetMultiplicity_ptttbar(size_t jet_n) const {
641      if (jet_n < 2)  return 0;
642      if (jet_n >= 2)               return 2;
643      return jet_n;
644    }
645
646    size_t TransformJetMultiplicity_mttbar(size_t jet_n) const {
647      if (jet_n == 0)  return 0;
648      if (jet_n == 1)  return 1;
649      if (jet_n >= 2)  return 2;
650      return jet_n;
651    }
652
653    /// @name Objects that are used by the event selection decisions
654    /// @{
655    map<string, Histo1DPtr> _h;
656    map<string, BinnedHistoPtr<string>> _d;
657    map<string, Histo1DGroupPtr> _h_multi;
658    /// @}
659
660  };
661
662
663  RIVET_DECLARE_PLUGIN(ATLAS_2019_I1750330);
664
665}