Rivet Analyses Reference

CMS_2018_I1620050

Measurement of normalized differential ttbar cross sections in the dilepton channel from pp collisions at 13 TeV
Experiment: CMS (LHC)
Inspire ID: 1620050
Status: VALIDATED
Authors:
  • Youn Jung Roh
  • Suyong Choi
  • Junghwan Goh
  • Dajeong Jeon
  • Jason S. H. Lee
References:
  • JHEP 1804 (2018) 060
  • DOI:10.1007/JHEP04(2018)060
  • arXiv: 1708.07638
  • https://www.hepdata.net/record/81686
  • CMS-TOP-16-007
Beams: p+ p+
Beam energies: (6500.0, 6500.0) GeV
Run details:
  • pp QCD interactions at $\sqrt{s} = 13$ TeV. Data collected by CMS during the year 2015. Selection of dilepton top pair candidate events at particle level.

Abstract: Normalized differential cross sections for top quark pair production are measured in the dilepton ($e^{+}e^{-}$, $\mu^{+}\mu^{-}$, and $\mu^{\mp}e^{\pm}$) decay channels in proton-proton collisions at a center-of-mass energy of 13TeV. The measurements are performed with data corresponding to an integrated luminosity of 2.1$fb^{-1}$ using the CMS detector at the LHC. The cross sections are measured differentially as a function of the kinematic properties of the leptons, jets from bottom quark hadronization, top quarks, and top quark pairs at the particle and parton levels. The results are compared to several Monte Carlo generators that implement calculations up to next-to-leading order in perturbative quantum chromodynamics interfaced with parton showering, and also to fixed-order theoretical calculations of top quark pair production up to next-to-next-to-leading order. Rivet: This analysis is to be run on $\text{t}\bar{\text{t}}$ Monte Carlo. The particle-level phase space is defined using the following definitions: \begin{description} \item[lepton]: an electron or muon with $p_\text{T}>30\,\text{GeV}$ and $|\eta|<2.4$, dressed within a cone of radius 0.1, \item[jet]: a jet is reconstructed with the anti-$k_t$ algorithm with a radius of 0.4, after removing the neutrinos and dressed leptons, with $p_\text{T]>30\,\text{GeV}$ and $|\eta|<2.4$, \item[b-jet]: a jet that contains a B-hadron. \end{description} A W boson at the particle level is defined by combining a dressed lepton and a neutrino. In each event, a pair of particle-level W bosons is chosen among the possible combinations such that the sum of the absolute values of the invariant mass differences with respect to the W boson mass is minimal. Similarly, a top quark at the particle level is defined by combining a particle-level W boson and a b jet. The combination of a W boson and a b jet with the minimum invariant mass difference from the correct top quark mass is selected.

Source code: CMS_2018_I1620050.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
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
#include "Rivet/Analysis.hh"
#include "Rivet/Projections/FinalState.hh"
#include "Rivet/Projections/FastJets.hh"
#include "Rivet/Projections/PartonicTops.hh"
#include "Rivet/Projections/DressedLeptons.hh"
#include "Rivet/Projections/IdentifiedFinalState.hh"
#include "Rivet/Projections/PromptFinalState.hh"
#include "Rivet/Projections/VetoedFinalState.hh"

namespace Rivet {


  /// Normalized dilepton ttbar differential cross-sections in pp collisions at 13 TeV
  class CMS_2018_I1620050 : public Analysis {
  public:

    RIVET_DEFAULT_ANALYSIS_CTOR(CMS_2018_I1620050);

