Rivet analyses referenceATLAS_2017_I1632756Photon + heavy flavour at 8 TeVExperiment: ATLAS (LHC) Inspire ID: 1632756 Status: VALIDATED Authors:
Beam energies: (4000.0, 4000.0) GeV Run details:
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// -*- C++ -*-
2#include "Rivet/Analysis.hh"
3#include "Rivet/Projections/FinalState.hh"
4#include "Rivet/Projections/PromptFinalState.hh"
5#include "Rivet/Projections/VisibleFinalState.hh"
6#include "Rivet/Projections/LeadingParticlesFinalState.hh"
7#include "Rivet/Projections/VetoedFinalState.hh"
8#include "Rivet/Projections/FastJets.hh"
9#include "Rivet/Projections/HeavyHadrons.hh"
10
11namespace Rivet {
12
13
14 /// @brief Measurement of prompt isolated photon + b/c-jet + X differential cross-sections
15 class ATLAS_2017_I1632756 : public Analysis {
16 public:
17
18 /// Constructor
19 RIVET_DEFAULT_ANALYSIS_CTOR(ATLAS_2017_I1632756);
20
21
22 /// Book histograms and initialise projections before the run
23 void init() {
24 // particles for photon isolation: no muons, no neutrinos
25 declare(VisibleFinalState(Cuts::abspid != PID::MUON), "caloParticles");
26
27 // Voronoi eta-phi tessellation with KT jets, for ambient energy density calculation
28 FastJets fj(FinalState(), JetAlg::KT, 0.5, JetMuons::NONE, JetInvisibles::NONE);
29 fj.useJetArea(new fastjet::AreaDefinition(fastjet::VoronoiAreaSpec()));
30 declare(fj, "KtJetsD05");
31
32 // Leading photon
33 LeadingParticlesFinalState photonfs(PromptFinalState(Cuts::abseta < 2.37 && Cuts::pT > 25*GeV));
34 photonfs.addParticleId(PID::PHOTON);
35 declare(photonfs, "LeadingPhoton");
36
37 // Jets
38 FastJets jetpro(FinalState(), JetAlg::ANTIKT, 0.4, JetMuons::DECAY, JetInvisibles::DECAY);
39 declare(jetpro, "Jets");
40
41 // Heavy hadrons
42 declare(HeavyHadrons(), "HeavyHadrons");
43
44 // Book the dsigma/dEt (in eta bins) histograms
45 // d02 and d03 are for photon+b; d04 and d05 are for photon+c
46 for (size_t i = 0; i < _eta_bins.size() - 1; ++i) {
47 if (fuzzyEquals(_eta_bins[i], 1.37)) continue; // This region is not used
48 int offset = i > 1? 1 : 2;
49 book(_h_Et_photonb[i], i + offset, 1, 1);
50 book(_h_Et_photonc[i], i + offset + 2, 1, 1);
51 }
52
53 }
54
55
56 /// Return eta bin for either dsigma/dET histogram (area_eta=false) or energy density correction (area_eta=true)
57 size_t _getEtaBin(double eta_w, bool area_eta) const {
58 const double eta = fabs(eta_w);
59 if (!area_eta) {
60 return binIndex(eta, _eta_bins);
61 } else {
62 return binIndex(eta, _eta_bins_areaoffset);
63 }
64 }
65
66
67 /// Perform the per-event analysis
68 void analyze(const Event& event) {
69
70 // Get the leading photon
71 const Particles& photons = apply<LeadingParticlesFinalState>(event, "LeadingPhoton").particlesByPt();
72 if (photons.empty()) vetoEvent;
73 const Particle& leadingPhoton = photons[0];
74
75 // Veto events with leading photon in ECAL crack
76 if (inRange(leadingPhoton.abseta(), 1.37, 1.56)) vetoEvent;
77
78 // Compute isolation energy in cone of radius .4 around photon (all particles except muons, neutrinos and leading photon)
79 FourMomentum mom_in_EtCone;
80 const Particles& fs = apply<FinalState>(event, "caloParticles").particles();
81 for (const Particle& p : fs) {
82 // increment if inside cone of 0.4
83 if (deltaR(leadingPhoton, p) < 0.4) mom_in_EtCone += p.momentum();
84 }
85 // Remove the photon energy from the isolation
86 mom_in_EtCone -= leadingPhoton.momentum();
87
88 // Get the area-filtered jet inputs for computing median energy density, etc.
89 vector<double> ptDensity;
90 vector< vector<double> > ptDensities(_eta_bins_areaoffset.size()-1);
91 const FastJets& fast_jets = apply<FastJets>(event, "KtJetsD05");
92 const auto clust_seq_area = fast_jets.clusterSeqArea();
93 for (const Jet& jet : fast_jets.jets()) {
94 const double area = clust_seq_area->area(jet);
95 if (area > 10e-4 && jet.abseta() < _eta_bins_areaoffset.back())
96 ptDensities.at( _getEtaBin(jet.abseta(), true) ).push_back(jet.pT()/area);
97 }
98
99 // Compute the median energy density, etc.
100 for (size_t b = 0; b < _eta_bins_areaoffset.size() - 1; ++b) {
101 const double ptmedian = (ptDensities[b].size() > 0) ? median(ptDensities[b]) : 0;
102 ptDensity.push_back(ptmedian);
103 }
104
105 // Compute the isolation energy correction (cone area*energy density)
106 const double etCone_area = PI * sqr(0.4);
107 const double correction = ptDensity[_getEtaBin(leadingPhoton.abseta(), true)] * etCone_area;
108
109 // Apply isolation cut on area-corrected value
110 // cut is Etiso < 4.8GeV + 4.2E-03 * Et_gamma.
111 if (mom_in_EtCone.Et() - correction > 4.8*GeV + 0.0042*leadingPhoton.Et()) vetoEvent;
112
113
114 // Get the leading jet
115 Jets jets = apply<FastJets>(event, "Jets").jetsByPt(Cuts::pT > 20*GeV);
116 idiscard(jets, deltaRLess(leadingPhoton, 0.4));
117 if (jets.empty()) vetoEvent;
118 const Jet& leadingJet = jets[0];
119
120 // Veto events with leading jet outside |y|<2.5
121 if (leadingJet.absrap() > 2.5) vetoEvent;
122
123 // Veto events with leading photon and leading jet with deltaR < 1.0
124 if (deltaR(leadingPhoton, leadingJet) < 1.0) vetoEvent;
125
126 // Veto events with leading jet not b-tagged (deltaR match with a B-hadron) nor c-tagged (deltaR match with a C-hadron)
127 const Particles& allBs = apply<HeavyHadrons>(event, "HeavyHadrons").bHadrons(Cuts::pT > 5*GeV);
128 bool bTagged = false;
129 for (const Particle& thisB : allBs) {
130 if(deltaR(thisB, leadingJet) < 0.3) {
131 bTagged = true;
132 break;
133 }
134 }
135
136 bool cTagged = false;
137 if (!bTagged) {
138 const Particles& allCs = apply<HeavyHadrons>(event, "HeavyHadrons").cHadrons(Cuts::pT > 5*GeV);
139 for (const Particle& thisC : allCs) {
140 if (deltaR(thisC, leadingJet) < 0.3) {
141 cTagged = true;
142 break;
143 }
144 }
145 if (!cTagged) vetoEvent;
146 }
147
148 // Fill histograms
149 const size_t eta_bin = _getEtaBin(leadingPhoton.abseta(), false);
150 if (bTagged) _h_Et_photonb[eta_bin]->fill(leadingPhoton.Et()/GeV);
151 if (cTagged) _h_Et_photonc[eta_bin]->fill(leadingPhoton.Et()/GeV);
152
153 }
154
155
156 /// Normalise histograms etc., after the run
157 void finalize() {
158 const double sf = crossSection() / (picobarn * sumOfWeights());
159 for (size_t i = 0; i < _eta_bins.size() - 1; ++i) {
160 if (fuzzyEquals(_eta_bins[i], 1.37)) continue; // This region is not used
161 scale(_h_Et_photonb[i], sf);
162 scale(_h_Et_photonc[i], sf);
163 }
164 }
165
166
167 private:
168
169 Histo1DPtr _h_Et_photonb[3];
170 Histo1DPtr _h_Et_photonc[3];
171
172 const vector<double> _eta_bins = { 0.00, 1.37, 1.56, 2.37 };
173 const vector<double> _eta_bins_areaoffset = { 0.0, 1.5, 3.0 };
174
175 };
176
177
178 RIVET_DECLARE_PLUGIN(ATLAS_2017_I1632756);
179
180
181}
|