Rivet Analyses Reference

ATLAS_2017_I1632756

Photon + heavy flavour at 8 TeV
Experiment: ATLAS (LHC)
Inspire ID: 1632756
Status: VALIDATED
Authors:
  • Sebastien Prince
References:Beams: p+ p+
Beam energies: (4000.0, 4000.0) GeV
Run details:
  • prompt photon + heavy-flavour jets

This Letter presents the measurement of differential cross sections of isolated prompt photons produced in association with a $b$-jet or a $c$-jet. These final states provide sensitivity to the heavy-flavour content of the proton and aspects related to the modelling of heavy-flavour quarks in perturbative QCD. The measurement uses proton-proton collision data at a centre-of-mass energy of 8 TeV recorded by the ATLAS detector at the LHC in 2012 corresponding to an integrated luminosity of up to 20.2 fb$^{-1}$. The differential cross sections are measured for each jet flavour with respect to the transverse energy of the leading photon in two photon pseudorapidity regions: $|\eta^\gamma| < 1.37$ and $1.56 < |\eta^\gamma| < 2.37$. The measurement covers photon transverse energies $25 < E^\gamma_\text{T} < 400$ GeV and $25 < E^\gamma_\text{T} < 350$ GeV respectively for the two $|\eta^\gamma|$ regions. For each jet flavour, the ratio of the cross sections in the two $|\eta^\gamma|$ regions is also measured. The measurement is corrected for detector effects and compared to leading-order and next-to-leading-order perturbative QCD calculations, based on various treatments and assumptions about the heavy-flavour content of the proton. Overall, the predictions agree well with the measurement, but some deviations are observed at high photon transverse energies. The total uncertainty in the measurement ranges between 13% and 66%, while the central $\gamma + b$ measurement exhibits the smallest uncertainty, ranging from 13% to 27%, which is comparable to the precision of the theoretical predictions.

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

namespace Rivet {


  /// @brief Measurement of prompt isolated photon + b/c-jet + X differential cross-sections
  class ATLAS_2017_I1632756 : public Analysis {
  public:

    /// Constructor
    RIVET_DEFAULT_ANALYSIS_CTOR(ATLAS_2017_I1632756);


    /// Book histograms and initialise projections before the run
    void init() {
      // particles for photon isolation: no muons, no neutrinos
      declare(VisibleFinalState(Cuts::abspid != PID::MUON), "caloParticles");

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

      // Leading photon
      LeadingParticlesFinalState photonfs(PromptFinalState(Cuts::abseta < 2.37 && Cuts::pT > 25*GeV));
      photonfs.addParticleId(PID::PHOTON);
      declare(photonfs, "LeadingPhoton");

      // Jets
      FastJets jetpro(FinalState(), FastJets::ANTIKT, 0.4, JetAlg::Muons::DECAY, JetAlg::Invisibles::DECAY);
      declare(jetpro, "Jets");

      // Heavy hadrons
      declare(HeavyHadrons(), "HeavyHadrons");

      // Book the dsigma/dEt (in eta bins) histograms
      // d02 and d03 are for photon+b; d04 and d05 are for photon+c
      for (size_t i = 0; i < _eta_bins.size() - 1; ++i) {
        if (fuzzyEquals(_eta_bins[i], 1.37)) continue; // This region is not used
        int offset = i > 1? 1 : 2;
        book(_h_Et_photonb[i], i + offset, 1, 1);
        book(_h_Et_photonc[i], i + offset + 2, 1, 1);
      }

    }


    /// Return eta bin for either dsigma/dET histogram (area_eta=false) or energy density correction (area_eta=true)
    size_t _getEtaBin(double eta_w, bool area_eta) const {
      const double eta = fabs(eta_w);
      if (!area_eta) {
        return binIndex(eta, _eta_bins);
      } else {
        return binIndex(eta, _eta_bins_areaoffset);
      }
    }


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

      // Get the leading photon
      const Particles& photons = apply<LeadingParticlesFinalState>(event, "LeadingPhoton").particlesByPt();
      if (photons.empty())  vetoEvent;
      const Particle& leadingPhoton = photons[0];

      // Veto events with leading photon in ECAL crack
      if (inRange(leadingPhoton.abseta(), 1.37, 1.56))  vetoEvent;