    void init() {

      // Parton level top quark to analyze dilepton channels only
      // Note: 2nd argument of PartonicTops to toggle tau->lepton channel (true to inclusive, false to exclusive)
      const bool acceptTauDecays = false;
      declare(PartonicTops(PartonicTops::DecayMode::MUON, acceptTauDecays), "PartonTopsToMuon"); // Partonic top decaying to mu
      declare(PartonicTops(PartonicTops::DecayMode::ELECTRON, acceptTauDecays), "PartonTopsToElectron"); // Partonic top decaying to e

      // Build particle level tops starting from FinalState
      const FinalState fs(Cuts::pT > 0. && Cuts::abseta < 6.);

      // Neutrinos
      IdentifiedFinalState neutrinos(fs);
      neutrinos.acceptNeutrinos();
      PromptFinalState prompt_neutrinos(neutrinos, true, true);
      declare(prompt_neutrinos, "Neutrinos");

      // Projection for electrons and muons
      Cut leptonCuts = Cuts::pt > 20*GeV && Cuts::abseta < 2.4;

      PromptFinalState fsLepton(fs);
      fsLepton.acceptMuonDecays(true);
      fsLepton.acceptTauDecays(true);
      SpecialDressedLeptons dressedLeptons(fsLepton, leptonCuts);
      declare(dressedLeptons, "DressedLeptons");

      // Projection for jets
      VetoedFinalState fs_jets(fs);
      fs_jets.addVetoOnThisFinalState(dressedLeptons);
      fs_jets.vetoNeutrinos();
      declare(FastJets(fs_jets, FastJets::ANTIKT, 0.4), "ak4jets");

      // Book hists
      book(_hist_lep_pt, "d01-x01-y01");
      book(_hist_jet_pt, "d02-x01-y01");
      book(_hist_top_pt, "d03-x01-y01");
      book(_hist_top_y, "d04-x01-y01");
      book(_hist_tt_pt, "d05-x01-y01");
      book(_hist_tt_y, "d06-x01-y01");
      book(_hist_tt_m, "d07-x01-y01");
      book(_hist_tt_dphi, "d08-x01-y01");
    }


    void analyze(const Event& event) {

      // Do the analysis only for the full-dleptonic channel
      const Particles partonTopsToMuon     = apply<ParticleFinder>(event, "PartonTopsToMuon").particles();
      //const Particles partonTopsToElectron = apply<ParticleFinder>(event, "PartonTopsToElectron").particles();
      Particles partonTopsToElectron;
      for (const Particle& x : apply<ParticleFinder>(event, "PartonTopsToElectron").particles() ) {
        bool isDuplicated = false;
        for (const Particle& y : partonTopsToMuon ) {
          if ( std::abs(x.pt()-y.pt()) < 0.01 and deltaR(x, y) < 0.01 ) {
            isDuplicated = true;
            break;
          }
        }
        if ( !isDuplicated ) partonTopsToElectron.push_back(x);
      }
      const int nPartonElectrons = partonTopsToElectron.size();
      const int nPartonMuons     = partonTopsToMuon.size();
      if ( nPartonElectrons+nPartonMuons != 2 ) vetoEvent;

      // Select leptons
      const vector<DressedLepton>& dressedLeptons = apply<SpecialDressedLeptons>(event, "DressedLeptons").dressedLeptons();
      if ( dressedLeptons.size() < 2 ) vetoEvent;
      sortByPt(dressedLeptons);

      const FourMomentum& lepton1 = dressedLeptons[0].momentum();
      const FourMomentum& lepton2 = dressedLeptons[1].momentum();
      const int channel = dressedLeptons[0].abspid() + dressedLeptons[1].abspid();
      if ( !((channel == 22 and nPartonElectrons == 2) or
             (channel == 24 and nPartonElectrons == 1 and nPartonMuons == 1) or
             (channel == 26 and nPartonMuons == 2)) ) vetoEvent;

      // Select neutrinos
      const Particles neutrinos = apply<PromptFinalState>(event, "Neutrinos").particlesByPt();
      if ( neutrinos.size() < 2 ) vetoEvent;

      // Select bjets
      const FastJets& fjJets = apply<FastJets>(event, "ak4jets");
      const Jets jets = fjJets.jetsByPt(Cuts::abseta < 2.4 && Cuts::pT > 30*GeV);
      const Jets bJets = filter_select(jets, hasBTag());
      // Jets bJets;
      // for ( Jets::const_iterator itjet = jets.begin(); itjet != jets.end() ; ++itjet) {
      //   if ( itjet->bTagged() ) { // Note: default b tagging algorithm is ghost association (see the Rivet Jet class reference manual)
      //     bJets.push_back(*itjet);
      //   }
      // }
      // sortByPt(bJets);
      // There should at least two b jets.
      if ( bJets.size() < 2 ) vetoEvent;

      // Construct particle level top
      FourMomentum nu1 = neutrinos[0].momentum();
      FourMomentum nu2 = neutrinos[1].momentum();
      if ( std::abs((lepton1+nu1).mass()-80.4) + std::abs((lepton2+nu2).mass()-80.4) >
           std::abs((lepton1+nu2).mass()-80.4) + std::abs((lepton2+nu1).mass()-80.4) ) {
        std::swap(nu1, nu2);
      }
      const FourMomentum w1 = lepton1 + nu1;
      const FourMomentum w2 = lepton2 + nu2;

      FourMomentum bjet1 = bJets[0].momentum();
      FourMomentum bjet2 = bJets[1].momentum();
      if ( std::abs((w1+bjet1).mass()-172.5) + std::abs((w2+bjet2).mass()-172.5) >
           std::abs((w1+bjet2).mass()-172.5) + std::abs((w2+bjet1).mass()-172.5) ) {
        std::swap(bjet1, bjet2);
      }
      const FourMomentum t1 = w1 + bjet1;
      const FourMomentum t2 = w2 + bjet2;
      const FourMomentum tt = t1+t2;

      _hist_lep_pt->fill(lepton1.pt());
      _hist_lep_pt->fill(lepton2.pt());
      _hist_jet_pt->fill(bjet1.pt());
      _hist_jet_pt->fill(bjet2.pt());

      _hist_top_pt->fill(t1.pt());
      _hist_top_pt->fill(t2.pt());
      _hist_top_y->fill(t1.rapidity());
      _hist_top_y->fill(t2.rapidity());

      _hist_tt_pt->fill(tt.pt());
      _hist_tt_y->fill(tt.rapidity());
      _hist_tt_m->fill(tt.mass());
      _hist_tt_dphi->fill(deltaPhi(t1.phi(), t2.phi()));
    }


