Rivet Analyses Reference

ATLAS_2016_I1426695

Track-based minimum bias at 8 TeV
Experiment: ATLAS (LHC)
Inspire ID: 1426695
Status: VALIDATED
Authors:
  • Roman Lysak
References:Beams: p+ p+
Beam energies: (4000.0, 4000.0) GeV
Run details:
  • $pp$ QCD interactions at 8 TeV. Diffractive events should be included. Multiple kinematic cuts should not be required.

This paper presents measurements of distributions of charged particles which are produced in proton-proton collisions at a centre-of-mass energy of $\sqrt{s} = 8$ TeV and recorded by the ATLAS detector at the LHC. A special dataset recorded in 2012 with a small number of interactions per beam crossing (below 0.004) and corresponding to an integrated luminosity of 160 $\mu\text{b}^{-1}$ was used. A minimum-bias trigger was utilised to select a data sample of more than 9 million collision events. The multiplicity, pseudorapidity, and transverse momentum distributions of charged particles are shown in different regions of kinematics and charged-particle multiplicity, including measurements of final states at high multiplicity. The results are corrected for detector effects and are compared to the predictions of various Monte Carlo event generator models which simulate the full hadronic final state.

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

namespace Rivet {


  /// @brief Track-based minimum bias at 8 TeV
  class ATLAS_2016_I1426695 : public Analysis {
  public:

    //phase space regions
    enum regionID {
      k_pt100_nch2 = 0,
      k_pt500_nch1,
      k_pt500_nch6,
      k_pt500_nch20,
      k_pt500_nch50,
      kNregions
    };


    /// Constructor
    RIVET_DEFAULT_ANALYSIS_CTOR(ATLAS_2016_I1426695);

    /// @name Analysis methods
    //@{

    /// Book histograms and initialise projections before the run
    void init() {

      for (int iR=0; iR < kNregions; ++iR)  {
        book(_sumW[iR], "_sumW" + to_str(iR))       ;
      }

      // Initialise and register projections
      declare(ChargedFinalState(Cuts::abseta < 2.5 && Cuts::pT > 100*MeV), "CFS_100");
      declare(ChargedFinalState(Cuts::abseta < 2.5 && Cuts::pT > 500*MeV), "CFS_500");

      // Book histograms
      for (int iR=0; iR < kNregions; ++iR)  {
        if (iR == k_pt100_nch2 || iR == k_pt500_nch1) {
          book(_hist_nch  [iR],  2 + iR, 1, 1);
          book(_hist_ptnch[iR], 14 + iR, 1, 1);
        }
        book(_hist_pt [iR], 4 + iR, 1, 1);
        book(_hist_eta[iR], 9 + iR, 1, 1);
      }
    }

    void fillPtEtaNch(const Particles particles, int nMin, int iRegion) {

      //skip if event fails multiplicity cut
      int nch =particles.size();
      if (nch < nMin)  return;

      _sumW[iRegion]->fill();

      // Fill nch
      if (iRegion == k_pt100_nch2 || iRegion == k_pt500_nch1) {
        _hist_nch[iRegion]->fill(nch);
      }

      for (const Particle &p : particles) {
      // Loop over particles, fill pT, eta and ptnch
        const double pt  = p.pT()/GeV;
        const double eta = p.eta();

        _hist_pt [iRegion]->fill(pt , 1.0/pt);
        _hist_eta[iRegion]->fill(eta);

        if (iRegion == k_pt100_nch2 || iRegion == k_pt500_nch1) {
          _hist_ptnch[iRegion]->fill(nch, pt);
        }
      } //end loop over particles
    }

    /// Perform the per-event analysis
    void analyze(const Event& event) {

      // Get charged particles, omitting some strange heavies
      const Cut& pcut = (
          (Cuts::abspid!=PID::SIGMAMINUS) && (Cuts::abspid!=PID::SIGMAPLUS) &&
          (Cuts::abspid!=PID::XIMINUS)    && (Cuts::abspid!=PID::OMEGAMINUS));
      const Particles& p_100 = apply<ChargedFinalState>(event, "CFS_100").particles(pcut);
      const Particles& p_500 = apply<ChargedFinalState>(event, "CFS_500").particles(pcut);

      fillPtEtaNch(p_100,  2, 0);
      fillPtEtaNch(p_500,  1, 1);
      fillPtEtaNch(p_500,  6, 2);
      fillPtEtaNch(p_500, 20, 3);
      fillPtEtaNch(p_500, 50, 4);
    }


    void finalize() {

      for (int iR = 0; iR < kNregions; ++iR)  {
        if (_sumW[iR]->val() > 0) {
          if (iR == k_pt100_nch2 || iR == k_pt500_nch1) {
            scale(_hist_nch[iR], 1.0/ *_sumW[iR]);
          }
          scale(_hist_pt [iR], 1.0/ dbl(*_sumW[iR])/TWOPI/5.);
          scale(_hist_eta[iR], 1.0/ *_sumW[iR]);
        }
      }
    }

    //@}


  private:

    CounterPtr _sumW[kNregions];

    /// @name Histograms
    Histo1DPtr   _hist_nch    [kNregions];
    Histo1DPtr   _hist_pt     [kNregions];
    Histo1DPtr   _hist_eta    [kNregions];
    Profile1DPtr _hist_ptnch  [kNregions];

  };


  RIVET_DECLARE_PLUGIN(ATLAS_2016_I1426695);

}