      // Compute isolation energy in cone of radius .4 around photon (all particles except muons, neutrinos and leading photon)
      FourMomentum mom_in_EtCone;
      const Particles& fs = apply<FinalState>(event, "caloParticles").particles();
      for (const Particle& p : fs) {
        // increment if inside cone of 0.4
        if (deltaR(leadingPhoton, p) < 0.4)  mom_in_EtCone += p.momentum();
      }
      // Remove the photon energy from the isolation
      mom_in_EtCone -= leadingPhoton.momentum();

      // Get the area-filtered jet inputs for computing median energy density, etc.
      vector<double> ptDensity;
      vector< vector<double> > ptDensities(_eta_bins_areaoffset.size()-1);
      const FastJets& fast_jets = apply<FastJets>(event, "KtJetsD05");
      const auto clust_seq_area = fast_jets.clusterSeqArea();
      for (const Jet& jet : fast_jets.jets()) {
        const double area = clust_seq_area->area(jet);
        if (area > 10e-4 && jet.abseta() < _eta_bins_areaoffset.back())
          ptDensities.at( _getEtaBin(jet.abseta(), true) ).push_back(jet.pT()/area);
      }

      // Compute the median energy density, etc.
      for (size_t b = 0; b < _eta_bins_areaoffset.size() - 1; ++b) {
        const double ptmedian = (ptDensities[b].size() > 0) ? median(ptDensities[b]) : 0;
        ptDensity.push_back(ptmedian);
      }

      // Compute the isolation energy correction (cone area*energy density)
      const double etCone_area = PI * sqr(0.4);
      const double correction = ptDensity[_getEtaBin(leadingPhoton.abseta(), true)] * etCone_area;

      // Apply isolation cut on area-corrected value
      // cut is Etiso < 4.8GeV + 4.2E-03 * Et_gamma.
      if (mom_in_EtCone.Et() - correction > 4.8*GeV + 0.0042*leadingPhoton.Et())  vetoEvent;


      // Get the leading jet
      Jets jets = apply<FastJets>(event, "Jets").jetsByPt(Cuts::pT > 20*GeV);
      ifilter_discard(jets, deltaRLess(leadingPhoton, 0.4));
      if (jets.empty())  vetoEvent;
      const Jet& leadingJet = jets[0];

      // Veto events with leading jet outside |y|<2.5
      if (leadingJet.absrap() > 2.5)  vetoEvent;

      // Veto events with leading photon and leading jet with deltaR < 1.0
      if (deltaR(leadingPhoton, leadingJet) < 1.0)  vetoEvent;

      // Veto events with leading jet not b-tagged (deltaR match with a B-hadron) nor c-tagged (deltaR match with a C-hadron)
      const Particles& allBs = apply<HeavyHadrons>(event, "HeavyHadrons").bHadrons(5*GeV);
      bool bTagged = false;
      for (const Particle& thisB : allBs) {
        if(deltaR(thisB, leadingJet) < 0.3) {
          bTagged = true;
          break;
        }
      }

      bool cTagged = false;
      if (!bTagged) {
        const Particles& allCs = apply<HeavyHadrons>(event, "HeavyHadrons").cHadrons(5*GeV);
        for (const Particle& thisC : allCs) {
          if (deltaR(thisC, leadingJet) < 0.3) {
            cTagged = true;
            break;
          }
        }
        if (!cTagged) vetoEvent;
      }

      // Fill histograms
      const size_t eta_bin = _getEtaBin(leadingPhoton.abseta(), false);
      if (bTagged) _h_Et_photonb[eta_bin]->fill(leadingPhoton.Et()/GeV);
      if (cTagged) _h_Et_photonc[eta_bin]->fill(leadingPhoton.Et()/GeV);

    }


    /// Normalise histograms etc., after the run
    void finalize() {
      const double sf = crossSection() / (picobarn * sumOfWeights());
      for (size_t i = 0; i < _eta_bins.size() - 1; ++i) {
        if (fuzzyEquals(_eta_bins[i], 1.37)) continue; // This region is not used
        scale(_h_Et_photonb[i], sf);
        scale(_h_Et_photonc[i], sf);
      }
    }


  private:

      Histo1DPtr _h_Et_photonb[3];
      Histo1DPtr _h_Et_photonc[3];

      const vector<double> _eta_bins = { 0.00, 1.37, 1.56, 2.37 };
      const vector<double> _eta_bins_areaoffset = { 0.0, 1.5, 3.0 };

  };


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


}