    /// Normalise histograms etc., after the run
    void finalize() {
      normalize(_hist_lep_pt);
      normalize(_hist_jet_pt);
      normalize(_hist_top_pt);
      normalize(_hist_top_y);
      normalize(_hist_tt_pt);
      normalize(_hist_tt_y);
      normalize(_hist_tt_m);
      normalize(_hist_tt_dphi);
    }


    /// @brief Special dressed-lepton finder
    ///
    /// Find dressed leptons by clustering all leptons and photons
    class SpecialDressedLeptons : public FinalState {
    public:
      /// The default constructor. May specify cuts
      SpecialDressedLeptons(const FinalState& fs, const Cut& cut)
        : FinalState(cut)
      {
        setName("SpecialDressedLeptons");
        IdentifiedFinalState ifs(fs);
        ifs.acceptIdPair(PID::PHOTON);
        ifs.acceptIdPair(PID::ELECTRON);
        ifs.acceptIdPair(PID::MUON);
        declare(ifs, "IFS");
        declare(FastJets(ifs, FastJets::ANTIKT, 0.1), "LeptonJets");
      }

      /// Clone on the heap.
      virtual unique_ptr<Projection> clone() const {
        return unique_ptr<Projection>(new SpecialDressedLeptons(*this));
      }

      /// Retrieve the dressed leptons
      const vector<DressedLepton>& dressedLeptons() const { return _clusteredLeptons; }

      /// Perform the calculation
      void project(const Event& e) {
        _theParticles.clear();
        _clusteredLeptons.clear();

        vector<DressedLepton> allClusteredLeptons;
        const Jets jets = applyProjection<FastJets>(e, "LeptonJets").jetsByPt(5*GeV);
        for (const Jet& jet : jets) {
          Particle lepCand;
          for (const Particle& cand : jet.particles()) {
            const int absPdgId = cand.abspid();
            if (absPdgId == PID::ELECTRON || absPdgId == PID::MUON) {
              if (cand.pt() > lepCand.pt()) lepCand = cand;
            }
          }

          if (!lepCand.isChargedLepton()) continue;

          DressedLepton lepton = DressedLepton(lepCand);
          for (const Particle& cand : jet.particles()) {
            if (isSame(cand, lepCand)) continue;
            lepton.addConstituent(cand, true);
          }
          allClusteredLeptons.push_back(lepton);
        }

        for (const DressedLepton& lepton : allClusteredLeptons) {
          if (accept(lepton)) {
            _clusteredLeptons.push_back(lepton);
            _theParticles.push_back(lepton.constituentLepton());
            _theParticles += lepton.constituentPhotons();
          }
        }
      }
    private:

      /// Container which stores the clustered lepton objects
      vector<DressedLepton> _clusteredLeptons;

    };

  private:

    Histo1DPtr _hist_lep_pt;
    Histo1DPtr _hist_jet_pt;
    Histo1DPtr _hist_top_pt;
    Histo1DPtr _hist_top_y;
    Histo1DPtr _hist_tt_pt;
    Histo1DPtr _hist_tt_y;
    Histo1DPtr _hist_tt_m;
    Histo1DPtr _hist_tt_dphi;

  };


  RIVET_DECLARE_PLUGIN(CMS_2018_I1620050);

}