rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

CMS_2010_I845323

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).

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.

Source code: CMS_2010_I845323.cc
  1// -*- C++ -*-
  2#include "Rivet/Analysis.hh"
  3#include "Rivet/Projections/ChargedFinalState.hh"
  4
  5namespace Rivet {
  6
  7
  8  /// Charged-particle pT and pseudorapidity spectra from pp collisions at 900 and 2360 GeV
  9  class CMS_2010_I845323 : public Analysis {
 10  public:
 11
 12    RIVET_DEFAULT_ANALYSIS_CTOR(CMS_2010_I845323);
 13
 14
 15    /// @{
 16
 17    void init() {
 18      ChargedFinalState cfs(Cuts::abseta < 2.5);
 19      declare(cfs, "CFS");
 20
 21      for (double eVal : allowedEnergies()) {
 22        const string en = toString(int(eVal));
 23        if (isCompatibleWithSqrtS(eVal))  _sqs = en;
 24        bool offset(en == "2360"s);
 25
 26        for (int d=0; d<3; ++d) {
 27          for (int y=0; y<4; ++y) {
 28            size_t bin = 4*d+y;
 29            book(_h[en+"dNch_dpT"+toString(bin)], d+(offset? 4 : 1), 1, y+1);
 30          }
 31        }
 32        book(_h[en+"dNch_dpT_all"], 7, 1, 1+offset);
 33        book(_h[en+"dNch_dEta"],    8, 1, 1+offset);
 34      }
 35      if (_sqs == "" && !merging()) {
 36        throw BeamError("Invalid beam energy for " + name() + "\n");
 37      }
 38    }
 39
 40
 41    void analyze(const Event& event) {
 42
 43      //charged particles
 44      const ChargedFinalState& charged = apply<ChargedFinalState>(event, "CFS");
 45
 46      for (const Particle& p : charged.particles()) {
 47        //selecting only charged hadrons
 48        if (! PID::isHadron(p.pid())) continue;
 49
 50        const double pT = p.pT();
 51        const double eta = p.eta();
 52
 53        // The data is actually a duplicated folded distribution. This should mimic it.
 54        _h[_sqs+"dNch_dEta"]->fill(eta, 0.5);
 55        _h[_sqs+"dNch_dEta"]->fill(-eta, 0.5);
 56        if (fabs(eta) < 2.4 && pT > 0.1*GeV) {
 57          if (pT < 4.0*GeV) {
 58            _h[_sqs+"dNch_dpT_all"]->fill(pT/GeV, 1.0/(pT/GeV));
 59            if (pT < 2.0*GeV) {
 60              const string suff = toString(int(fabs(eta)/0.2));
 61              _h[_sqs+"dNch_dpT"+suff]->fill(pT/GeV);
 62            }
 63          }
 64        }
 65      }
 66    }
 67
 68
 69    void finalize() {
 70      const double normfac = 1.0/sumOfWeights(); // Normalizing to unit eta is automatic
 71      // The pT distributions in bins of eta must be normalized to unit eta.  This is a factor of 2
 72      // for the |eta| times 0.2 (eta range).
 73      // The pT distributions over all eta are normalized to unit eta (2.0*2.4) and by 1/2*pi*pT.
 74      // The 1/pT part is taken care of in the filling.  The 1/2pi is taken care of here.
 75      const double normpT = normfac/(2.0*0.2);
 76      const double normpTall = normfac/(2.0*M_PI*2.0*2.4);
 77
 78      for (auto& item : _h) {
 79        if (item.first.find("_dEta") != string::npos) {
 80          scale(item.second, normfac);
 81        }
 82        else if (item.first.find("_all") != string::npos) {
 83          scale(item.second, normpTall);
 84        }
 85        else {
 86          scale(item.second, normpT);
 87        }
 88      }
 89    }
 90
 91    /// @}
 92
 93
 94  private:
 95
 96    /// @{
 97    map<string, Histo1DPtr> _h;
 98
 99    string _sqs = "";
100    /// @}
101
102  };
103
104
105
106  RIVET_DECLARE_ALIASED_PLUGIN(CMS_2010_I845323, CMS_2010_S8547297);
107
108}