rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

ATLAS_2011_I921594

Inclusive isolated prompt photon analysis with full 2010 LHC data
Experiment: ATLAS (LHC 7TeV)
Inspire ID: 921594
Status: VALIDATED
Authors:
  • Giovanni Marchiori
References: Beams: p+ p+
Beam energies: (3500.0, 3500.0) GeV
Run details:
  • Inclusive photon $+ X$ events (primary $\gamma$+jet events) at $\sqrt{s} = 7$ TeV.

A measurement of the differential cross-section for the inclusive production of isolated prompt photons in $pp$ collisions at a center-of-mass energy $\sqrt{s} = 7$ TeV is presented. The measurement covers the pseudorapidity ranges $|\eta|<1.37$ and $1.52<|\eta|<2.37$ in the transverse energy range $45<E_T<400$ GeV. The results are based on an integrated luminosity of 35 pb$^{-1}$, collected with the ATLAS detector at the LHC. The yields of the signal photons are measured using a data-driven technique, based on the observed distribution of the hadronic energy in a narrow cone around the photon candidate and the photon selection criteria. The results are compared with next-to-leading order perturbative QCD calculations and found to be in good agreement over four orders of magnitude in cross-section.

Source code: ATLAS_2011_I921594.cc
  1// -*- C++ -*-
  2#include "Rivet/Analysis.hh"
  3#include "Rivet/Projections/FinalState.hh"
  4#include "Rivet/Projections/LeadingParticlesFinalState.hh"
  5#include "Rivet/Projections/FastJets.hh"
  6
  7namespace Rivet {
  8
  9
 10  /// Inclusive isolated prompt photon analysis with full 2010 LHC data
 11  class ATLAS_2011_I921594 : public Analysis {
 12  public:
 13
 14    /// Constructor
 15    RIVET_DEFAULT_ANALYSIS_CTOR(ATLAS_2011_I921594);
 16
 17
 18    /// Book histograms and initialise projections before the run
 19    void init() {
 20
 21      FinalState fs;
 22      declare(fs, "FS");
 23
 24      // Consider the final state jets for the energy density calculation
 25      FastJets fj(fs, JetAlg::KT, 0.5);
 26      fj.useJetArea(new fastjet::AreaDefinition(fastjet::VoronoiAreaSpec()));
 27      declare(fj, "KtJetsD05");
 28
 29      // Consider the leading pt photon with |eta|<2.37 and pT>45 GeV
 30      LeadingParticlesFinalState photonfs(FinalState((Cuts::etaIn(-2.37, 2.37) && Cuts::pT >=  45*GeV)));
 31      photonfs.addParticleId(PID::PHOTON);
 32      declare(photonfs, "LeadingPhoton");
 33
 34      // Book the dsigma/dEt (in eta bins) histograms
 35      size_t d = 1;
 36      for (size_t i = 0; i < _eta_bins.size()-1; ++i) {
 37        if (fuzzyEquals(_eta_bins[i], 1.37)) continue; // skip this bin
 38        book(_h_Et_photon[i], d, 1, 1);
 39        ++d;
 40      }
 41    }
 42
 43
 44    /// Return eta bin for either dsigma/dET histogram (area_eta=false) or energy density correction (area_eta=true)
 45    size_t _getEtaBin(double eta, bool area_eta) const {
 46      const double aeta = fabs(eta);
 47      return (!area_eta) ? binIndex(aeta, _eta_bins) : binIndex(aeta, _eta_bins_areaoffset);
 48    }
 49
 50
 51    /// Perform the per-event analysis
 52    void analyze(const Event& event) {
 53      // Retrieve leading photon
 54      const Particles& photons = apply<LeadingParticlesFinalState>(event, "LeadingPhoton").particles();
 55      if (photons.size() != 1) vetoEvent;
 56      const Particle& leadingPhoton = photons[0];
 57
 58      // Veto events with photon in ECAL crack
 59      if (inRange(leadingPhoton.abseta(), 1.37, 1.52)) vetoEvent;
 60
 61      // Compute isolation energy in cone of radius .4 around photon (all particles)
 62      FourMomentum mom_in_EtCone;
 63      Particles fs = apply<FinalState>(event, "FS").particles();
 64      for (const Particle& p : fs) {
 65        // Check if it's outside the cone of 0.4
 66        if (deltaR(leadingPhoton, p) >= 0.4) continue;
 67        // Don't count particles in the 5x7 central core
 68        if (deltaEta(leadingPhoton, p) < .025*5.0*0.5 &&
 69            deltaPhi(leadingPhoton, p) < (PI/128.)*7.0*0.5) continue;
 70        // Increment isolation energy
 71        mom_in_EtCone += p.momentum();
 72      }
 73
 74      // Get the area-filtered jet inputs for computing median energy density, etc.
 75      vector< vector<double> > ptDensities(_eta_bins_areaoffset.size()-1);
 76      FastJets fast_jets = apply<FastJets>(event, "KtJetsD05");
 77      const shared_ptr<fastjet::ClusterSequenceArea> clust_seq_area = fast_jets.clusterSeqArea();
 78      for (const Jet& jet : fast_jets.jets()) {
 79        const double area = clust_seq_area->area(jet); //< Implicit call to .pseudojet()
 80        if (area > 1e-4 && jet.abseta() < _eta_bins_areaoffset.back())
 81          ptDensities.at( _getEtaBin(jet.abseta(), true) ).push_back(jet.pT()/area);
 82      }
 83
 84      // Compute the median energy density, etc.
 85      vector<double> ptDensity;
 86      for (size_t b = 0; b < _eta_bins_areaoffset.size()-1; b++) {
 87        ptDensity += ptDensities[b].empty() ? 0 : median(ptDensities[b]);
 88      }
 89
 90      // Compute the isolation energy correction (cone area*energy density)
 91      const double ETCONE_AREA = M_PI*sqr(0.4) - (7.0*.025)*(5.0*PI/128.);
 92      const double correction = ptDensity[_getEtaBin(leadingPhoton.abseta(), true)] * ETCONE_AREA;
 93
 94      // Apply isolation cut on area-corrected value
 95      if (mom_in_EtCone.Et() - correction > 4*GeV) vetoEvent;
 96
 97      // Fill histograms
 98      const size_t eta_bin = _getEtaBin(leadingPhoton.abseta(), false);
 99      _h_Et_photon[eta_bin]->fill(leadingPhoton.Et());
100    }
101
102
103    /// Normalise histograms etc., after the run
104    void finalize() {
105      for (size_t i = 0; i < _eta_bins.size()-1; i++) {
106        if (fuzzyEquals(_eta_bins[i], 1.37)) continue;
107        scale(_h_Et_photon[i], crossSection()/picobarn/sumOfWeights());
108      }
109    }
110
111
112  private:
113
114    Histo1DPtr _h_Et_photon[5];
115
116    const vector<double> _eta_bins = {0.00, 0.60, 1.37, 1.52, 1.81, 2.37};
117    const vector<double> _eta_bins_areaoffset = {0.0, 1.5, 3.0};
118
119  };
120
121
122  RIVET_DECLARE_PLUGIN(ATLAS_2011_I921594);
123
124}