rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

CMS_2017_I1518399

Differential $t\bar{t}$ production cross-section as a function of the leading jet mass for boosted top quarks at 8 TeV
Experiment: CMS (LHC)
Inspire ID: 1518399
Status: VALIDATED
Authors:
  • Torben Dreyer
  • Roman Kogler
  • Johannes Haller
References:
  • TOP-15-015
  • arXiv: 1703.06330
  • Submitted to Eur. Phys. J. C.
Beams: p+ p+
Beam energies: (4000.0, 4000.0) GeV
Run details:
  • $t\bar{t}$ events at $\sqrt{s}=8 \text{TeV}$. Boosted topology requires high statistics in a single run or custom merging of the YODA files.

Measurement of the differential and normalised differential $t\bar{t}$ production cross section as a function of the leading jet mass. The measurement was performed in the lepton+jets channel and the fiducial measurement phase space is enriched with events where the leading jet includes all decay products of a hadronically decaying top quark. This analysis is to be run on ${t\bar{t}}$ simulation. The objects are defined as follows - Lepton: electron or muon originating from the decay of a W boson - Jets: jets are clustered from final state particles using the Cambridge-Aachen algorithm with a distance parameter of 1.2. If a lepton overlaps with a jet (distance in the eta-phi plane smaller 1.2) its four momentum is subtracted from the four momentum of the jet. Definition of the particle level phase space: - only $t\bar{t}$ decays with one top quark decaying hadronically and one top quark decaying leptonically including an electron or muon from the $W$ decay. - lepton $\pt > 45 \text{GeV}$ and $|\eta| < 2.1$ - at least one jet with $\pt > 400 \text{GeV}$ and $|\eta| < 2.5$ - a second jet with $\pt > 150 \text{GeV}$ and $|\eta| < 2.5$ - veto on additional jets with $\pt > 150 \text{GeV}$ and $|\eta| < 2.5$ - distance between the second jet and the lepton in the eta-phi plane $\Delta R(jet2, lepton)$ smaller 1.2 - invariant mass of the leading jet larger than the invariant mass of the combination of the second jet and lepton The differential and normalised differential cross section was measured as a function of the leading jet mass.

Source code: CMS_2017_I1518399.cc
  1// -*- C++ -*-
  2#include "Rivet/Analysis.hh"
  3#include "Rivet/Projections/PartonicTops.hh"
  4#include "Rivet/Projections/FinalState.hh"
  5#include "Rivet/Projections/VetoedFinalState.hh"
  6#include "Rivet/Projections/PromptFinalState.hh"
  7#include "Rivet/Projections/FastJets.hh"
  8#include "Rivet/Projections/LeptonFinder.hh"
  9#include "Rivet/Projections/ChargedLeptons.hh"
 10
 11namespace Rivet {
 12
 13
 14  /// Leading-jet mass for boosted top quarks at 8 TeV
 15  class CMS_2017_I1518399 : public Analysis {
 16  public:
 17
 18    /// Constructor
 19    RIVET_DEFAULT_ANALYSIS_CTOR(CMS_2017_I1518399);
 20
 21
 22    /// @name Analysis methods
 23    /// @{
 24
 25    /// Book histograms and initialise projections before the run
 26    void init() {
 27
 28      // Dressed leptons
 29      IdentifiedFinalState photons(PID::PHOTON);
 30      ChargedLeptons charged_leptons;
 31      PromptFinalState prompt_leptons(charged_leptons);
 32      Cut leptonCuts = Cuts::pT > 45*GeV && Cuts::abseta < 2.1;
 33      LeptonFinder dressed_leptons(prompt_leptons, photons, 0.1, leptonCuts);
 34      declare(dressed_leptons, "LeptonFinder");
 35
 36      // Jets
 37      VetoedFinalState fs_jets;
 38      fs_jets.vetoNeutrinos();
 39      declare(FastJets(fs_jets, JetAlg::CAM, 1.2), "JetsCA12");
 40
 41      // Partonic top for decay channel definition
 42      declare(PartonicTops(TopDecay::E_MU, PromptEMuFromTau::NO), "LeptonicTops");
 43      declare(PartonicTops(TopDecay::HADRONIC), "HadronicTops");
 44
 45      // Main histograms
 46      book(_hist_mass     , "d01-x01-y01");
 47      book(_hist_mass_norm, "d02-x01-y01");
 48
 49    }
 50
 51
 52    /// Perform the per-event analysis
 53    void analyze(const Event& event) {
 54
 55      // Decay mode check
 56      const Particles& leptonicTops = apply<PartonicTops>(event, "LeptonicTops").particlesByPt();
 57      const Particles& hadronicTops = apply<PartonicTops>(event, "HadronicTops").particlesByPt();
 58      if (leptonicTops.size() != 1 || hadronicTops.size() != 1) vetoEvent;
 59
 60      // Get the leptons
 61      const LeptonFinder& dressed_leptons = apply<LeptonFinder>(event, "LeptonFinder");
 62
 63      // Leading dressed lepton
 64      const DressedLeptons leptons = dressed_leptons.dressedLeptons();
 65      if (leptons.empty()) vetoEvent;
 66      Particle lepton;
 67      for (const Particle& l : leptons) {
 68        if (l.pT() > lepton.pT()) lepton = l;
 69      }
 70
 71      // Get the jets
 72      const Jets& psjetsCA12 = apply<FastJets>(event, "JetsCA12").jetsByPt(Cuts::pT > 50*GeV);
 73
 74      // Subtract the lepton four vector from a jet in case of overlap and clean jets
 75      Jets cleanedJets;
 76      for (Jet jet : psjetsCA12) {
 77        if (deltaR(jet, lepton) < 1.2 )
 78          jet = Jet(jet.momentum()-lepton.momentum(), jet.particles(), jet.tags());
 79        if (jet.abseta() < 2.5) cleanedJets.push_back(jet);
 80      }
 81      std::sort(cleanedJets.begin(), cleanedJets.end(), cmpMomByPt);
 82
 83      // Jet pT cuts
 84      if (cleanedJets.size() < 2) vetoEvent;
 85      if (cleanedJets.at(0).pT() < 400*GeV) vetoEvent;
 86      if (cleanedJets.at(1).pT() < 150*GeV) vetoEvent;
 87
 88      // Jet veto
 89      if (cleanedJets.size() > 2 && cleanedJets.at(2).pT() > 150*GeV) vetoEvent;
 90
 91      // Small distance between 2nd jet and lepton
 92      if (deltaR(cleanedJets.at(1), lepton) > 1.2) vetoEvent;
 93
 94      // m(jet1) > m(jet2 +lepton)
 95      FourMomentum secondJetLepton = cleanedJets.at(1).momentum() + lepton.momentum();
 96      if (cleanedJets.at(0).mass() < secondJetLepton.mass()) vetoEvent;
 97
 98      // Fill histograms
 99      _hist_mass->fill(cleanedJets.at(0).mass());
100      _hist_mass_norm->fill(cleanedJets.at(0).mass());
101    }
102
103
104    /// Normalise histograms etc., after the run
105    void finalize() {
106      const double sf = crossSection()/picobarn * 1000 / sumOfWeights();
107      scale(_hist_mass, sf);
108      normalize(_hist_mass_norm, 1.0, false);
109    }
110
111    /// @}
112
113
114  private:
115
116    // Histograms
117    Histo1DPtr _hist_mass, _hist_mass_norm;
118
119  };
120
121
122  RIVET_DECLARE_PLUGIN(CMS_2017_I1518399);
123
124}