rivet is hosted by Hepforge, IPPP Durham

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: 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// -*- C++ -*-
 2#include "Rivet/Analysis.hh"
 3#include "Rivet/Tools/Logging.hh"
 4#include "Rivet/Projections/FinalState.hh"
 5#include "Rivet/Projections/FastJets.hh"
 6
 7namespace Rivet {
 8
 9
10  /// @brief Inclusive b-jet production in pp collisions at 7 TeV
11  class CMS_2012_I1089835 : public Analysis {
12  public:
13
14    /// Constructor
15    RIVET_DEFAULT_ANALYSIS_CTOR(CMS_2012_I1089835);
16
17
18    /// @name Analysis methods
19    /// @{
20
21    /// Book histograms and initialise projections before the run
22    void init() {
23
24      const FinalState cnfs;
25      declare(cnfs, "FS");
26      declare(FastJets(cnfs, JetAlg::ANTIKT, 0.5), "Jets");
27
28      book(_h_dsigdpty05, 4, 1, 1);
29      book(_h_dsigdpty10, 5, 1, 1);
30      book(_h_dsigdpty15, 6, 1, 1);
31      book(_h_dsigdpty20, 7, 1, 1);
32      book(_h_dsigdpty22, 8, 1, 1);
33
34    }
35
36
37    /// Perform the per-event analysis
38    void analyze(const Event& event) {
39
40      const FastJets& fastjets = apply<FastJets>(event, "Jets");
41      const Jets jets = fastjets.jetsByPt(Cuts::pT > 10*GeV);
42
43      for (const Jet& j : jets) {
44
45        const double ptB= j.pT();
46        const double yB= j.rapidity();
47
48        if (j.bTagged()) {
49
50          if( fabs(yB) < 0.5) { _h_dsigdpty05->fill( ptB, 1.0 );}
51          else if( fabs(yB) >= 0.5 && fabs(yB) < 1.0) { _h_dsigdpty10->fill( ptB, 1.0 );}
52          else if( fabs(yB) >= 1.0 && fabs(yB) < 1.5) { _h_dsigdpty15->fill( ptB, 1.0 );}
53          else if( fabs(yB) >= 1.5 && fabs(yB) < 2.0) { _h_dsigdpty20->fill( ptB, 1.0 );}
54          else if( fabs(yB) >= 2.0 && fabs(yB) < 2.2) { _h_dsigdpty22->fill( ptB, 1.0 );}
55        }
56      }
57    }
58
59    /// Normalise histograms etc., after the run
60    void finalize() {
61
62      double invlumi = crossSection()/picobarn/sumOfWeights();
63
64      scale(_h_dsigdpty05, invlumi);
65      scale(_h_dsigdpty10, invlumi);
66      scale(_h_dsigdpty15, invlumi);
67      scale(_h_dsigdpty20, invlumi);
68      scale(_h_dsigdpty22, invlumi/0.4);
69
70    }
71
72    /// @}
73
74
75  private:
76
77    /// @name Histograms
78    /// @{
79    Histo1DPtr _h_dsigdpty05;
80    Histo1DPtr _h_dsigdpty10;
81    Histo1DPtr _h_dsigdpty15;
82    Histo1DPtr _h_dsigdpty20;
83    Histo1DPtr _h_dsigdpty22;
84    /// @}
85
86  };
87
88
89  RIVET_DECLARE_PLUGIN(CMS_2012_I1089835);
90
91}