Rivet Analyses Reference

ATLAS_2012_I1125575

Studies of the underlying event at 7 TeV with track-jets
Experiment: ATLAS (LHC)
Inspire ID: 1125575
Status: VALIDATED
Authors:
  • Kiran Joshi
References:Beams: p+ p+
Beam energies: (3500.0, 3500.0) GeV
Run details:
  • Minimum bias events at 7 TeV.

Distributions sensitive to the underlying event are studied in events containing one or more charged particles. Jets are reconstructed using the anti-$k_t$ algorithm with radius parameter $R$ varying between 0.2 and 1.0. Distributions of the charged-particle multiplicity, the scalar sum of the transverse momentum of charged particles, and the average charged-particle pT are measured as functions of $p_\perp^\text{jet}$ in regions transverse to and opposite the leading jet for $4 \text{GeV} < p_\perp^\text{jet} < 100 \text{GeV}$. In addition, the $R$-dependence of the mean values of these observables is studied.

Source code: ATLAS_2012_I1125575.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
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
// -*- C++ -*-
#include "Rivet/Analysis.hh"
#include "Rivet/Projections/ChargedFinalState.hh"
#include "Rivet/Projections/FastJets.hh"
#include "Rivet/Tools/BinnedHistogram.hh"

namespace Rivet {


  /// ATLAS charged particle jet underlying event and jet radius dependence
  class ATLAS_2012_I1125575 : public Analysis {
  public:

    /// @name Constructors etc.
    //@{

    /// Constructor
    ATLAS_2012_I1125575()
      : Analysis("ATLAS_2012_I1125575")
    {    }

    //@}


    /// @name Analysis methods
    //@{

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

      const ChargedFinalState jet_input((Cuts::etaIn(-2.5, 2.5) && Cuts::pT >=  0.5*GeV));
      declare(jet_input, "JET_INPUT");

      const ChargedFinalState track_input((Cuts::etaIn(-1.5, 1.5) && Cuts::pT >=  0.5*GeV));
      declare(track_input, "TRACK_INPUT");

      const FastJets jets02(jet_input, FastJets::ANTIKT, 0.2);
      declare(jets02, "JETS_02");

      const FastJets jets04(jet_input, FastJets::ANTIKT, 0.4);
      declare(jets04, "JETS_04");

      const FastJets jets06(jet_input, FastJets::ANTIKT, 0.6);
      declare(jets06, "JETS_06");

      const FastJets jets08(jet_input, FastJets::ANTIKT, 0.8);
      declare(jets08, "JETS_08");

      const FastJets jets10(jet_input, FastJets::ANTIKT, 1.0);
      declare(jets10, "JETS_10");

      // Mean number of tracks
      initializeProfiles(_h_meanNch, 1);

      // Mean of the average track pT in each region
      initializeProfiles(_h_meanPtAvg, 2);

      // Mean of the scalar sum of track pT in each region
      initializeProfiles(_h_meanPtSum, 3);

      // Distribution of Nch, in bins of leading track-jet pT
      initializeHistograms(_h_Nch, 4);

      // Distribution of average track-jet pT, in bins of leading track-jet pT
      initializeHistograms(_h_PtAvg, 5);

      // Distribution of sum of track-jet pT, in bins of leading track-jet pT
      initializeHistograms(_h_PtSum, 6);

      for (int i = 0; i < 5; ++i)
        book(_nEvents[i], "nEvents_"+to_str(i));
    }


    void initializeProfiles(Profile1DPtr plots[5][2], int distribution) {
      for (int i = 0; i < 5; ++i) {
        for (int j = 0; j < 2; ++j) {
          book(plots[i][j] ,distribution, i+1, j+1);
        }
      }
    }


    void initializeHistograms(BinnedHistogram plots[5][2], int distribution) {
      Scatter2D refscatter = refData(1, 1, 1);
      for (int i = 0; i < 5; ++i) {
        for (int y = 0; y < 2; ++y) {
          for (size_t j = 0; j < refscatter.numPoints(); ++j) {
            int histogram_number = ((j+1)*2)-((y+1)%2);
            double low_edge = refscatter.point(j).xMin();
            double high_edge = refscatter.point(j).xMax();
            Histo1DPtr tmp;
            plots[i][y].add(low_edge, high_edge, book(tmp, distribution, i+1, histogram_number));
          }
        }
      }
    }


