rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

CMS_2012_I1107658

Measurement of the underlying event activity in the Drell-Yan process at a centre-of-mass energy of 7 TeV
Experiment: CMS (LHC)
Inspire ID: 1107658
Status: VALIDATED
Authors:
  • Sunil Bansal
References: Beams: p+ p+
Beam energies: (3500.0, 3500.0) GeV
Run details:
  • Drell-Yan events with $Z/\gamma^* \to \mu\mu$. $m(\mu,\mu) > 20$ GeV

A measurement of the underlying event activity using Drell-Yan events using muonic final state. The production of charged particles with pseudorapidity $|\eta| < 2$ and transverse momentum $p_\perp > 0.5\,\text{GeV}/c$ is studied in towards, transverse and away region w.r.t. to the direction of di-muon system. The UE activity is measured in terms of of a particle density and an energy density. The particle density is computed as the average number of primary charged particles per unit pseudorapidity and per unit azimuth. The energy density is expressed in terms of the average of the scalar sum of the transverse momenta of primary charged particles per unit pseudorapidity and azimuth. The ratio of the energy and particle density is also reported in 3 regions. UE activity is studied as a function of invariant mass of muon pair ($M_{\mu\mu}$) by limiting the ISR contribution by requiring transverse momentum of muon pair $p_\perp(\mu\mu) < 5\,\text{GeV}/c$. The $p_\perp(\mu\mu)$ dependence is studied for the events having $M_{\mu\mu}$ in window of 81--101 GeV/$c$. The normalized charged particle multiplicity and $p_\perp$ spectrum of the charged particles in three regions also been reported for events having $M_{\mu\mu}$ in window of 81--101 GeV/$c$. Multiplicity and $p_\perp$ spectra in the transverse region are also reported, for events having $p_\perp(\mu\mu) < 5\,\text{GeV}/c$.

