Rivet Analyses Reference

CMS_2011_S8884919

Measurement of the NSD charged particle multiplicity at $\sqrt{s} = 0.9$, 2.36, and 7 TeV
Experiment: CMS (LHC)
Inspire ID: 879315
Status: VALIDATED
Authors:
  • Romain Rougny
References:Beams: p+ p+
Beam energies: (450.0, 450.0); (1180.0, 1180.0); (3500.0, 3500.0) GeV
Run details:
  • Non-single-diffractive (NSD) events only. Should include double-diffractive (DD) events and non-diffractive (ND) events but NOT single-diffractive (SD) events. For example, in Pythia6 the SD processes to be turned off are 92 and 93 and in Pythia8 the SD processes are 103 and 104 (also called SoftQCD:singleDiffractive). Beam energy must be specified as analysis option "ENERGY" when rivet-merge'ing samples.

Measurements of primary charged hadron multiplicity distributions are presented for non-single-diffractive events in proton-proton collisions at centre-of-mass energies of $\sqrt{s} = 0.9$, 2.36, and 7 TeV, in five pseudorapidity ranges from $|\eta| < 0.5$ to $|\eta| < 2.4$. The data were collected with the minimum-bias trigger of the CMS experiment during the LHC commissioning runs in 2009 and the 7 TeV run in 2010. The average transverse momentum as a function of the multiplicity is also presented. The measurement of higher-order moments of the multiplicity distribution confirms the violation of Koba-Nielsen-Olesen scaling that has been observed at lower energies. Beam energy must be specified (in GeV) as analysis option "ENERGY" when rivet-merging samples.

Source code: CMS_2011_S8884919.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
// -*- C++ -*-
#include "Rivet/Analysis.hh"
#include "Rivet/Projections/ChargedFinalState.hh"
#include "Rivet/Projections/Beam.hh"
using namespace std;

namespace Rivet {


  /// Measurement of the NSD charged particle multiplicity at 0.9, 2.36, and 7 TeV
  class CMS_2011_S8884919 : public Analysis {
  public:

    RIVET_DEFAULT_ANALYSIS_CTOR(CMS_2011_S8884919);


    void init() {
      ChargedFinalState cfs((Cuts::etaIn(-2.4, 2.4)));
      declare(cfs, "CFS");

      // eta bins
      _etabins.push_back(0.5);
      _etabins.push_back(1.0);
      _etabins.push_back(1.5);
      _etabins.push_back(2.0);
      _etabins.push_back(2.4) ;

      if (isCompatibleWithSqrtS(900)) {
        for (size_t ietabin=0; ietabin < _etabins.size(); ietabin++) {
          _h_dNch_dn.push_back( Histo1DPtr() );
          book( _h_dNch_dn.back(), 2 + ietabin, 1, 1);
        }
        book(_h_dNch_dn_pt500_eta24 ,20, 1, 1);
        book(_h_dmpt_dNch_eta24 ,23, 1, 1);
      }

      if (isCompatibleWithSqrtS(2360)) {
        for (size_t ietabin=0; ietabin < _etabins.size(); ietabin++) {
          _h_dNch_dn.push_back( Histo1DPtr() );
          book(_h_dNch_dn.back(), 7 + ietabin, 1, 1);
        }
        book(_h_dNch_dn_pt500_eta24 ,21, 1, 1);
        book(_h_dmpt_dNch_eta24 ,24, 1, 1);
      }

      if (isCompatibleWithSqrtS(7000)) {
        for (size_t ietabin=0; ietabin < _etabins.size(); ietabin++) {
          _h_dNch_dn.push_back( Histo1DPtr() );
          book(_h_dNch_dn.back(), 12 + ietabin, 1, 1);
        }
        book(_h_dNch_dn_pt500_eta24 ,22, 1, 1);
        book(_h_dmpt_dNch_eta24 ,25, 1, 1);
      }
    }


    void analyze(const Event& event) {

      // Get the charged particles
      const ChargedFinalState& charged = apply<ChargedFinalState>(event, "CFS");

      // Resetting the multiplicity for the event to 0;
      vector<int> _nch_in_Evt;
      vector<int> _nch_in_Evt_pt500;
      _nch_in_Evt.assign(_etabins.size(), 0);
      _nch_in_Evt_pt500.assign(_etabins.size(), 0);
      double sumpt = 0;

      // Loop over particles in event
      for (const Particle& p : charged.particles()) {
        // Selecting only charged hadrons
        if (! PID::isHadron(p.pid())) continue;

        double pT = p.pT();
        double eta = p.eta();
        sumpt += pT;
        for (size_t ietabin = _etabins.size(); ietabin > 0; --ietabin) {
          if (fabs(eta) > _etabins[ietabin-1]) break;
          ++_nch_in_Evt[ietabin-1];
          if (pT > 0.5/GeV) ++_nch_in_Evt_pt500[ietabin-1];
        }
      }

      // Filling multiplicity-dependent histogramms
      for (size_t ietabin = 0; ietabin < _etabins.size(); ietabin++) {
        _h_dNch_dn[ietabin]->fill(_nch_in_Evt[ietabin]);
      }

      // Do only if eta bins are the needed ones
      if (_etabins[4] == 2.4 && _etabins[0] == 0.5) {
        if (_nch_in_Evt[4] != 0) {
          _h_dmpt_dNch_eta24->fill(_nch_in_Evt[4], sumpt/GeV / _nch_in_Evt[4]);
        }
        _h_dNch_dn_pt500_eta24->fill(_nch_in_Evt_pt500[4]);
      } else {
        MSG_WARNING("You changed the number of eta bins, but forgot to propagate it everywhere !!");
      }
    }


    void finalize() {
      for (size_t ietabin = 0; ietabin < _etabins.size(); ietabin++){
        normalize(_h_dNch_dn[ietabin]);
      }
      normalize(_h_dNch_dn_pt500_eta24);
    }


  private:

    /// @{
    vector<Histo1DPtr> _h_dNch_dn;
    Histo1DPtr _h_dNch_dn_pt500_eta24;
    Profile1DPtr _h_dmpt_dNch_eta24;
    /// @}

    vector<double> _etabins;

  };



  RIVET_DECLARE_ALIASED_PLUGIN(CMS_2011_S8884919, CMS_2011_I879315);

}