Rivet Analyses Reference

ATLAS_2012_I1183818

Pseudorapidity dependence of the total transverse energy at 7 TeV
Experiment: ATLAS (LHC 7TeV)
Inspire ID: 1183818
Status: VALIDATED
Authors:
  • Robindra Prabhu
  • Peter Wijeratne
  • Roman Lysak
References:Beams: p+ p+
Beam energies: (3500.0, 3500.0) GeV
Run details:
  • pp QCD interactions at 7 TeV, min bias and di-jet events

The measurement of the sum of the transverse energy of particles as a function of particle pseudorapidity, eta, in proton-proton collisions at a centre-of-mass energy, $\sqrt{s} = 7 \text{TeV}$ using the ATLAS detector at the Large Hadron Collider. The measurements are performed in the region $|\eta| < 4.8$ for two event classes: those requiring the presence of particles with a low transverse momentum and those requiring particles with a significant transverse momentum (dijet events where both jets have $E_T > 20$ GeV). In the second dataset measurements are made in the region transverse to the hard scatter.

Source code: ATLAS_2012_I1183818.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/FinalState.hh"
#include "Rivet/Projections/ChargedFinalState.hh"
#include "Rivet/Projections/FastJets.hh"


/// @author Peter Wijeratne <paw@hep.ucl.ac.uk>
/// @author Robindra Prabhu <prabhu@cern.ch>
namespace Rivet {

  // A very basic analysis sensitive to ET flow in minbias and dijet events
  class ATLAS_2012_I1183818 : public Analysis {

  public:

    ATLAS_2012_I1183818()
      : Analysis("ATLAS_2012_I1183818")
    {}


  public:

    void init() {

      const FinalState cnfs((Cuts::etaIn(-4.8, 4.8)));
      const ChargedFinalState cfs((Cuts::etaIn(-2.5, 2.5) && Cuts::pT >=  250*MeV));
      declare(cnfs, "FS");
      declare(cfs, "CFS");

      const FastJets jetsAntiKt4(cnfs, FastJets::ANTIKT, 0.4);
      declare(jetsAntiKt4, "AntiKt4Jets");

      // ------- MINBIAS HISTOGRAMS --------
      //
      // MB event counter
      book(m_chargedEvents, "m_chargedEvents_MB");

      book(_h_ETflowEta ,1, 1, 1);
      book(_h_SumETbin1 ,3, 1, 1);
      book(_h_SumETbin2 ,4, 1, 1);
      book(_h_SumETbin3 ,5, 1, 1);
      book(_h_SumETbin4 ,6, 1, 1);
      book(_h_SumETbin5 ,7, 1, 1);
      book(_h_SumETbin6 ,8, 1, 1);

      // ------- DIJET HISTOGRAMS --------
      //
      // Dijet event counter
      book(m_events_dijets, "m_chargedEvents_dijets");

      // sumET
      book(_h_transETflowEta , 2, 1, 1);
      book(_h_transSumETbin1 , 9, 1, 1);
      book(_h_transSumETbin2 ,10, 1, 1);
      book(_h_transSumETbin3 ,11, 1, 1);
      book(_h_transSumETbin4 ,12, 1, 1);
      book(_h_transSumETbin5 ,13, 1, 1);
      book(_h_transSumETbin6 ,14, 1, 1);


    }