Source code: CMS_2012_I1107658.cc
  1// -*- C++ -*-
  2#include "Rivet/Analysis.hh"
  3#include "Rivet/Projections/FinalState.hh"
  4#include "Rivet/Projections/DileptonFinder.hh"
  5#include "Rivet/Projections/FinalState.hh"
  6#include "Rivet/Projections/ChargedFinalState.hh"
  7#include "Rivet/Projections/VetoedFinalState.hh"
  8
  9namespace Rivet {
 10
 11
 12  /// Underlying event activity in the Drell-Yan process at 7 TeV
 13  class CMS_2012_I1107658 : public Analysis {
 14  public:
 15
 16    /// Constructor
 17    RIVET_DEFAULT_ANALYSIS_CTOR(CMS_2012_I1107658);
 18
 19
 20    /// Initialization
 21    void init() {
 22
 23      /// @note Using a bare muon Z (but with a clustering radius!?)
 24      Cut cut = Cuts::abseta < 2.4 && Cuts::pT > 20*GeV;
 25      DileptonFinder zfinder(91.2*GeV, 0.2, cut && Cuts::abspid == PID::MUON, Cuts::massIn(4*GeV, 140*GeV));
 26      declare(zfinder, "DileptonFinder");
 27
 28      ChargedFinalState nonmuons(Cuts::abseta < 2 && Cuts::pT > 500*MeV && Cuts::abspid != PID::MUON);
 29      declare(nonmuons, "nonmuons");
 30
 31      book(_h_Nchg_towards_pTmumu                 ,1, 1, 1);
 32      book(_h_Nchg_transverse_pTmumu              ,2, 1, 1);
 33      book(_h_Nchg_away_pTmumu                    ,3, 1, 1);
 34      book(_h_pTsum_towards_pTmumu                ,4, 1, 1);
 35      book(_h_pTsum_transverse_pTmumu             ,5, 1, 1);
 36      book(_h_pTsum_away_pTmumu                   ,6, 1, 1);
 37      book(_h_avgpT_towards_pTmumu                ,7, 1, 1);
 38      book(_h_avgpT_transverse_pTmumu             ,8, 1, 1);
 39      book(_h_avgpT_away_pTmumu                   ,9, 1, 1);
 40      book(_h_Nchg_towards_plus_transverse_Mmumu  ,10, 1, 1);
 41      book(_h_pTsum_towards_plus_transverse_Mmumu ,11, 1, 1);
 42      book(_h_avgpT_towards_plus_transverse_Mmumu ,12, 1, 1);
 43      book(_h_Nchg_towards_zmass_81_101           ,13, 1, 1);
 44      book(_h_Nchg_transverse_zmass_81_101        ,14, 1, 1);
 45      book(_h_Nchg_away_zmass_81_101              ,15, 1, 1);
 46      book(_h_pT_towards_zmass_81_101             ,16, 1, 1);
 47      book(_h_pT_transverse_zmass_81_101          ,17, 1, 1);
 48      book(_h_pT_away_zmass_81_101                ,18, 1, 1);
 49      book(_h_Nchg_transverse_zpt_5               ,19, 1, 1);
 50      book(_h_pT_transverse_zpt_5                 ,20, 1, 1);
 51    }
 52
 53
 54    /// Perform the per-event analysis
 55    void analyze(const Event& event) {
 56      const DileptonFinder& zfinder = apply<DileptonFinder>(event, "DileptonFinder");
 57
 58      if (zfinder.bosons().size() != 1) vetoEvent;
 59
 60      double Zpt = zfinder.bosons()[0].pT()/GeV;
 61      double Zphi = zfinder.bosons()[0].phi();
 62      double Zmass = zfinder.bosons()[0].mass()/GeV;
 63
 64      Particles particles = apply<VetoedFinalState>(event, "nonmuons").particles();
 65
 66      int nTowards = 0;
 67      int nTransverse = 0;
 68      int nAway = 0;
 69      double ptSumTowards = 0;
 70      double ptSumTransverse = 0;
 71      double ptSumAway = 0;
 72
 73      for (const Particle& p : particles) {
 74        double dphi = fabs(deltaPhi(Zphi, p.phi()));
 75        double pT = p.pT();
 76
 77        if ( dphi < M_PI/3 ) {
 78          nTowards++;
 79          ptSumTowards += pT;
 80          if (Zmass > 81. && Zmass < 101.) _h_pT_towards_zmass_81_101->fill(pT);
 81        } else if ( dphi < 2.*M_PI/3 ) {
 82          nTransverse++;
 83          ptSumTransverse += pT;
 84          if (Zmass > 81. && Zmass < 101.) _h_pT_transverse_zmass_81_101->fill(pT);
 85          if (Zpt < 5.) _h_pT_transverse_zpt_5->fill(pT);
 86        } else {
 87          nAway++;
 88          ptSumAway += pT;
 89          if (Zmass > 81. && Zmass < 101.) _h_pT_away_zmass_81_101->fill(pT);
 90        }
 91
 92      } // Loop over particles
 93
 94
 95      const double area = 8./3.*M_PI;
 96      if (Zmass > 81. && Zmass < 101.) {
 97        _h_Nchg_towards_pTmumu->         fill(Zpt, 1./area * nTowards);
 98        _h_Nchg_transverse_pTmumu->      fill(Zpt, 1./area * nTransverse);
 99        _h_Nchg_away_pTmumu->            fill(Zpt, 1./area * nAway);
100        _h_pTsum_towards_pTmumu->        fill(Zpt, 1./area * ptSumTowards);
101        _h_pTsum_transverse_pTmumu->     fill(Zpt, 1./area * ptSumTransverse);
102        _h_pTsum_away_pTmumu->           fill(Zpt, 1./area * ptSumAway);
103        if (nTowards > 0)    _h_avgpT_towards_pTmumu->    fill(Zpt, ptSumTowards/nTowards);
104        if (nTransverse > 0) _h_avgpT_transverse_pTmumu-> fill(Zpt, ptSumTransverse/nTransverse);
105        if (nAway > 0)       _h_avgpT_away_pTmumu->       fill(Zpt, ptSumAway/nAway);
106        _h_Nchg_towards_zmass_81_101->   fill(nTowards);
107        _h_Nchg_transverse_zmass_81_101->fill(nTransverse);
108        _h_Nchg_away_zmass_81_101->      fill(nAway);
109      }
110
111      if (Zpt < 5.) {
112        _h_Nchg_towards_plus_transverse_Mmumu->fill(Zmass, (nTowards + nTransverse)/(2.*area));
113        _h_pTsum_towards_plus_transverse_Mmumu->fill(Zmass, (ptSumTowards + ptSumTransverse)/(2.*area));
114        if ((nTowards + nTransverse) > 0) _h_avgpT_towards_plus_transverse_Mmumu->fill(Zmass, (ptSumTowards + ptSumTransverse)/(nTowards + nTransverse));
115        _h_Nchg_transverse_zpt_5->fill(nTransverse);
116      }
117
118    }
119
120
121    /// Normalise histograms etc., after the run
122    void finalize() {
123      scale(_h_pT_towards_zmass_81_101,    safediv(1, _h_Nchg_towards_zmass_81_101->integral(), 0));
124      scale(_h_pT_transverse_zmass_81_101, safediv(1, _h_Nchg_transverse_zmass_81_101->integral(), 0));
125      scale(_h_pT_away_zmass_81_101,       safediv(1, _h_Nchg_away_zmass_81_101->integral(), 0));
126      scale(_h_pT_transverse_zpt_5,        safediv(1, _h_Nchg_transverse_zpt_5->integral(), 0));
127      normalize(_h_Nchg_towards_zmass_81_101);
128      normalize(_h_Nchg_transverse_zmass_81_101);
129      normalize(_h_Nchg_away_zmass_81_101);
130      normalize(_h_Nchg_transverse_zpt_5);
131    }
132
133
134  private:
135
136
137    /// @name Histogram objects
138    /// @{
139    Profile1DPtr _h_Nchg_towards_pTmumu;
140    Profile1DPtr _h_Nchg_transverse_pTmumu;
141    Profile1DPtr _h_Nchg_away_pTmumu;
142    Profile1DPtr _h_pTsum_towards_pTmumu;
143    Profile1DPtr _h_pTsum_transverse_pTmumu;
144    Profile1DPtr _h_pTsum_away_pTmumu;
145    Profile1DPtr _h_avgpT_towards_pTmumu;
146    Profile1DPtr _h_avgpT_transverse_pTmumu;
147    Profile1DPtr _h_avgpT_away_pTmumu;
148    Profile1DPtr _h_Nchg_towards_plus_transverse_Mmumu;
149    Profile1DPtr _h_pTsum_towards_plus_transverse_Mmumu;
150    Profile1DPtr _h_avgpT_towards_plus_transverse_Mmumu;
151    Histo1DPtr _h_Nchg_towards_zmass_81_101;
152    Histo1DPtr _h_Nchg_transverse_zmass_81_101;
153    Histo1DPtr _h_Nchg_away_zmass_81_101;
154    Histo1DPtr _h_pT_towards_zmass_81_101;
155    Histo1DPtr _h_pT_transverse_zmass_81_101;
156    Histo1DPtr _h_pT_away_zmass_81_101;
157    Histo1DPtr _h_Nchg_transverse_zpt_5;
158    Histo1DPtr _h_pT_transverse_zpt_5;
159    /// @}
160
161  };
162
163
164  RIVET_DECLARE_PLUGIN(CMS_2012_I1107658);
165
166}