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:
  • CMS-QCD-11-012
  • CERN-PH-EP-2012-085
  • arXiv: 1204.1411
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
  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
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
// -*- C++ -*-
#include "Rivet/Analysis.hh"
#include "Rivet/Projections/FinalState.hh"
#include "Rivet/Projections/ZFinder.hh"
#include "Rivet/Projections/FinalState.hh"
#include "Rivet/Projections/ChargedFinalState.hh"
#include "Rivet/Projections/VetoedFinalState.hh"

namespace Rivet {


  /// Underlying event activity in the Drell-Yan process at 7 TeV
  class CMS_2012_I1107658 : public Analysis {
  public:

    /// Constructor
    CMS_2012_I1107658()
      : Analysis("CMS_2012_I1107658")
    {   }


    /// Initialization
    void init() {

      /// @note Using a bare muon Z (but with a clustering radius!?)
      Cut cut = Cuts::abseta < 2.4 && Cuts::pT > 20*GeV;
      ZFinder zfinder(FinalState(), cut, PID::MUON, 4*GeV, 140*GeV, 0.2, ZFinder::ClusterPhotons::NONE);
      declare(zfinder, "ZFinder");

      ChargedFinalState cfs((Cuts::etaIn(-2, 2) && Cuts::pT >=  500*MeV));
      VetoedFinalState nonmuons(cfs);
      nonmuons.addVetoPairId(PID::MUON);
      declare(nonmuons, "nonmuons");

      book(_h_Nchg_towards_pTmumu                 ,1, 1, 1);
      book(_h_Nchg_transverse_pTmumu              ,2, 1, 1);
      book(_h_Nchg_away_pTmumu                    ,3, 1, 1);
      book(_h_pTsum_towards_pTmumu                ,4, 1, 1);
      book(_h_pTsum_transverse_pTmumu             ,5, 1, 1);
      book(_h_pTsum_away_pTmumu                   ,6, 1, 1);
      book(_h_avgpT_towards_pTmumu                ,7, 1, 1);
      book(_h_avgpT_transverse_pTmumu             ,8, 1, 1);
      book(_h_avgpT_away_pTmumu                   ,9, 1, 1);
      book(_h_Nchg_towards_plus_transverse_Mmumu  ,10, 1, 1);
      book(_h_pTsum_towards_plus_transverse_Mmumu ,11, 1, 1);
      book(_h_avgpT_towards_plus_transverse_Mmumu ,12, 1, 1);
      book(_h_Nchg_towards_zmass_81_101           ,13, 1, 1);
      book(_h_Nchg_transverse_zmass_81_101        ,14, 1, 1);
      book(_h_Nchg_away_zmass_81_101              ,15, 1, 1);
      book(_h_pT_towards_zmass_81_101             ,16, 1, 1);
      book(_h_pT_transverse_zmass_81_101          ,17, 1, 1);
      book(_h_pT_away_zmass_81_101                ,18, 1, 1);
      book(_h_Nchg_transverse_zpt_5               ,19, 1, 1);
      book(_h_pT_transverse_zpt_5                 ,20, 1, 1);
    }


