Rivet Analyses Reference

CMS_2012_I1089835

Inclusive b-jet production in pp collisions at 7 TeV
Experiment: CMS (LHC)
Inspire ID: 1089835
Status: VALIDATED
Authors:
  • cms-pag-conveners-bph@cern.ch
  • Paolo Gunnellini
References:
  • CMS-BPH-11-022
  • JHEP 1204 (2012) 084
  • DOI:10.1007/JHEP04(2012)084
  • arXiv: 1202.4617
Beams: p+ p+
Beam energies: (3500.0, 3500.0) GeV
Run details:
  • QCD events with transverse momentum greater than 10 GeV

The inclusive b-jet production cross section in pp collisions at a center-of-mass energy of 7 TeV is measured using data collected by the CMS experiment at the LHC. The cross section is presented as a function of the jet transverse momentum in the range 18 $<$ pT $<$ 200 GeV for several rapidity intervals. The results are also given as the ratio of the b-jet production cross section to the inclusive jet production cross section. The measurement is performed with two different analyses, which differ in their trigger selection and b-jet identification: a jet analysis that selects events with a b jet using a sample corresponding to an integrated luminosity of 34 inverse picobarns, and a muon analysis requiring a b jet with a muon based on an integrated luminosity of 3 inverse picobarns. In both approaches the b jets are identified by requiring a secondary vertex. The results from the two methods are in agreement with each other and with next-to-leading order calculations, as well as with predictions based on the PYTHIA event generator.

Source code: CMS_2012_I1089835.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
// -*- C++ -*-
#include "Rivet/Analysis.hh"
#include "Rivet/Tools/Logging.hh"
#include "Rivet/Projections/FinalState.hh"
#include "Rivet/Projections/FastJets.hh"

namespace Rivet {
  
  /// @brief Inclusive b-jet production in pp collisions at 7 TeV
  class CMS_2012_I1089835 : public Analysis {
  public:
    
    /// @name Constructors etc.
    //@{
    
    /// Constructor
    CMS_2012_I1089835()
      : Analysis("CMS_2012_I1089835")
    {        
      
    }
    
    //@}
    
    
  public:
    
    /// @name Analysis methods
    //@{
    
    /// Book histograms and initialise projections before the run
    void init() {
      
      const FinalState cnfs;
      declare(cnfs, "FS");
      declare(FastJets(cnfs, FastJets::ANTIKT, 0.5), "Jets");
      
      book(_h_dsigdpty05, 4, 1, 1);
      book(_h_dsigdpty10, 5, 1, 1);
      book(_h_dsigdpty15, 6, 1, 1);
      book(_h_dsigdpty20, 7, 1, 1);
      book(_h_dsigdpty22, 8, 1, 1);
      
    }
    
    
    /// Perform the per-event analysis
    void analyze(const Event& event) {
      
      const FastJets& fastjets = apply<FastJets>(event, "Jets"); 
      const Jets jets = fastjets.jetsByPt(10.);
      
      for (const Jet& j : jets) {
        
        const double ptB= j.pT();
        const double yB= j.rapidity();
        
        if (j.bTagged()) {
          
          if( fabs(yB) < 0.5) { _h_dsigdpty05->fill( ptB, 1.0 );}
          else if( fabs(yB) >= 0.5 && fabs(yB) < 1.0) { _h_dsigdpty10->fill( ptB, 1.0 );}
          else if( fabs(yB) >= 1.0 && fabs(yB) < 1.5) { _h_dsigdpty15->fill( ptB, 1.0 );}
          else if( fabs(yB) >= 1.5 && fabs(yB) < 2.0) { _h_dsigdpty20->fill( ptB, 1.0 );}
          else if( fabs(yB) >= 2.0 && fabs(yB) < 2.2) { _h_dsigdpty22->fill( ptB, 1.0 );}
        }
      }
    }
    
    /// Normalise histograms etc., after the run
    void finalize() {
      
      double invlumi = crossSection()/picobarn/sumOfWeights();
      
      scale(_h_dsigdpty05, invlumi); 
      scale(_h_dsigdpty10, invlumi); 
      scale(_h_dsigdpty15, invlumi); 
      scale(_h_dsigdpty20, invlumi); 
      scale(_h_dsigdpty22, invlumi/0.4); 
      
    }
    
    //@}
    
    
  private:
    
    /// @name Histograms
    //@{
    Histo1DPtr _h_dsigdpty05;
    Histo1DPtr _h_dsigdpty10;
    Histo1DPtr _h_dsigdpty15;
    Histo1DPtr _h_dsigdpty20;
    Histo1DPtr _h_dsigdpty22;
    //@}
    
    
  };
  
  
  
  // The hook for the plugin system
  DECLARE_RIVET_PLUGIN(CMS_2012_I1089835);

}