rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

CMS_2010_I855299

Charged-particle transverse momentum and pseudorapidity spectra from proton-proton collisions at 7000 GeV
Experiment: CMS (LHC)
Inspire ID: 855299
Status: VALIDATED
Authors:
  • A. Knutsson
References: Beams: p+ p+
Beam energies: (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).

Charged particle spectra are measured in proton-proton collisions at center-of-mass energies 7000 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 ($p_\perp$) spectra from 0.1 to 2 GeV in bins of pseudorapidity ($\eta$) and the $p_\perp$ spectrum from 0.1 to 6 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 $p_\perp$. 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_I855299.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 7000 GeV
 9  class CMS_2010_I855299 : public Analysis {
10  public:
11
12    RIVET_DEFAULT_ANALYSIS_CTOR(CMS_2010_I855299);
13
14
15    void init() {
16      ChargedFinalState cfs((Cuts::etaIn(-2.5, 2.5)));
17      declare(cfs, "CFS");
18
19      for (int d=1; d<=3; d++) {
20        for (int y=1; y<=4; y++) {
21          _h_dNch_dpT.push_back(Histo1DPtr());
22          book(_h_dNch_dpT.back(), d, 1, y);
23        }
24      }
25
26      book(_h_dNch_dpT_all ,4, 1, 1);
27      book(_h_dNch_dEta ,5, 1, 1);
28    }
29
30
31    void analyze(const Event& event) {
32
33      //charged particles
34      const ChargedFinalState& charged = apply<ChargedFinalState>(event, "CFS");
35
36      for (const Particle& p : charged.particles()) {
37        //selecting only charged hadrons
38        if (! PID::isHadron(p.pid())) continue;
39
40        const double pT = p.pT();
41        const double eta = p.eta();
42
43        // The data is actually a duplicated folded distribution.  This should mimic it.
44        _h_dNch_dEta->fill(eta, 0.5);
45        _h_dNch_dEta->fill(-eta, 0.5);
46        if (fabs(eta) < 2.4 && pT > 0.1*GeV) {
47          if (pT < 6.0*GeV) {
48            _h_dNch_dpT_all->fill(pT/GeV, 1.0/(pT/GeV));
49            if (pT < 2.0*GeV) {
50              int ietabin = int(fabs(eta)/0.2);
51              _h_dNch_dpT[ietabin]->fill(pT/GeV);
52            }
53          }
54        }
55      }
56    }
57
58
59    void finalize() {
60      const double normfac = 1.0/sumOfWeights(); // Normalizing to unit eta is automatic
61      // The pT distributions in bins of eta must be normalized to unit eta.  This is a factor of 2
62      // for the |eta| times 0.2 (eta range).
63      // The pT distributions over all eta are normalized to unit eta (2.0*2.4) and by 1/2*pi*pT.
64      // The 1/pT part is taken care of in the filling.  The 1/2pi is taken care of here.
65      const double normpT = normfac/(2.0*0.2);
66      const double normpTall = normfac/(2.0*M_PI*2.0*2.4);
67
68      for (size_t ietabin=0; ietabin < _h_dNch_dpT.size(); ietabin++){
69        scale(_h_dNch_dpT[ietabin], normpT);
70      }
71      scale(_h_dNch_dpT_all, normpTall);
72      scale(_h_dNch_dEta, normfac);
73    }
74
75
76  private:
77
78    /// @{
79    std::vector<Histo1DPtr> _h_dNch_dpT;
80    Histo1DPtr _h_dNch_dpT_all;
81    Histo1DPtr _h_dNch_dEta;
82    /// @}
83
84  };
85
86
87
88  RIVET_DECLARE_ALIASED_PLUGIN(CMS_2010_I855299, CMS_2010_S8656010);
89
90}