rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

CMS_2016_I1486238

Studies of 2 $b$-jet + 2 jet production in proton-proton collisions at 7 TeV
Experiment: CMS (LHC)
Inspire ID: 1486238
Status: VALIDATED
Authors:
  • P. Gunnellini
References:
  • CMS-FSQ-13-010
  • CERN-EP-2016-191
  • arXiv: 1609.03489
  • Phys.Rev. D94 (2016) no.11, 112005
Beams: p+ p+
Beam energies: (3500.0, 3500.0) GeV
Run details:
  • Hard QCD events with p_\perp cut at generator level of 15 GeV. Needs 2 $b$-jets with $p_\perp > 20$ GeV in $|\eta|$ $<$ 2.4; 2 jets with $p_\perp > 20$ GeV in $|\eta|$ $<$ 4.7

Measurements are presented of the cross section for the production of at least four jets, of which at least two originate from $b$ quarks, in proton-proton collisions. Data collected with the CMS detector at the LHC at a center-of-mass energy of 7 TeV are used, corresponding to an integrated luminosity of 3 inverse picobarns. The cross section is measured as a function of the jet transverse momentum for p_\perp 20 GeV, and of the jet pseudorapidity for $|\eta|<2.4$ ($b$ jets) and 4.7 (untagged jets). The correlations in azimuthal angle and p_\perp between the jets are also studied. The $\eta$ and p_\perp distributions of the four jets and the correlations between them are well reproduced by event generators that combine perturbative QCD calculations at next-to-leading-order accuracy with contributions from parton showers and multiparton interactions.