    /// Perform the per-event analysis
    void analyze(const Event& event) {
      vector<Jets*> all_jets;
      Jets jets_02 = apply<FastJets>(event, "JETS_02").jetsByPt(Cuts::pT > 4*GeV && Cuts::abseta < 1.5);
      all_jets.push_back(&jets_02);
      Jets jets_04 = apply<FastJets>(event, "JETS_04").jetsByPt(Cuts::pT > 4*GeV && Cuts::abseta < 1.5);
      all_jets.push_back(&jets_04);
      Jets jets_06 = apply<FastJets>(event, "JETS_06").jetsByPt(Cuts::pT > 4*GeV && Cuts::abseta < 1.5);
      all_jets.push_back(&jets_06);
      Jets jets_08 = apply<FastJets>(event, "JETS_08").jetsByPt(Cuts::pT > 4*GeV && Cuts::abseta < 1.5);
      all_jets.push_back(&jets_08);
      Jets jets_10 = apply<FastJets>(event, "JETS_10").jetsByPt(Cuts::pT > 4*GeV && Cuts::abseta < 1.5);
      all_jets.push_back(&jets_10);

      // Count the number of tracks in the away and transverse regions, for each set of jets
      double n_ch[5][2] = { {0,0}, {0,0}, {0,0}, {0,0}, {0,0} };
      // Also add up the sum pT
      double sumpt[5][2] = { {0,0}, {0,0}, {0,0}, {0,0}, {0,0} };
      // ptmean = sumpt / n_ch
      double ptavg[5][2] = { {0,0}, {0,0}, {0,0}, {0,0}, {0,0} };
      // lead jet pT defines which bin we want to fill
      double lead_jet_pts[5] = {0.0};

      // Loop over each of the jet radii:
      for (int i = 0; i < 5; ++i) {
        if (all_jets[i]->size() < 1) continue;

        // Find the lead jet pT
        lead_jet_pts[i] = all_jets[i]->at(0).pT();

        // Loop over each of the charged particles
        const Particles& tracks = apply<ChargedFinalState>(event, "TRACK_INPUT").particlesByPt();
        for(const Particle& t : tracks) {

          // Get the delta-phi between the track and the leading jet
          double dphi = deltaPhi(all_jets[i]->at(0), t);

          // Find out which region this puts it in.
          // 0 = away region, 1 = transverse region, 2 = toward region
          int region = region_index(dphi);

          // If the track is in the toward region, ignore it.
          if (region == 2) continue;

          // Otherwise, increment the relevant counters
          ++n_ch[i][region];
          sumpt[i][region] += t.pT();

        }
        // Calculate the pT_avg for the away and transverse regions.
        // (And make sure we don't try to divide by zero.)
        ptavg[i][0] = (n_ch[i][0] == 0 ? 0.0 : sumpt[i][0] / n_ch[i][0]);
        ptavg[i][1] = (n_ch[i][1] == 0 ? 0.0 : sumpt[i][1] / n_ch[i][1]);

        _nEvents[i]->fill();
      }

      fillProfiles(_h_meanNch,    n_ch, lead_jet_pts, 1.0 / (2*PI));
      fillProfiles(_h_meanPtAvg, ptavg, lead_jet_pts, 1.0);
      fillProfiles(_h_meanPtSum, sumpt, lead_jet_pts, 1.0 / (2*PI));

      fillHistograms(_h_Nch,    n_ch, lead_jet_pts);
      fillHistograms(_h_PtAvg, ptavg, lead_jet_pts);
      fillHistograms(_h_PtSum, sumpt, lead_jet_pts);
    }


    void fillProfiles(Profile1DPtr plots[5][2], double var[5][2], double lead_pt[5], double scale) {
      for (int i=0; i<5; ++i) {
        double pt = lead_pt[i];
        for (int j=0; j<2; ++j) {
          double v = var[i][j];
          plots[i][j]->fill(pt, v*scale);
        }
      }
    }


    void fillHistograms(BinnedHistogram plots[5][2], double var[5][2], double lead_pt[5]) {
      for (int i=0; i<5; ++i) {
        double pt = lead_pt[i];
        for (int j=0; j<2; ++j) {
          double v = var[i][j];
          plots[i][j].fill(pt, v);
        }
      }
    }


    int region_index(double dphi) {
      assert(inRange(dphi, 0.0, PI, CLOSED, CLOSED));
      if (dphi < PI/3.0) return 2;
      if (dphi < 2*PI/3.0) return 1;
      return 0;
    }


    /// Normalise histograms etc., after the run
    void finalize() {
      finalizeHistograms(_h_Nch);
      finalizeHistograms(_h_PtAvg);
      finalizeHistograms(_h_PtSum);
    }


    void finalizeHistograms(BinnedHistogram plots[5][2]) {
      for (int i = 0; i < 5; ++i) {
        for (int j = 0; j < 2; ++j) {
          vector<Histo1DPtr> histos = plots[i][j].histos();
          for(Histo1DPtr h : histos) {
            scale(h, 1.0/ *_nEvents[i]);
          }
        }
      }
    }

    //@}


  private:

    // Data members like post-cuts event weight counters go here
    CounterPtr _nEvents[5];

    Profile1DPtr _h_meanNch[5][2];
    Profile1DPtr _h_meanPtAvg[5][2];
    Profile1DPtr _h_meanPtSum[5][2];

    BinnedHistogram _h_Nch[5][2];
    BinnedHistogram _h_PtAvg[5][2];
    BinnedHistogram _h_PtSum[5][2];

  };



  // The hook for the plugin system
  RIVET_DECLARE_PLUGIN(ATLAS_2012_I1125575);

}