Rivet Analyses Reference

ATLAS_2012_I1093738

Isolated prompt photon + jet cross-section
Experiment: ATLAS (LHC)
Inspire ID: 1093738
Status: VALIDATED
Authors:
  • Giovanni Marchiori
References:Beams: p+ p+
Beam energies: (3500.0, 3500.0) GeV
Run details:
  • Inclusive photon+jet+X events at $\sqrt{s} = 7$ TeV.

A measurement of the production cross section for isolated photons in association with jets in $pp$ collisions at $\sqrt{s} = 7$ TeV. Photons with $|\eta|<1.37$ and $E_T>25$ GeV and jets with $|y|<4.4$ and $p_T>20$ GeV are selected. The differential cross section as a function of the photon transverse energy is measured, for three leading jet rapidity configurations, separately for the cases where the photon and jet rapidities have the same or the opposite sign. The measurement uses 37 pb$^{-1}$ of integrated luminosity collected with the ATLAS detector.

Source code: ATLAS_2012_I1093738.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
// -*- C++ -*-
#include "Rivet/Analysis.hh"
#include "Rivet/Projections/FinalState.hh"
#include "Rivet/Projections/LeadingParticlesFinalState.hh"
#include "Rivet/Projections/VetoedFinalState.hh"
#include "Rivet/Projections/FastJets.hh"

namespace Rivet {


  /// @brief Measurement of isolated gamma + jet + X differential cross-sections
  ///
  /// Inclusive isolated gamma + jet cross-sections, differential in pT(gamma), for
  /// various photon and jet rapidity configurations.
  ///
  /// @author Giovanni Marchiori
  class ATLAS_2012_I1093738 : public Analysis {
  public:

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


    // Book histograms and initialise projections before the run
    void init() {
      // Final state
      FinalState fs;
      declare(fs, "FS");

      // Voronoi eta-phi tessellation with KT jets, for ambient energy density calculation
      FastJets fj(fs, FastJets::KT, 0.5);
      fj.useJetArea(new fastjet::AreaDefinition(fastjet::VoronoiAreaSpec()));
      declare(fj, "KtJetsD05");

      // Leading photon
      LeadingParticlesFinalState photonfs(FinalState((Cuts::etaIn(-1.37, 1.37) && Cuts::pT >=  25.0*GeV)));
      photonfs.addParticleId(PID::PHOTON);
      declare(photonfs, "LeadingPhoton");

      // FS excluding the leading photon
      VetoedFinalState vfs(fs);
      vfs.addVetoOnThisFinalState(photonfs);
      declare(vfs, "JetFS");

      // Jets
      FastJets jetpro(vfs, FastJets::ANTIKT, 0.4);
      jetpro.useInvisibles();
      declare(jetpro, "Jets");

      book(_h_phbarrel_jetcentral_SS ,1, 1, 1);
      book(_h_phbarrel_jetmedium_SS  ,2, 1, 1);
      book(_h_phbarrel_jetforward_SS ,3, 1, 1);

      book(_h_phbarrel_jetcentral_OS ,4, 1, 1);
      book(_h_phbarrel_jetmedium_OS  ,5, 1, 1);
      book(_h_phbarrel_jetforward_OS ,6, 1, 1);
    }


    int getEtaBin(double eta, int what) const {
      const double aeta = fabs(eta);
      if (what == 0) return binIndex(aeta, _eta_bins_ph);
      if (what == 1) return binIndex(aeta, _eta_bins_jet);
      return binIndex(aeta, _eta_bins_areaoffset);
    }


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

      // Get the photon
      const FinalState& photonfs = apply<FinalState>(event, "LeadingPhoton");
      if (photonfs.particles().size() < 1) vetoEvent;
      const FourMomentum photon = photonfs.particles().front().momentum();

      // Get the jet
      Jets jets = apply<FastJets>(event, "Jets").jetsByPt(20.0*GeV);
      if (jets.empty()) vetoEvent;
      FourMomentum leadingJet = jets[0].momentum();

      // Require jet separated from photon
      if (deltaR(photon, leadingJet) < 1.0) vetoEvent;

      // Veto if leading jet is outside plotted rapidity regions
      if (leadingJet.absrap() > 4.4) vetoEvent;

      // Compute the jet pT densities
      vector< vector<double> > ptDensities(_eta_bins_areaoffset.size()-1);
      FastJets fastjets = apply<FastJets>(event, "KtJetsD05");
      const shared_ptr<fastjet::ClusterSequenceArea> clust_seq_area = fastjets.clusterSeqArea();
      for (const Jet& jet : fastjets.jets()) {
        const double area = clust_seq_area->area(jet); //< Implicit call to pseudojet()
        if (area > 1e-4 && jet.abseta() < _eta_bins_areaoffset.back()) {
          ptDensities.at(getEtaBin(jet.abseta(), 2)) += jet.pT()/area;
        }
      }

