rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

CMS_2015_I1385107

Underlying event measurement with leading jets at $\sqrt{s} = 2.76$ \text{TeV}
Experiment: CMS (LHC)
Inspire ID: 1385107
Status: VALIDATED
Authors:
  • Wei Yang Wang
  • Xavier Janssen
No references listed
Beams: p+ p+
Beam energies: (1380.0, 1380.0) GeV
Run details:
  • Requires inclusive inelastic events (non-diffractive and inelastic diffractive). The profile plots require large statistics.

A measurement of the underlying event (UE) activity in proton-proton collisions is performed using events with charged-particle jets produced in the central pseudorapidity region ($|\eta|^\text{jet} < 2$) and with transverse momentum $1 \leq p_T^\text{jet} < 100$ \text{GeV}. The analysis uses a data sample collected at a centre-of-mass energy of 2.76 \text{TeV} with the CMS experiment at the LHC. The UE activity is measured as a function of $p_T^\text{jet}$ in terms of the average multiplicity and scalar sum of transverse momenta of charged particles, with $|\eta| < 2$ and $p_T > 0.5$ \text{GeV}, in the azimuthal region transverse to the highest-$p_T$ jet direction. By further dividing the transverse region into two regions of smaller and larger activity, various components of the UE activity are separated.

Source code: CMS_2015_I1385107.cc
  1// -*- C++ -*-
  2#include "Rivet/Analysis.hh"
  3#include "Rivet/Projections/FinalState.hh"
  4#include "Rivet/Projections/ChargedFinalState.hh"
  5#include "Rivet/Projections/FastJets.hh"
  6
  7namespace Rivet {
  8
  9
 10  /// CMS UE charged particles vs. leading jet at 2.76 TeV
 11  class CMS_2015_I1385107 : public Analysis {
 12  public:
 13    /// Constructor
 14    CMS_2015_I1385107() : Analysis("CMS_2015_I1385107"),
 15                          ETACUT(2.0),
 16                          AREATOT(2*ETACUT * 2*M_PI),
 17                          AREA3(AREATOT / 3.),
 18                          AREA6(AREATOT / 6.)
 19    {   }
 20
 21
 22    /// Book histograms and initialise projections before the run
 23    void init() {
 24
 25      const ChargedFinalState cfs(Cuts::abseta < 2 && Cuts::pT > 500*MeV);
 26      declare(cfs, "CFS");
 27
 28      const ChargedFinalState cfsforjet(Cuts::abseta < 2.5 && Cuts::pT > 500*MeV);
 29      const FastJets jetpro(cfsforjet, JetAlg::SISCONE, 0.5);
 30      declare(jetpro, "Jets");
 31
 32      book(_h_Nch_TransAVE_vs_pT ,1, 1, 1); // Nch vs. pT_max      (TransAVE)
 33      book(_h_Sum_TransAVE_vs_pT ,2, 1, 1); // sum(pT) vs. pT_max  (TransAVE)
 34      book(_h_Nch_TransMAX_vs_pT ,3, 1, 1); // Nch vs. pT_max      (TransMAX)
 35      book(_h_Sum_TransMAX_vs_pT ,4, 1, 1); // sum(pT) vs. pT_max  (TransMAX)
 36      book(_h_Nch_TransMIN_vs_pT ,5, 1, 1); // Nch vs. pT_max      (TransMIN)
 37      book(_h_Sum_TransMIN_vs_pT ,6, 1, 1); // sum(pT) vs. pT_max  (TransMIN)
 38      book(_h_Nch_TransDIF_vs_pT ,7, 1, 1); // Nch vs. pT_max      (TransDIF)
 39      book(_h_Sum_TransDIF_vs_pT ,8, 1, 1); // sum(pT) vs. pT_max  (TransDIF)
 40    }
 41
 42
 43    /// Local definition of a signed dphi, for use in differentating L and R trans regions
 44    double signedDeltaPhi(double jetphi, double partphi) {
 45      double delta = partphi - jetphi;
 46      while (delta <= -PI) delta += 2 * PI;
 47      while (delta > PI) delta -= 2 * PI;
 48      return delta;
 49    }
 50
 51    /// Perform the per-event analysis
 52    void analyze(const Event& event) {
 53
 54      // Find the lead jet, applying a restriction that the jets must be within |eta| < 2.
 55      FourMomentum p_lead;
 56      for (const Jet& j : apply<FastJets>(event, "Jets").jetsByPt(Cuts::abseta < 2.0 && Cuts::pT > 1*GeV)) {
 57        p_lead = j.momentum();
 58        break;
 59      }
 60      if (p_lead.isZero()) vetoEvent;
 61      const double phi_lead = p_lead.phi();
 62      const double pT_lead  = p_lead.pT();
 63
 64      // Loop on charged particles and separate Left and Right transverse regions
 65      Particles particles = apply<ChargedFinalState>(event, "CFS").particlesByPt();
 66      int nch_TransLeft = 0, nch_TransRight = 0;
 67      double ptSum_TransLeft = 0., ptSum_TransRight = 0.;
 68      for (const Particle& p : particles) {
 69        const double dphi = signedDeltaPhi(phi_lead, p.momentum().phi());
 70        if (!inRange(fabs(dphi), PI/3, 2*PI/3.)) continue; //< only fill trans regions
 71        if (dphi < 0) {  // Transverse Right region
 72          nch_TransRight += 1;
 73          ptSum_TransRight += p.pT() / GeV;
 74        } else if (dphi > 0) {  // Transverse Left region
 75          nch_TransLeft += 1;
 76          ptSum_TransLeft += p.pT() / GeV;
 77        }
 78      }
 79
 80      // Translate to min and max (+sum and diff) Transverse regions
 81      const int nch_TransMIN = std::min(nch_TransLeft, nch_TransRight);
 82      const int nch_TransMAX = std::max(nch_TransLeft, nch_TransRight);
 83      const int nch_TransSUM = nch_TransMAX + nch_TransMIN;
 84      const int nch_TransDIF = nch_TransMAX - nch_TransMIN;
 85      //
 86      const double ptSum_TransMIN = std::min(ptSum_TransLeft, ptSum_TransRight);
 87      const double ptSum_TransMAX = std::max(ptSum_TransLeft, ptSum_TransRight);
 88      const double ptSum_TransSUM = ptSum_TransMAX + ptSum_TransMIN;
 89      const double ptSum_TransDIF = ptSum_TransMAX - ptSum_TransMIN;
 90
 91      // Fill profiles
 92      _h_Nch_TransMIN_vs_pT->fill(pT_lead/GeV, 1/AREA6 * nch_TransMIN);
 93      _h_Sum_TransMIN_vs_pT->fill(pT_lead/GeV, 1/AREA6 * ptSum_TransMIN);
 94      //
 95      _h_Nch_TransMAX_vs_pT->fill(pT_lead/GeV, 1/AREA6 * nch_TransMAX);
 96      _h_Sum_TransMAX_vs_pT->fill(pT_lead/GeV, 1/AREA6 * ptSum_TransMAX);
 97      //
 98      _h_Nch_TransAVE_vs_pT->fill(pT_lead/GeV, 1/AREA3 * nch_TransSUM);
 99      _h_Sum_TransAVE_vs_pT->fill(pT_lead/GeV, 1/AREA3 * ptSum_TransSUM);
100      //
101      _h_Nch_TransDIF_vs_pT->fill(pT_lead/GeV, 1/AREA6 * nch_TransDIF);
102      _h_Sum_TransDIF_vs_pT->fill(pT_lead/GeV, 1/AREA6 * ptSum_TransDIF);
103    }
104
105
106  private:
107
108    // Data members like post-cuts event weight counters go here
109    const double ETACUT, AREATOT, AREA3, AREA6;
110
111    /// Histograms
112    Profile1DPtr _h_Nch_TransAVE_vs_pT, _h_Sum_TransAVE_vs_pT;
113    Profile1DPtr _h_Nch_TransDIF_vs_pT, _h_Sum_TransDIF_vs_pT;
114    Profile1DPtr _h_Nch_TransMIN_vs_pT, _h_Sum_TransMIN_vs_pT;
115    Profile1DPtr _h_Nch_TransMAX_vs_pT, _h_Sum_TransMAX_vs_pT;
116
117  };
118
119
120  RIVET_DECLARE_PLUGIN(CMS_2015_I1385107);
121
122}