Rivet Analyses Reference

CMS_2010_S8547297

Charged-particle pT and pseudorapidity spectra from pp collisions at 900 and 2360 GeV.
Experiment: CMS (LHC)
Inspire ID: 845323
Status: VALIDATED
Authors:
  • A. Knutsson
References:Beams: p+ p+
Beam energies: (450.0, 450.0); (1180.0, 1180.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. Examples, 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.

Charged particle spectra are measured in proton-proton collisions at center-of-mass energies 900 and 2360 GeV. The spectra are normalized to all non-single-diffractive (NSD) events using corrections for trigger and selection efficiency, acceptance, and branching ratios. There are transverse-momentum (pT) spectra from 0.1 to 2 GeV in bins of pseudorapidity (eta) and pT spectra from 0.1 to 4 GeV for |eta|<2.4. The eta spectra come from the average of three methods and cover |eta|<2.5 and are corrected to include all pT. The data were corrected according to the SD/DD/ND content of the CMS trigger, as predicted by PYTHIA6. The uncertainties connected with correct or incorrect modelling of diffraction were included in the systematic errors. Beam energy must be specified (in GeV) as analysis option "ENERGY" when rivet-merging samples.

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

namespace Rivet {


  /// Charged-particle pT and pseudorapidity spectra from pp collisions at 900 and 2360 GeV
  class CMS_2010_S8547297 : public Analysis {
  public:

    RIVET_DEFAULT_ANALYSIS_CTOR(CMS_2010_S8547297);


    /// @{

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

      if (isCompatibleWithSqrtS(900)) {
        for (int d=1; d<=3; d++) {
          for (int y=1; y<=4; y++) {
            _h_dNch_dpT.push_back(Histo1DPtr());
            book(_h_dNch_dpT.back(), d, 1, y);
          }
        }
        book(_h_dNch_dpT_all ,7, 1, 1);
        book(_h_dNch_dEta ,8, 1, 1);
      } else if (isCompatibleWithSqrtS(2360)) {
        for (int d=4; d<=6; d++) {
          for (int y=1; y<=4; y++) {
            _h_dNch_dpT.push_back(Histo1DPtr());
            book(_h_dNch_dpT.back(), d, 1, y);
          }
        }
        book(_h_dNch_dpT_all ,7, 1, 2);
        book(_h_dNch_dEta ,8, 1, 2);
      }
    }


    void analyze(const Event& event) {
      const double weight = 1.0;

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

      for (const Particle& p : charged.particles()) {
        //selecting only charged hadrons
        if (! PID::isHadron(p.pid())) continue;

        const double pT = p.pT();
        const double eta = p.eta();

        // The data is actually a duplicated folded distribution.  This should mimic it.
        _h_dNch_dEta->fill(eta, 0.5*weight);
        _h_dNch_dEta->fill(-eta, 0.5*weight);
        if (fabs(eta) < 2.4 && pT > 0.1*GeV) {
          if (pT < 4.0*GeV) {
            _h_dNch_dpT_all->fill(pT/GeV, weight/(pT/GeV));
            if (pT < 2.0*GeV) {
              int ietabin = int(fabs(eta)/0.2);
              _h_dNch_dpT[ietabin]->fill(pT/GeV, weight);
            }
          }
        }
      }
    }


    void finalize() {
      const double normfac = 1.0/sumOfWeights(); // Normalizing to unit eta is automatic
      // The pT distributions in bins of eta must be normalized to unit eta.  This is a factor of 2
      // for the |eta| times 0.2 (eta range).
      // The pT distributions over all eta are normalized to unit eta (2.0*2.4) and by 1/2*pi*pT.
      // The 1/pT part is taken care of in the filling.  The 1/2pi is taken care of here.
      const double normpT = normfac/(2.0*0.2);
      const double normpTall = normfac/(2.0*M_PI*2.0*2.4);

      for (size_t ietabin=0; ietabin < _h_dNch_dpT.size(); ietabin++){
        scale(_h_dNch_dpT[ietabin], normpT);
      }
      scale(_h_dNch_dpT_all, normpTall);
      scale(_h_dNch_dEta, normfac);
    }

    /// @}


  private:

    /// @{
    std::vector<Histo1DPtr> _h_dNch_dpT;
    Histo1DPtr _h_dNch_dpT_all;
    Histo1DPtr _h_dNch_dEta;
    /// @}

  };



  RIVET_DECLARE_ALIASED_PLUGIN(CMS_2010_S8547297, CMS_2010_I845323);

}