Rivet analyses referenceATLAS_2016_I1467230Track-based minimum bias with low-pT tracks at 13 TeVExperiment: ATLAS (LHC) Inspire ID: 1467230 Status: VALIDATED Authors:
Beam energies: (6500.0, 6500.0) GeV Run details:
Measurements of distributions of charged particles produced in proton-proton collisions with a centre-of-mass energy of 13 TeV are presented. The data were recorded by the ATLAS detector at the LHC and correspond to an integrated luminosity of 151 $\mu$b$^{-1}$. The particles are required to have a transverse momentum greater than 100 MeV and an absolute pseudorapidity less than 2.5. The charged-particle multiplicity, its dependence on transverse momentum and pseudorapidity and the dependence of the mean transverse momentum on multiplicity are measured in events containing at least two charged particles satisfying the above kinematic criteria. The results are corrected for detector effects and compared to the predictions from several Monte Carlo event generators. Source code: ATLAS_2016_I1467230.cc 1// -*- C++ -*-
2#include "Rivet/Analysis.hh"
3#include "Rivet/Projections/ChargedFinalState.hh"
4
5namespace Rivet {
6
7
8 /// @brief ATLAS 13 TeV minimum bias analysis for low-pT tracks
9 class ATLAS_2016_I1467230 : public Analysis {
10 public:
11
12 /// Particle types included
13 enum PartTypes {
14 k_NoStrange,
15 k_AllCharged,
16 kNPartTypes
17 };
18
19 /// Phase space regions
20 enum regionID {
21 k_pt100_nch2_eta25,
22 kNregions
23 };
24
25
26 /// Default constructor
27 RIVET_DEFAULT_ANALYSIS_CTOR(ATLAS_2016_I1467230);
28
29
30 /// Initialization, called once before running
31 void init() {
32
33 for (int iT = 0; iT < kNPartTypes; ++iT) {
34 for (int iR = 0; iR < kNregions; ++iR) {
35 book(_sumW[iT][iR], "_sumW" + to_str(iT) + to_str(iR));
36 }
37 }
38
39 // Initialize and register projections
40 declare(ChargedFinalState(Cuts::abseta < 2.5 && Cuts::pT > 100*MeV), "CFS100_25");
41
42 for (int iT = 0; iT < kNPartTypes; ++iT) {
43 for (int iR = 0; iR < kNregions; ++iR) {
44 book(_hist_nch [iT][iR], 1, iR + 1, iT + 1);
45 book(_hist_pt [iT][iR], 2, iR + 1, iT + 1);
46 book(_hist_eta [iT][iR], 3, iR + 1, iT + 1);
47 book(_hist_ptnch[iT][iR], 4, iR + 1, iT + 1);
48 }
49 }
50
51 }
52
53
54 /// Fill histograms for the given particle selection and phase-space region
55 void fillPtEtaNch(const Particles& particles, int ptype, int iRegion) {
56
57 // Skip if event fails multiplicity cut
58 const size_t nch = particles.size();
59 if (nch < 2) return;
60
61 _sumW[ptype][iRegion]->fill();
62
63 // Fill nch
64 _hist_nch[ptype][iRegion]->fill(nch);
65
66 // Loop over particles, fill pT, eta and ptnch
67 for (const Particle& p : particles) {
68 const double pt = p.pT()/GeV;
69 const double eta = p.eta();
70 _hist_pt [ptype][iRegion]->fill(pt , 1.0/pt);
71 _hist_eta [ptype][iRegion]->fill(eta);
72 _hist_ptnch[ptype][iRegion]->fill(nch, pt);
73 }
74 }
75
76
77 /// Per-event analysis
78 void analyze(const Event& event) {
79
80 // Get all charged particles
81 const ChargedFinalState& cfs = apply<ChargedFinalState>(event, "CFS100_25");
82 const Particles& pall = cfs.particles();
83
84 // Get charged particles, filtered to omit charged strange baryons
85 const Cut& pcut = Cuts::abspid != PID::SIGMAMINUS && Cuts::abspid != PID::SIGMAPLUS && Cuts::abspid != PID::XIMINUS && Cuts::abspid != PID::OMEGAMINUS;
86 const Particles& pnostrange = cfs.particles(pcut);
87
88 // Fill all histograms
89 for (int iR = 0; iR < kNregions; ++iR) {
90 fillPtEtaNch(pall, k_AllCharged, iR);
91 fillPtEtaNch(pnostrange, k_NoStrange, iR);
92 }
93
94 }
95
96
97 /// Post-run data manipulation
98 void finalize() {
99
100 // Scale all histograms
101 for (int iT = 0; iT < kNPartTypes; ++iT) {
102 for (int iR = 0; iR < kNregions; ++iR) {
103 if (_sumW[iT][iR]->val() > 0) {
104 scale(_hist_nch[iT][iR], 1.0/ *_sumW[iT][iR]);
105 scale(_hist_pt [iT][iR], 1.0/ dbl(*_sumW[iT][iR])/TWOPI/5.);
106 scale(_hist_eta[iT][iR], 1.0/ *_sumW[iT][iR]);
107 }
108 }
109 }
110
111 }
112
113
114 private:
115
116 /// Weight sums
117 CounterPtr _sumW[kNPartTypes][kNregions];
118
119 /// @name Histogram arrays
120 /// @{
121 Histo1DPtr _hist_nch [kNPartTypes][kNregions];
122 Histo1DPtr _hist_pt [kNPartTypes][kNregions];
123 Histo1DPtr _hist_eta [kNPartTypes][kNregions];
124 Profile1DPtr _hist_ptnch [kNPartTypes][kNregions];
125 /// @}
126
127 };
128
129
130 RIVET_DECLARE_PLUGIN(ATLAS_2016_I1467230);
131
132}
|