rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

LHCB_2018_I1662483

Measurement of the forward top pair production cross-section in the dilepton channel
Experiment: LHCB (LHC)
Inspire ID: 1662483
Status: VALIDATED
Authors:
  • Stephen Farry
References: Beams: p+ p+
Beam energies: (6500.0, 6500.0) GeV
Run details:
  • Proton-proton collisions at c.o.m. energy of 13 TeV producing hard QCD events containing a top--anti-top quark pair decaying into dilepton channel accompanied by a b-jet in the final state.

Measurement of the top pair cross-section in the dilepton channel, where a muon, electron and one $b$-jet are required to be present in the event. A fiducial cross-section is measuremed where the leptons and jets are required to be within the pseudorapidity acceptance of the LHCb acceptance ($2<\eta<4.5$), and to have a transverse momentum greater than 20 GeV. Additionally leptons are required to be separated from each other by a radial distance ($\Delta R$) of 0.1, and from the $b$-jet by a distance of 0.5.

Source code: LHCB_2018_I1662483.cc
  1// -*- C++ -*-
  2#include "Rivet/Analysis.hh"
  3#include "Rivet/Projections/FinalState.hh"
  4#include "Rivet/Projections/FinalPartons.hh"
  5#include "Rivet/Projections/FastJets.hh"
  6#include "Rivet/Projections/LeptonFinder.hh"
  7#include "Rivet/Projections/MissingMomentum.hh"
  8#include "Rivet/Projections/PromptFinalState.hh"
  9
 10namespace Rivet {
 11
 12
 13  /// @brief Measurement of forward top pair production in the mu e channel at LHCb
 14  ///
 15  /// Both leptons with pt > 20 GeV and in pseudorapidity range 2.0 - 4.5.
 16  /// Leptons are separated from each other by a dR of 0.1.
 17  /// Leading jet is required to be a b-jet with pt > 20 GeV and pseudorapidity in range 2.2 - 4.2.
 18  /// The jet is separated from both leptons by dR of 0.5.
 19  ///
 20  class LHCB_2018_I1662483 : public Analysis {
 21  public:
 22
 23    /// Default constructor
 24    RIVET_DEFAULT_ANALYSIS_CTOR(LHCB_2018_I1662483);
 25
 26
 27    /// @name Analysis methods
 28    /// @{
 29
 30    /// Book histograms and initialise projections before the run
 31    void init() {
 32
 33      // Initialise and register projections
 34
 35      // the basic final-state projection:
 36      // all final-state particles
 37      const FinalState fs;
 38
 39      //get the final state b quarks for b-tagging
 40      FinalPartons bquarks(Cuts::abspid == 5 && Cuts::pT > 5*GeV);
 41      declare(bquarks, "bquarks");
 42
 43      // the final-state particles declared above are clustered using FastJet with
 44      // the anti-kT algorithm and a jet-radius parameter 0.5
 45      // muons and neutrinos are excluded from the clustering
 46      FastJets jetfs(fs, JetAlg::ANTIKT, 0.5, JetMuons::NONE, JetInvisibles::NONE);
 47      declare(jetfs, "jets");
 48
 49      // FinalState of bare muons and electrons in the event
 50      Cut lepton_cuts = Cuts::pT>=20*GeV && Cuts::etaIn(2.0, 4.5);
 51      FinalState fsmuons(Cuts::abspid == PID::MUON && lepton_cuts );
 52      FinalState fselectrons(Cuts::abspid == PID::ELECTRON && lepton_cuts);
 53
 54      declare(fsmuons, "muons");
 55      declare(fselectrons, "electrons");
 56
 57      // Inclusive cross-section measurement
 58      book(_h_fiducial_xsect ,1,1,1);
 59    }
 60
 61
 62    /// Perform the per-event analysis
 63    void analyze(const Event& event) {
 64
 65      // retrieve leptons, sorted by pT
 66      const FinalState& muons = apply<FinalState>(event, "muons");
 67      const FinalState& electrons = apply<FinalState>(event, "electrons");
 68      const FinalPartons& bquarks = apply<FinalPartons>(event, "bquarks");
 69
 70      // retrieve clustered jets, sorted by pT, with a minimum pT cut
 71      Jets jets = apply<FastJets>(event, "jets").jetsByPt(Cuts::pT > 20*GeV && Cuts::etaIn(2.2,4.2));
 72
 73      if (jets.empty()) vetoEvent;
 74      bool pass = false;
 75      for (const auto& m : muons.particles()) {
 76        for (const auto& e : electrons.particles()) {
 77          if (deltaR(m.momentum(), e.momentum()) < 0.1 ) continue;
 78          vector<Jet> lepton_jets;
 79          for (const auto& j : jets ){
 80            if (deltaR(j.momentum(), m.momentum()) > 0.5 &&
 81                deltaR(j.momentum(), e.momentum()) > 0.5)
 82              lepton_jets.push_back(j);
 83            if (lepton_jets.size() > 0){
 84              for (const auto& b : bquarks.particles()) {
 85                if ( deltaR(b.momentum(), lepton_jets.at(0).momentum()) < 0.5 )  {
 86                  pass = true;
 87                }
 88              }
 89            }
 90          }
 91        }
 92      }
 93
 94      // veto event if it doesn't pass our selection
 95      if (!pass) vetoEvent;
 96      // append to cross-section
 97      _h_fiducial_xsect->fill(13000);
 98    }
 99
100
101    /// Normalise histograms etc., after the run
102    void finalize() {
103      scale(_h_fiducial_xsect, crossSection()/femtobarn/sumOfWeights()); // norm to cross section
104    }
105
106
107    /// Histogram
108    BinnedHistoPtr<int> _h_fiducial_xsect;
109
110  };
111
112
113  RIVET_DECLARE_PLUGIN(LHCB_2018_I1662483);
114
115}