    void analyze(const Event& event) {

      const FinalState& cfs = apply<FinalState>(event, "CFS");

      bool isCharged = false;
      if (cfs.size() >= 2) {  // event selection: > 2 charged particles with pT > 250.MeV and |eta| < 2.5
        isCharged = true;
        m_chargedEvents->fill();
      }

      const FinalState& cnfs = apply<FinalState>(event, "FS");

      Particles particles;
      for( const Particle& p : cnfs.particles() ) {
        // enforce truth selection representing detected particle sensitivity
        double pp = p.p3().mod();
        if (PID::charge3(p.pid()) != 0 && pp < 0.5*GeV) continue;
        if (PID::charge3(p.pid()) == 0 && pp < 0.2*GeV) continue;

        particles.push_back(p);
      }


      // get jets
      const FastJets& jetsAntiKt4 = apply<FastJets>(event, "AntiKt4Jets");
      const Jets& jets = jetsAntiKt4.jetsByPt(20.0*GeV);

      // initialise sumET variables
      double sumETbin1 = 0;
      double sumETbin2 = 0;
      double sumETbin3 = 0;
      double sumETbin4 = 0;
      double sumETbin5 = 0;
      double sumETbin6 = 0;

      // if (passes event selection)
      if (isCharged) {

        for( const Particle& p : particles ) {

          ///calculate variables
          double ET = p.Et()/GeV;
          double eta = p.abseta();

          // fill histograms
          _h_ETflowEta->fill(eta, ET);

          if      (eta <  0.8) sumETbin1 += ET;
          else if (eta <  1.6) sumETbin2 += ET;
          else if (eta <  2.4) sumETbin3 += ET;
          else if (eta <  3.2) sumETbin4 += ET;
          else if (eta <  4.0) sumETbin5 += ET;
          else if (eta <= 4.8) sumETbin6 += ET;

        } // end of for

        _h_SumETbin1->fill(sumETbin1);
        _h_SumETbin2->fill(sumETbin2);
        _h_SumETbin3->fill(sumETbin3);
        _h_SumETbin4->fill(sumETbin4);
        _h_SumETbin5->fill(sumETbin5);
        _h_SumETbin6->fill(sumETbin6);
      }

      // --- do dijet analysis ---

      if ( jets.size() >= 2                       && // require at least two jets
           jets[0].Et() >= 20.*GeV     && // require two leading jets to pass ET cuts
           jets[1].Et() >= 20.*GeV     &&
           fabs(jets[0].eta()) < 2.5   && // require leading jets to be central
           fabs(jets[1].eta()) < 2.5   &&
           deltaPhi(jets[0], jets[1]) > 2.5       && // require back-to-back topology
           jets[1].Et()/jets[0].Et() >= 0.5) { //require ET-balance

        // found an event that satisfies dijet selection, now fill histograms...
        // initialise dijet sumET variables
        double trans_sumET_bin1 = 0.;
        double trans_sumET_bin2 = 0.;
        double trans_sumET_bin3 = 0.;
        double trans_sumET_bin4 = 0.;
        double trans_sumET_bin5 = 0.;
        double trans_sumET_bin6 = 0.;

        m_events_dijets->fill();

        // loop over all particles and check their relation to leading jet
        for( const Particle& particle : particles ) {

          // calculate variables
          double dPhi = deltaPhi( jets[0], particle.momentum() );
          double ET   = particle.Et()/GeV;
          double eta  = fabs(particle.eta());

          // Transverse region
          if ( dPhi > 1./3.*M_PI && dPhi < 2./3.*M_PI ) {
            _h_transETflowEta->fill( eta, ET );
            if      (eta <  0.8) { trans_sumET_bin1 += ET; }
            else if (eta <  1.6) { trans_sumET_bin2 += ET; }
            else if (eta <  2.4) { trans_sumET_bin3 += ET; }
            else if (eta <  3.2) { trans_sumET_bin4 += ET; }
            else if (eta <  4.0) { trans_sumET_bin5 += ET; }
            else if (eta <= 4.8) { trans_sumET_bin6 += ET; }
          }

        } // end loop over particles

        _h_transSumETbin1->fill( trans_sumET_bin1);
        _h_transSumETbin2->fill( trans_sumET_bin2);
        _h_transSumETbin3->fill( trans_sumET_bin3);
        _h_transSumETbin4->fill( trans_sumET_bin4);
        _h_transSumETbin5->fill( trans_sumET_bin5);
        _h_transSumETbin6->fill( trans_sumET_bin6);
      } // end of dijet selection cuts

    }


    void finalize() {
      /// several scale factors here:
      /// 1. nEvents (m_chargedEvents)
      /// 2. phase-space (2*M_PI)
      /// 3. double binning due to symmetrisation (2)
      scale( _h_ETflowEta, 1./m_chargedEvents->val()/(4.*M_PI) );
      scale( _h_SumETbin1, 1./m_chargedEvents->val() );
      scale( _h_SumETbin2, 1./m_chargedEvents->val() );
      scale( _h_SumETbin3, 1./m_chargedEvents->val() );
      scale( _h_SumETbin4, 1./m_chargedEvents->val() );
      scale( _h_SumETbin5, 1./m_chargedEvents->val() );
      scale( _h_SumETbin6, 1./m_chargedEvents->val() );

      //Dijet analysis

      // Dijet scale factors:
      //1. number of events passing dijet selection
      //2. phase-space: 1. / 2/3*M_PI
      //3. double binning due to symmetrisation in |eta| plots : 1/2
      scale( _h_transETflowEta, 1./m_events_dijets->val() * 1./(4./3.*M_PI) );
      scale( _h_transSumETbin1, 1./m_events_dijets->val() );
      scale( _h_transSumETbin2, 1./m_events_dijets->val() );
      scale( _h_transSumETbin3, 1./m_events_dijets->val() );
      scale( _h_transSumETbin4, 1./m_events_dijets->val() );
      scale( _h_transSumETbin5, 1./m_events_dijets->val() );
      scale( _h_transSumETbin6, 1./m_events_dijets->val() );
    }

  private:

    // Event counts
    CounterPtr m_chargedEvents;
    CounterPtr m_events_dijets;

    // Minbias-analysis: variable + histograms
    Histo1DPtr _h_ETflowEta;
    Histo1DPtr _h_SumETbin1;
    Histo1DPtr _h_SumETbin2;
    Histo1DPtr _h_SumETbin3;
    Histo1DPtr _h_SumETbin4;
    Histo1DPtr _h_SumETbin5;
    Histo1DPtr _h_SumETbin6;

    // Transverse region
    Histo1DPtr _h_transETflowEta;
    Histo1DPtr _h_transSumETbin1;
    Histo1DPtr _h_transSumETbin2;
    Histo1DPtr _h_transSumETbin3;
    Histo1DPtr _h_transSumETbin4;
    Histo1DPtr _h_transSumETbin5;
    Histo1DPtr _h_transSumETbin6;

  };


  RIVET_DECLARE_PLUGIN(ATLAS_2012_I1183818);

}