|
Rivet analyses reference
ATLAS_2016_I1467230
Track-based minimum bias with low-pT tracks at 13 TeV
Experiment: ATLAS (LHC)
Inspire ID: 1467230
Status: VALIDATED
Authors:
References:
Beams: p+ p+
Beam energies: (6500.0, 6500.0) GeV
Run details:
- $pp$ QCD interactions at 13 TeV. Diffractive events should be included. Multiple kinematic cuts should not be required.
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
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
126
127
128
129
130
131
132
133 | // -*- C++ -*-
#include "Rivet/Analysis.hh"
#include "Rivet/Projections/ChargedFinalState.hh"
namespace Rivet {
/// @brief ATLAS 13 TeV minimum bias analysis for low-pT tracks
class ATLAS_2016_I1467230 : public Analysis {
public:
/// Particle types included
enum PartTypes {
k_NoStrange,
k_AllCharged,
kNPartTypes
};
/// Phase space regions
enum regionID {
k_pt100_nch2_eta25,
kNregions
};
/// Default constructor
RIVET_DEFAULT_ANALYSIS_CTOR(ATLAS_2016_I1467230);
/// Initialization, called once before running
void init() {
for (int iT = 0; iT < kNPartTypes; ++iT) {
for (int iR = 0; iR < kNregions; ++iR) {
book(_sumW[iT][iR], "_sumW" + to_str(iT) + to_str(iR));
}
}
// Initialize and register projections
declare(ChargedFinalState(Cuts::abseta < 2.5 && Cuts::pT > 100*MeV), "CFS100_25");
for (int iT = 0; iT < kNPartTypes; ++iT) {
for (int iR = 0; iR < kNregions; ++iR) {
book(_hist_nch [iT][iR], 1, iR + 1, iT + 1);
book(_hist_pt [iT][iR], 2, iR + 1, iT + 1);
book(_hist_eta [iT][iR], 3, iR + 1, iT + 1);
book(_hist_ptnch[iT][iR], 4, iR + 1, iT + 1);
}
}
}
/// Fill histograms for the given particle selection and phase-space region
void fillPtEtaNch(const Particles& particles, int ptype, int iRegion) {
// Skip if event fails multiplicity cut
const size_t nch = particles.size();
if (nch < 2) return;
_sumW[ptype][iRegion]->fill();
// Fill nch
_hist_nch[ptype][iRegion]->fill(nch);
// Loop over particles, fill pT, eta and ptnch
for (const Particle& p : particles) {
const double pt = p.pT()/GeV;
const double eta = p.eta();
_hist_pt [ptype][iRegion]->fill(pt , 1.0/pt);
_hist_eta [ptype][iRegion]->fill(eta);
_hist_ptnch[ptype][iRegion]->fill(nch, pt);
}
}
/// Per-event analysis
void analyze(const Event& event) {
// Get all charged particles
const ChargedFinalState& cfs = apply<ChargedFinalState>(event, "CFS100_25");
const Particles& pall = cfs.particles();
// Get charged particles, filtered to omit charged strange baryons
const Cut& pcut = Cuts::abspid != PID::SIGMAMINUS && Cuts::abspid != PID::SIGMAPLUS && Cuts::abspid != PID::XIMINUS && Cuts::abspid != PID::OMEGAMINUS;
const Particles& pnostrange = cfs.particles(pcut);
// Fill all histograms
for (int iR = 0; iR < kNregions; ++iR) {
fillPtEtaNch(pall, k_AllCharged, iR);
fillPtEtaNch(pnostrange, k_NoStrange, iR);
}
}
/// Post-run data manipulation
void finalize() {
// Scale all histograms
for (int iT = 0; iT < kNPartTypes; ++iT) {
for (int iR = 0; iR < kNregions; ++iR) {
if (_sumW[iT][iR]->val() > 0) {
scale(_hist_nch[iT][iR], 1.0/ *_sumW[iT][iR]);
scale(_hist_pt [iT][iR], 1.0/ dbl(*_sumW[iT][iR])/TWOPI/5.);
scale(_hist_eta[iT][iR], 1.0/ *_sumW[iT][iR]);
}
}
}
}
private:
/// Weight sums
CounterPtr _sumW[kNPartTypes][kNregions];
/// @name Histogram arrays
//@{
Histo1DPtr _hist_nch [kNPartTypes][kNregions];
Histo1DPtr _hist_pt [kNPartTypes][kNregions];
Histo1DPtr _hist_eta [kNPartTypes][kNregions];
Profile1DPtr _hist_ptnch [kNPartTypes][kNregions];
//@}
};
// The hook for the plugin system
RIVET_DECLARE_PLUGIN(ATLAS_2016_I1467230);
}
|
|