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:
  • Eur.Phys.J. C76 (2016) 670
  • DOI:10.1140/epjc/s10052-016-4521-y
  • arXiv: 1607.08430
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
  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
// -*- C++ -*-
#include "Rivet/Analysis.hh"
#include "Rivet/Projections/FastJets.hh"
#include "Rivet/Projections/HeavyHadrons.hh"
#include "Rivet/Tools/BinnedHistogram.hh"

namespace Rivet {


  /// @brief Measurement of the bbar dijet dijet cross section at 7 TeV
  class ATLAS_2016_I1478355 : public Analysis {
  public:

    /// Constructor
    RIVET_DEFAULT_ANALYSIS_CTOR(ATLAS_2016_I1478355);

    /// @name Analysis methods
    ///@{

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

      // Initialise and register projections
      FinalState fs(Cuts::abseta < 3.2);
      FastJets fj(fs, FastJets::ANTIKT, 0.4);
      fj.useInvisibles();
      declare(fj, "Jets");
      declare(HeavyHadrons(Cuts::abseta < 3.2 && Cuts::pT > 5*GeV), "BHadrons");

      book(_h["m_bb"], 1, 1, 1);
      book(_h["Delta_phi"], 2, 1, 1);
      book(_h["y_diff"], 3, 1, 1);
      book(_h["Delta_R"], 4, 1, 1);
      book(_h["pT_bb"], 5, 1, 1);
      book(_h["y_B"], 6, 1, 1);
    }


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

      // Retrieve clustered jets, sorted by pT, with a minimum pT cut
      const Jets& jets = apply<FastJets>(event, "Jets").jetsByPt(20*GeV);
      const Particles& bHadrons = apply<HeavyHadrons>(event, "BHadrons").bHadrons();

      if (jets.size() < 1)  vetoEvent;
      if (jets[0].pT() > 270*GeV){

        Jets foundBjets;
        for (const Jet& j : jets){
          for (const Particle& b : bHadrons){
            if ((deltaR(j, b) < 0.3) && (j.pT() > 20*GeV) && (fabs(j.eta()) < 2.5)) foundBjets.push_back(j);
          }
        }

        Jet firstBjet;
        Jet secondBjet;
        bool bJetPair = false;
        for (const Jet& b1 : foundBjets){
          for (const Jet& b2 : foundBjets){
            if(deltaR(b1, b2) > 0.4){
              bJetPair = true;
              firstBjet = b1;
              secondBjet = b2;
              break;
            }
          }
          if (bJetPair) break;
        }

        if (bJetPair) {
          const double mass = (firstBjet.momentum() + secondBjet.momentum()).mass();
          _h["m_bb"]->fill(mass/GeV);
          const double Delta_phi = deltaPhi(firstBjet.phi(), secondBjet.phi());
          _h["Delta_phi"]->fill(Delta_phi);
          const double y_diff = fabs(firstBjet.rapidity() - secondBjet.rapidity())/2;
          _h["y_diff"]->fill(y_diff);
          const double Delta_R = deltaR(firstBjet.momentum(), secondBjet.momentum());
          _h["Delta_R"]->fill(Delta_R);
          const double pT = (firstBjet.momentum() + secondBjet.momentum()).pT();
          _h["pT_bb"]->fill(pT/GeV);
          const double y_B = (firstBjet.rapidity() + secondBjet.rapidity())/2;
          _h["y_B"]->fill(y_B);
        }
      }
    }

    /// Normalise histograms etc., after the run
    void finalize() {

      const double xsec = crossSectionPerEvent()/picobarn; 
      scale(_h, xsec);
    }

    ///@}

    /// @name Histograms
    ///@{
     map<string, Histo1DPtr> _h;
    ///@}

  };

  DECLARE_RIVET_PLUGIN(ATLAS_2016_I1478355);
}