rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

ATLAS_2016_I1478355

Measurement of the bbar dijet dijet cross section at 7 TeV
Experiment: ATLAS (LHC)
Inspire ID: 1478355
Status: VALIDATED
Authors:
  • Lucia Keszeghova
  • Deepak Kar
References: Beams: p+ p+
Beam energies: (3500.0, 3500.0) GeV
Run details:
  • Dijet production, high pThatmin

The dijet production cross section for jets containing a $b$-hadron ($b$-jets) has been measured in proton-proton collisions with a centre-of-mass energy of $\sqrt{s} = 7$ TeV, using the ATLAS detector at the LHC. The data used correspond to an integrated luminosity of 4.2 fb$^{-1}$. The cross section is measured for events with two identified $b$-jets with a transverse momentum $p_\mathrm{T} > 20$GeV and a minimum separation in the $\eta$-$\phi$ plane of $\Delta R$ = 0.4. At least one of the jets in the event is required to have $p_\mathrm{T} > 270$GeV. The cross section is measured differentially as a function of dijet invariant mass, dijet transverse momentum, boost of the dijet system, and the rapidity difference, azimuthal angle and angular distance between the $b$-jets. The results are compared to different predictions of leading order and next-to-leading order perturbative quantum chromodynamics matrix elements supplemented with models for parton-showers and hadronization.

Source code: ATLAS_2016_I1478355.cc
  1// -*- C++ -*-
  2#include "Rivet/Analysis.hh"
  3#include "Rivet/Projections/FastJets.hh"
  4#include "Rivet/Projections/HeavyHadrons.hh"
  5
  6namespace Rivet {
  7
  8
  9  /// @brief Measurement of the bbar dijet dijet cross section at 7 TeV
 10  class ATLAS_2016_I1478355 : public Analysis {
 11  public:
 12
 13    /// Constructor
 14    RIVET_DEFAULT_ANALYSIS_CTOR(ATLAS_2016_I1478355);
 15
 16    /// @name Analysis methods
 17    ///@{
 18
 19    /// Book histograms and initialise projections before the run
 20    void init() {
 21
 22      // Initialise and register projections
 23      FinalState fs(Cuts::abseta < 3.2);
 24      FastJets fj(fs, JetAlg::ANTIKT, 0.4);
 25      fj.useInvisibles();
 26      declare(fj, "Jets");
 27      declare(HeavyHadrons(Cuts::abseta < 3.2 && Cuts::pT > 5*GeV), "BHadrons");
 28
 29      book(_h["m_bb"], 1, 1, 1);
 30      book(_h["Delta_phi"], 2, 1, 1);
 31      book(_h["y_diff"], 3, 1, 1);
 32      book(_h["Delta_R"], 4, 1, 1);
 33      book(_h["pT_bb"], 5, 1, 1);
 34      book(_h["y_B"], 6, 1, 1);
 35    }
 36
 37
 38    /// Perform the per-event analysis
 39    void analyze(const Event& event) {
 40
 41      // Retrieve clustered jets, sorted by pT, with a minimum pT cut
 42      const Jets& jets = apply<FastJets>(event, "Jets").jetsByPt(Cuts::pT > 20*GeV);
 43      const Particles& bHadrons = apply<HeavyHadrons>(event, "BHadrons").bHadrons();
 44
 45      if (jets.size() < 1)  vetoEvent;
 46      if (jets[0].pT() > 270*GeV){
 47
 48        Jets foundBjets;
 49        for (const Jet& j : jets){
 50          for (const Particle& b : bHadrons){
 51            if ((deltaR(j, b) < 0.3) && (j.pT() > 20*GeV) && (fabs(j.eta()) < 2.5)) foundBjets.push_back(j);
 52          }
 53        }
 54
 55        Jet firstBjet;
 56        Jet secondBjet;
 57        bool bJetPair = false;
 58        for (const Jet& b1 : foundBjets){
 59          for (const Jet& b2 : foundBjets){
 60            if(deltaR(b1, b2) > 0.4){
 61              bJetPair = true;
 62              firstBjet = b1;
 63              secondBjet = b2;
 64              break;
 65            }
 66          }
 67          if (bJetPair) break;
 68        }
 69
 70        if (bJetPair) {
 71          const double mass = (firstBjet.momentum() + secondBjet.momentum()).mass();
 72          _h["m_bb"]->fill(mass/GeV);
 73          const double Delta_phi = deltaPhi(firstBjet.phi(), secondBjet.phi());
 74          _h["Delta_phi"]->fill(Delta_phi);
 75          const double y_diff = fabs(firstBjet.rapidity() - secondBjet.rapidity())/2;
 76          _h["y_diff"]->fill(y_diff);
 77          const double Delta_R = deltaR(firstBjet.momentum(), secondBjet.momentum());
 78          _h["Delta_R"]->fill(Delta_R);
 79          const double pT = (firstBjet.momentum() + secondBjet.momentum()).pT();
 80          _h["pT_bb"]->fill(pT/GeV);
 81          const double y_B = (firstBjet.rapidity() + secondBjet.rapidity())/2;
 82          _h["y_B"]->fill(y_B);
 83        }
 84      }
 85    }
 86
 87    /// Normalise histograms etc., after the run
 88    void finalize() {
 89
 90      const double xsec = crossSectionPerEvent()/picobarn; 
 91      scale(_h, xsec);
 92    }
 93
 94    ///@}
 95
 96    /// @name Histograms
 97    ///@{
 98     map<string, Histo1DPtr> _h;
 99    ///@}
100
101  };
102
103  RIVET_DECLARE_PLUGIN(ATLAS_2016_I1478355);
104}