      // Compute the median event energy density
      /// @todo This looks equivalent to median(ptDensities[b]) -- isn't SKIPNHARDJETS meant to be used as an offset?
      const unsigned int SKIPNHARDJETS = 0;
      vector<double> ptDensity;
      for (size_t b = 0; b < _eta_bins_areaoffset.size()-1; b++) {
        double median = 0.0;
        if (ptDensities[b].size() > SKIPNHARDJETS) {
          std::sort(ptDensities[b].begin(), ptDensities[b].end());
          const int nDens = ptDensities[b].size() - SKIPNHARDJETS;
          if (nDens % 2 == 0) {
            median = (ptDensities[b][nDens/2]+ptDensities[b][(nDens-2)/2])/2;
          } else {
            median = ptDensities[b][(nDens-1)/2];
          }
        }
        ptDensity.push_back(median);
      }

      // Compute photon isolation with a standard ET cone
      const Particles fs = apply<FinalState>(event, "JetFS").particles();
      FourMomentum mom_in_EtCone;
      const double ISO_DR = 0.4;
      const double CLUSTER_ETA_WIDTH = 0.25*5.0;
      const double CLUSTER_PHI_WIDTH = (PI/128.)*7.0;
      for (const Particle& p : fs) {
        // Check if it's in the cone of .4
        if (deltaR(photon, p) >= ISO_DR) continue;
        // Check if it's in the 5x7 central core
        if (fabs(deltaEta(photon, p)) < CLUSTER_ETA_WIDTH*0.5 &&
            fabs(deltaPhi(photon, p)) < CLUSTER_PHI_WIDTH*0.5) continue;
        // Increment sum
        mom_in_EtCone += p.momentum();
      }

      // Figure out the correction (area*density)
      const double ETCONE_AREA = PI*ISO_DR*ISO_DR - CLUSTER_ETA_WIDTH*CLUSTER_PHI_WIDTH;
      const double correction = ptDensity[getEtaBin(photon.abseta(),2)] * ETCONE_AREA;

      // Require photon to be isolated
      if (mom_in_EtCone.Et()-correction > 4.0*GeV) vetoEvent;

      const int photon_jet_sign = sign( leadingJet.rapidity() * photon.rapidity() );

      // Fill histos
      const double abs_jet_rapidity = fabs(leadingJet.rapidity());
      const double photon_pt = photon.pT()/GeV;
      const double abs_photon_eta = fabs(photon.eta());
      if (abs_photon_eta < 1.37) {
        if (abs_jet_rapidity < 1.2) {
          if (photon_jet_sign >= 1) {
            _h_phbarrel_jetcentral_SS->fill(photon_pt);
          } else {
            _h_phbarrel_jetcentral_OS->fill(photon_pt);
          }
        } else if (abs_jet_rapidity < 2.8) {
          if (photon_jet_sign >= 1) {
            _h_phbarrel_jetmedium_SS->fill(photon_pt);
          } else {
            _h_phbarrel_jetmedium_OS->fill(photon_pt);
          }
        } else if (abs_jet_rapidity < 4.4) {
          if (photon_jet_sign >= 1) {
            _h_phbarrel_jetforward_SS->fill(photon_pt);
          } else {
            _h_phbarrel_jetforward_OS->fill(photon_pt);
          }
        }
      }

    }


    /// Normalise histograms etc., after the run
    void finalize() {
      scale(_h_phbarrel_jetcentral_SS, crossSection()/sumOfWeights());
      scale(_h_phbarrel_jetcentral_OS, crossSection()/sumOfWeights());
      scale(_h_phbarrel_jetmedium_SS, crossSection()/sumOfWeights());
      scale(_h_phbarrel_jetmedium_OS, crossSection()/sumOfWeights());
      scale(_h_phbarrel_jetforward_SS, crossSection()/sumOfWeights());
      scale(_h_phbarrel_jetforward_OS, crossSection()/sumOfWeights());
    }


  private:

    Histo1DPtr _h_phbarrel_jetcentral_SS, _h_phbarrel_jetmedium_SS, _h_phbarrel_jetforward_SS;
    Histo1DPtr _h_phbarrel_jetcentral_OS, _h_phbarrel_jetmedium_OS, _h_phbarrel_jetforward_OS;

    const vector<double> _eta_bins_ph = {0.0, 1.37, 1.52, 2.37};
    const vector<double> _eta_bins_jet = {0.0, 1.2, 2.8, 4.4};
    const vector<double> _eta_bins_areaoffset = {0.0, 1.5, 3.0};

  };


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


}