    /// Perform the per-event analysis
    void analyze(const Event& event) {
      const double weight = 1.0;
      const ZFinder& zfinder = apply<ZFinder>(event, "ZFinder");

      if (zfinder.bosons().size() != 1) vetoEvent;

      double Zpt = zfinder.bosons()[0].pT()/GeV;
      double Zphi = zfinder.bosons()[0].phi();
      double Zmass = zfinder.bosons()[0].mass()/GeV;

      Particles particles = apply<VetoedFinalState>(event, "nonmuons").particles();

      int nTowards = 0;
      int nTransverse = 0;
      int nAway = 0;
      double ptSumTowards = 0;
      double ptSumTransverse = 0;
      double ptSumAway = 0;

      for (const Particle& p : particles) {
        double dphi = fabs(deltaPhi(Zphi, p.phi()));
        double pT = p.pT();

        if ( dphi < M_PI/3 ) {
          nTowards++;
          ptSumTowards += pT;
          if (Zmass > 81. && Zmass < 101.) _h_pT_towards_zmass_81_101->fill(pT, weight);
        } else if ( dphi < 2.*M_PI/3 ) {
          nTransverse++;
          ptSumTransverse += pT;
          if (Zmass > 81. && Zmass < 101.) _h_pT_transverse_zmass_81_101->fill(pT, weight);
          if (Zpt < 5.) _h_pT_transverse_zpt_5->fill(pT, weight);
        } else {
          nAway++;
          ptSumAway += pT;
          if (Zmass > 81. && Zmass < 101.) _h_pT_away_zmass_81_101->fill(pT, weight);
        }

      } // Loop over particles


      const double area = 8./3.*M_PI;
      if (Zmass > 81. && Zmass < 101.) {
        _h_Nchg_towards_pTmumu->         fill(Zpt, 1./area * nTowards, weight);
        _h_Nchg_transverse_pTmumu->      fill(Zpt, 1./area * nTransverse, weight);
        _h_Nchg_away_pTmumu->            fill(Zpt, 1./area * nAway, weight);
        _h_pTsum_towards_pTmumu->        fill(Zpt, 1./area * ptSumTowards, weight);
        _h_pTsum_transverse_pTmumu->     fill(Zpt, 1./area * ptSumTransverse, weight);
        _h_pTsum_away_pTmumu->           fill(Zpt, 1./area * ptSumAway, weight);
        if (nTowards > 0)    _h_avgpT_towards_pTmumu->    fill(Zpt, ptSumTowards/nTowards, weight);
        if (nTransverse > 0) _h_avgpT_transverse_pTmumu-> fill(Zpt, ptSumTransverse/nTransverse, weight);
        if (nAway > 0)       _h_avgpT_away_pTmumu->       fill(Zpt, ptSumAway/nAway, weight);
        _h_Nchg_towards_zmass_81_101->   fill(nTowards, weight);
        _h_Nchg_transverse_zmass_81_101->fill(nTransverse, weight);
        _h_Nchg_away_zmass_81_101->      fill(nAway, weight);
      }

      if (Zpt < 5.) {
        _h_Nchg_towards_plus_transverse_Mmumu->fill(Zmass, (nTowards + nTransverse)/(2.*area), weight);
        _h_pTsum_towards_plus_transverse_Mmumu->fill(Zmass, (ptSumTowards + ptSumTransverse)/(2.*area), weight);
        if ((nTowards + nTransverse) > 0) _h_avgpT_towards_plus_transverse_Mmumu->fill(Zmass, (ptSumTowards + ptSumTransverse)/(nTowards + nTransverse), weight);
        _h_Nchg_transverse_zpt_5->fill(nTransverse, weight);
      }

    }


    /// Normalise histograms etc., after the run
    void finalize() {
      scale(_h_pT_towards_zmass_81_101,    safediv(1, _h_Nchg_towards_zmass_81_101->integral(), 0));
      scale(_h_pT_transverse_zmass_81_101, safediv(1, _h_Nchg_transverse_zmass_81_101->integral(), 0));
      scale(_h_pT_away_zmass_81_101,       safediv(1, _h_Nchg_away_zmass_81_101->integral(), 0));
      scale(_h_pT_transverse_zpt_5,        safediv(1, _h_Nchg_transverse_zpt_5->integral(), 0));
      normalize(_h_Nchg_towards_zmass_81_101);
      normalize(_h_Nchg_transverse_zmass_81_101);
      normalize(_h_Nchg_away_zmass_81_101);
      normalize(_h_Nchg_transverse_zpt_5);
    }


  private:


    /// @name Histogram objects
    //@{
    Profile1DPtr _h_Nchg_towards_pTmumu;
    Profile1DPtr _h_Nchg_transverse_pTmumu;
    Profile1DPtr _h_Nchg_away_pTmumu;
    Profile1DPtr _h_pTsum_towards_pTmumu;
    Profile1DPtr _h_pTsum_transverse_pTmumu;
    Profile1DPtr _h_pTsum_away_pTmumu;
    Profile1DPtr _h_avgpT_towards_pTmumu;
    Profile1DPtr _h_avgpT_transverse_pTmumu;
    Profile1DPtr _h_avgpT_away_pTmumu;
    Profile1DPtr _h_Nchg_towards_plus_transverse_Mmumu;
    Profile1DPtr _h_pTsum_towards_plus_transverse_Mmumu;
    Profile1DPtr _h_avgpT_towards_plus_transverse_Mmumu;
    Histo1DPtr _h_Nchg_towards_zmass_81_101;
    Histo1DPtr _h_Nchg_transverse_zmass_81_101;
    Histo1DPtr _h_Nchg_away_zmass_81_101;
    Histo1DPtr _h_pT_towards_zmass_81_101;
    Histo1DPtr _h_pT_transverse_zmass_81_101;
    Histo1DPtr _h_pT_away_zmass_81_101;
    Histo1DPtr _h_Nchg_transverse_zpt_5;
    Histo1DPtr _h_pT_transverse_zpt_5;
    //@}

  };


  // Hook for the plugin system
  RIVET_DECLARE_PLUGIN(CMS_2012_I1107658);

}