Rivet Analyses Reference

CMS_2015_I1380605

Per-event yield of the highest transverse momentum charged particle and charged-particle jet
Experiment: CMS (LHC)
Inspire ID: 1380605
Status: VALIDATED
Authors:
  • Anastasia Grebenyuk
  • Hannes Jung
References:Beams: p+ p+
Beam energies: (4000.0, 4000.0) GeV
Run details:
  • $sqrt{s}=8$ TeV, minbias events, no p_\perp cut on partons

The per-event yield of the highest transverse momentum charged particle and charged-particle jet, integrated above a given $p_{T min}$ threshold starting at $p_{T min}=0.8$ and $1$ GeV, respectively, is studied in $pp$ collisions at $\sqrt{s} = 8$ \text{TeV}. The particles and the jets are measured in the pseudorapidity ranges $|\eta| < 2.4$ and $1.9$, respectively.

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

namespace Rivet {


  /// Per-event yield of the highest transverse momentum charged particle and charged-particle jet
  class CMS_2015_I1380605 : public Analysis {
  public:

    RIVET_DEFAULT_ANALYSIS_CTOR(CMS_2015_I1380605);


    /// @name Analysis methods
    //@{

    /// Book histograms and initialise projections before the run
    void init() {
      const ChargedFinalState cfs((Cuts::etaIn(-7., 7.)));
      declare(cfs, "CFS");
      declare(FastJets(cfs, FastJets::ANTIKT, 0.5), "Jets");

      book(_h_tracks, 1, 1, 1);
      book(_h_jets  , 2, 1, 1);
      book(_ntracks, "ntracks");
    }


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

      // Veto events without forward activity on both sides
      const ChargedFinalState& cfs = apply<ChargedFinalState>(event, "CFS");
      const size_t count_plus  = cfs.particles(Cuts::eta > 5.3  && Cuts::eta < 6.5).size();
      const size_t count_minus = cfs.particles(Cuts::eta < -5.3 && Cuts::eta > -6.5).size();
      if (!(count_plus > 0 || count_minus > 0)) vetoEvent; //< @todo Should this really also veto the jet analysis?
      /// @warning Needs migration to an AO Counter
      /// @note This isn't the number of tracks, it's the sum of event weights passing the veto
      _ntracks->fill();

      // Do track analysis here
      // Find pttrackmax
      double track_ptmax = 0;
      for (const Particle& p : cfs.particles(Cuts::abseta < 2.4)) track_ptmax = max(track_ptmax, p.pT());
      // Fill track analysis histograms
      for (size_t i = 0; i < _h_tracks->numBins(); ++i) {
        const double binlimitlow_t = _h_tracks->bin(i).xMin();
        const double weightbw_t = _h_tracks->bin(i).xWidth();
        const double xbin_t = _h_tracks->bin(i).xMid();
        if (track_ptmax > binlimitlow_t) _h_tracks -> fill(xbin_t, weightbw_t);
      }

      // Do jet analysis here
      const Jets jetsdeta = apply<FastJets>(event,"Jets").jets(Cuts::pT > 1*GeV && Cuts::pT < 60*GeV && Cuts::abseta < 1.9);
      // Find ptjetmax
      double jet_ptmax = 0;
      for (const Jet& j : jetsdeta) jet_ptmax = max(jet_ptmax, j.pT());
      // Fill jet analysis histograms
      for (size_t i = 0; i < _h_jets->numBins(); ++i) {
        const double binlimitlow_j = _h_jets->bin(i).xMin();
        const double weightbw_j = _h_jets->bin(i).xWidth();
        const double xbin_j = _h_jets->bin(i).xMid();
        if (jet_ptmax > binlimitlow_j) _h_jets -> fill(xbin_j, weightbw_j);
      }
    }


    /// Normalise histograms etc., after the run
    void finalize() {
      const double norm_t0 = _h_tracks->bin(7).height()/2.056170e-03;
      //const double norm_t1 = _h_tracks->bin(7).sumW()/2.056170e-03;
      const double norm_j0 = _h_jets->bin(13).height()/3.575290e-03;
      //const double norm_j1 = _h_jets->bin(13).sumW()/3.575290e-03;
      if (norm_t0 > 0 ) scale(_h_tracks, 1./ norm_t0);
      if (norm_j0 > 0 ) scale(_h_jets, 1./ norm_j0);
    }

    //@}


  private:

    /// @name Histograms
    //@{
    Histo1DPtr _h_tracks, _h_jets;
    CounterPtr _ntracks;
    //@}

  };


  // The hook for the plugin system
  RIVET_DECLARE_PLUGIN(CMS_2015_I1380605);

}