Source code: CMS_2016_I1486238.cc
  1// -*- C++ -*-
  2#include "Rivet/Analysis.hh"
  3#include "Rivet/Projections/FinalState.hh"
  4#include "Rivet/Projections/FastJets.hh"
  5
  6#define I_KNOW_THE_INITIAL_QUARKS_PROJECTION_IS_DODGY_BUT_NEED_TO_USE_IT
  7#include "Rivet/Projections/InitialQuarks.hh"
  8
  9namespace Rivet {
 10
 11
 12  /// Studies of 2 b-jet + 2 jet production in proton-proton collisions at 7 TeV
 13  class CMS_2016_I1486238 : public Analysis {
 14  public:
 15
 16    /// Constructor
 17    RIVET_DEFAULT_ANALYSIS_CTOR(CMS_2016_I1486238);
 18
 19
 20    /// @name Analysis methods
 21    //@{
 22
 23    /// Book histograms and initialise projections before the run
 24    void init() {
 25
 26      FastJets akt(FinalState(), FastJets::ANTIKT, 0.5);
 27      declare(akt, "antikT");
 28
 29      book(_h_Deltaphi_newway, 1,1,1);
 30      book(_h_deltaphiafterlight, 9,1,1);
 31      book(_h_SumPLight, 5,1,1);
 32
 33      book(_h_LeadingBJetpt, 11,1,1);
 34      book(_h_SubleadingBJetpt, 15,1,1);
 35      book(_h_LeadingLightJetpt, 13,1,1);
 36      book(_h_SubleadingLightJetpt, 17,1,1);
 37
 38      book(_h_LeadingBJeteta, 10,1,1);
 39      book(_h_SubleadingBJeteta, 14,1,1);
 40      book(_h_LeadingLightJeteta, 12,1,1);
 41      book(_h_SubleadingLightJeteta, 16,1,1);
 42    }
 43
 44
 45    /// Perform the per-event analysis
 46    void analyze(const Event& event) {
 47      const Jets& jets = apply<JetAlg>(event, "antikT").jetsByPt(Cuts::absrap < 4.7 && Cuts::pT > 20*GeV);
 48      if (jets.size() < 4) vetoEvent;
 49
 50      // Initial quarks
 51      /// @note Quark-level tagging...
 52      Particles bquarks;
 53      for (ConstGenParticlePtr p : HepMCUtils::particles(event.genEvent())) {
 54        if (abs(p->pdg_id()) == PID::BQUARK) bquarks += Particle(p);
 55      }
 56      Jets bjets, ljets;
 57      for (const Jet& j : jets) {
 58        const bool btag = any(bquarks, deltaRLess(j, 0.3));
 59        // for (const Particle& b : bquarks) if (deltaR(j, b) < 0.3) btag = true;
 60        (btag && j.abseta() < 2.4 ? bjets : ljets).push_back(j);
 61      }
 62
 63      // Fill histograms
 64      const double weight = 1.0;
 65      if (bjets.size() >= 2 && ljets.size() >= 2) {
 66        _h_LeadingBJetpt->fill(bjets[0].pT()/GeV, weight);
 67        _h_SubleadingBJetpt->fill(bjets[1].pT()/GeV, weight);
 68        _h_LeadingLightJetpt->fill(ljets[0].pT()/GeV, weight);
 69        _h_SubleadingLightJetpt->fill(ljets[1].pT()/GeV, weight);
 70        //
 71        _h_LeadingBJeteta->fill(bjets[0].eta(), weight);
 72        _h_SubleadingBJeteta->fill(bjets[1].eta(), weight);
 73        _h_LeadingLightJeteta->fill(ljets[0].eta(), weight);
 74        _h_SubleadingLightJeteta->fill(ljets[1].eta(), weight);
 75
 76        const double lightdphi = deltaPhi(ljets[0], ljets[1]);
 77        _h_deltaphiafterlight->fill(lightdphi, weight);
 78
 79        const double vecsumlightjets = sqrt(sqr(ljets[0].px()+ljets[1].px()) + sqr(ljets[0].py()+ljets[1].py())); //< @todo Just (lj0+lj1).pT()? Or use add_quad
 80        const double term2 = vecsumlightjets/(sqrt(sqr(ljets[0].px()) + sqr(ljets[0].py())) + sqrt(sqr(ljets[1].px()) + sqr(ljets[1].py()))); //< @todo lj0.pT() + lj1.pT()? Or add_quad
 81        _h_SumPLight->fill(term2, weight);
 82
 83        const double pxBsyst2 = bjets[0].px()+bjets[1].px(); // @todo (bj0+bj1).px()
 84        const double pyBsyst2 = bjets[0].py()+bjets[1].py(); // @todo (bj0+bj1).py()
 85        const double pxJetssyst2 = ljets[0].px()+ljets[1].px(); // @todo (lj0+lj1).px()
 86        const double pyJetssyst2 = ljets[0].py()+ljets[1].py(); // @todo (lj0+lj1).py()
 87        const double modulusB2 = sqrt(sqr(pxBsyst2)+sqr(pyBsyst2)); //< @todo add_quad
 88        const double modulusJets2 = sqrt(sqr(pxJetssyst2)+sqr(pyJetssyst2)); //< @todo add_quad
 89        const double cosphiBsyst2 = pxBsyst2/modulusB2;
 90        const double cosphiJetssyst2 = pxJetssyst2/modulusJets2;
 91        const double phiBsyst2 = ((pyBsyst2 > 0) ? 1 : -1) * acos(cosphiBsyst2); //< @todo sign(pyBsyst2)
 92        const double phiJetssyst2 = sign(pyJetssyst2) * acos(cosphiJetssyst2);
 93        const double Dphi2 = deltaPhi(phiBsyst2, phiJetssyst2);
 94        _h_Deltaphi_newway->fill(Dphi2,weight);
 95      }
 96    }
 97
 98
 99    /// Normalise histograms etc., after the run
100    void finalize() {
101      const double invlumi = crossSection()/picobarn/sumOfWeights();
102      normalize(_h_SumPLight);  normalize(_h_deltaphiafterlight); normalize(_h_Deltaphi_newway);
103      scale(_h_LeadingLightJetpt, invlumi); scale(_h_SubleadingLightJetpt, invlumi); 
104      scale(_h_LeadingBJetpt, invlumi); scale(_h_SubleadingBJetpt, invlumi);
105      scale(_h_LeadingLightJeteta, invlumi); scale(_h_SubleadingLightJeteta, invlumi);
106      scale(_h_LeadingBJeteta, invlumi); scale(_h_SubleadingBJeteta, invlumi);
107    }
108
109    //@}
110
111
112  private:
113
114    /// @name Histograms
115    //@{
116
117    Histo1DPtr _h_deltaphiafterlight, _h_Deltaphi_newway, _h_SumPLight;
118    Histo1DPtr _h_LeadingBJetpt, _h_SubleadingBJetpt, _h_LeadingLightJetpt, _h_SubleadingLightJetpt;
119    Histo1DPtr _h_LeadingBJeteta, _h_SubleadingBJeteta, _h_LeadingLightJeteta, _h_SubleadingLightJeteta;
120
121  };
122
123
124  // Hook for the plugin system
125  RIVET_DECLARE_PLUGIN(CMS_2016_I1486238);
126
127}