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 ${\rm 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
  2
  3
  4
  5
  6
  7
  8
  9
 10
 11
 12
 13
 14
 15
 16
 17
 18
 19
 20
 21
 22
 23
 24
 25
 26
 27
 28
 29
 30
 31
 32
 33
 34
 35
 36
 37
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
// -*- C++ -*-
#include "Rivet/Analysis.hh"
#include "Rivet/Projections/PartonicTops.hh"
#include "Rivet/Projections/FinalState.hh"
#include "Rivet/Projections/VetoedFinalState.hh"
#include "Rivet/Projections/PromptFinalState.hh"
#include "Rivet/Projections/FastJets.hh"
#include "Rivet/Projections/DressedLeptons.hh"
#include "Rivet/Projections/ChargedLeptons.hh"

namespace Rivet {


  /// Leading jet mass for boosted top quarks at 8 TeV
  class CMS_2017_I1518399 : public Analysis {
  public:

    /// Constructor
    RIVET_DEFAULT_ANALYSIS_CTOR(CMS_2017_I1518399);


    /// @name Analysis methods
    //@{

    /// Book histograms and initialise projections before the run
    void init() {

      // Dressed leptons
      IdentifiedFinalState photons(PID::PHOTON);
      ChargedLeptons charged_leptons;
      PromptFinalState prompt_leptons(charged_leptons);
      Cut leptonCuts = Cuts::pT > 45*GeV && Cuts::abseta < 2.1;
      DressedLeptons dressed_leptons(photons, prompt_leptons, 0.1, leptonCuts);
      declare(dressed_leptons, "DressedLeptons");

      // Jets
      VetoedFinalState fs_jets;
      fs_jets.vetoNeutrinos();
      declare(FastJets(fs_jets, FastJets::CAM, 1.2), "JetsCA12");

      // Partonic top for decay channel definition
      declare(PartonicTops(PartonicTops::DecayMode::E_MU, false), "LeptonicTops");
      declare(PartonicTops(PartonicTops::DecayMode::HADRONIC), "HadronicTops");

      // Main histograms
      book(_hist_mass     , "d01-x01-y01");
      book(_hist_mass_norm, "d02-x01-y01");

    }


    /// Perform the per-event analysis
    void analyze(const Event& event) {

      // Decay mode check
      const Particles& leptonicTops = apply<PartonicTops>(event, "LeptonicTops").particlesByPt();
      const Particles& hadronicTops = apply<PartonicTops>(event, "HadronicTops").particlesByPt();
      if (leptonicTops.size() != 1 || hadronicTops.size() != 1) vetoEvent;

      // Get the leptons
      const DressedLeptons& dressed_leptons = apply<DressedLeptons>(event, "DressedLeptons");

      // Leading dressed lepton
      const vector<DressedLepton> leptons = dressed_leptons.dressedLeptons();
      if (leptons.empty()) vetoEvent;
      Particle lepton;
      for (const Particle& l : leptons) {
        if (l.pT() > lepton.pT()) lepton = l;
      }

      // Get the jets
      const Jets& psjetsCA12 = applyProjection<FastJets>(event, "JetsCA12").jetsByPt(Cuts::pT > 50*GeV);

      // Subtract the lepton four vector from a jet in case of overlap and clean jets
      Jets cleanedJets;
      for (Jet jet : psjetsCA12) {
        if (deltaR(jet, lepton) < 1.2 )
          jet = Jet(jet.momentum()-lepton.momentum(), jet.particles(), jet.tags());
        if (jet.abseta() < 2.5) cleanedJets.push_back(jet);
      }
      std::sort(cleanedJets.begin(), cleanedJets.end(), cmpMomByPt);

      // Jet pT cuts
      if (cleanedJets.size() < 2) vetoEvent;
      if (cleanedJets.at(0).pT() < 400*GeV) vetoEvent;
      if (cleanedJets.at(1).pT() < 150*GeV) vetoEvent;

      // Jet veto
      if (cleanedJets.size() > 2 && cleanedJets.at(2).pT() > 150*GeV) vetoEvent;

      // Small distance between 2nd jet and lepton
      if (deltaR(cleanedJets.at(1), lepton) > 1.2) vetoEvent;

      // m(jet1) > m(jet2 +lepton)
      FourMomentum secondJetLepton = cleanedJets.at(1).momentum() + lepton.momentum();
      if (cleanedJets.at(0).mass() < secondJetLepton.mass()) vetoEvent;

      // Fill histograms
      const double weight = 1.0;
      _hist_mass->fill(cleanedJets.at(0).mass(), weight);
      _hist_mass_norm->fill(cleanedJets.at(0).mass(), weight);
    }


    /// Normalise histograms etc., after the run
    void finalize() {
      const double sf = crossSection() * 1000 / sumOfWeights();
      scale(_hist_mass, sf);
      normalize(_hist_mass_norm, 1.0, false);
    }

    //@}


  private:

    // Histograms
    Histo1DPtr _hist_mass, _hist_mass_norm;

  };


  // The hook for the plugin system
  RIVET_DECLARE_PLUGIN(CMS_2017_I1518399);


}