rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

ATLAS_2024_I2768921

ttbar + gamma at 13 TeV
Experiment: ATLAS (LHC)
Inspire ID: 2768921
Status: VALIDATED
Authors:
  • Carmen Diez Pardos
  • Stefanie Mueller
References: Beams: p+ p+
Beam energies: (6500.0, 6500.0) GeV
Run details:
  • pp -> ttbar + photon at 13TeV

Inclusive and differential cross-sections are measured at particle level for the associated production of a top quark pair and a photon ($t\bar{t}\gamma$). The analysis is performed using an integrated luminosity of 140 fb$^{-1}$ of proton-proton collisions at a centre-of-mass energy of 13 TeV collected by the ATLAS detector. The measurements are performed in the single-lepton and dilepton top quark pair decay channels focusing on $t\bar{t}\gamma$ topologies where the photon is radiated from an initial-state parton or one of the top quarks. The absolute and normalised differential cross-sections are measured for several variables characterising the photon, lepton and jet kinematics as well as the angular separation between those objects. The observables are found to be in good agreement with the Monte Carlo predictions. The photon transverse momentum differential distribution is used to set limits on effective field theory parameters related to the electroweak dipole moments of the top quark. The combined limits using the photon and the $Z$ boson transverse momentum measured in $t\bar{t}$ production in associations with a $Z$ boson are also set.

Source code: ATLAS_2024_I2768921.cc
  1// -*- C++ -*-
  2#include "Rivet/Analysis.hh"
  3#include "Rivet/Projections/FinalState.hh"
  4#include "Rivet/Projections/ChargedFinalState.hh"
  5#include "Rivet/Projections/FastJets.hh"
  6#include "Rivet/Projections/VetoedFinalState.hh"
  7#include "Rivet/Projections/InvisibleFinalState.hh"
  8#include "Rivet/Projections/LeptonFinder.hh"
  9#include "Rivet/Projections/PromptFinalState.hh"
 10
 11namespace Rivet {
 12
 13  /// @brief ttbar + gamma at 13 TeV
 14  class ATLAS_2024_I2768921 : public Analysis {
 15  public:
 16
 17
 18    // Constructor
 19    RIVET_DEFAULT_ANALYSIS_CTOR(ATLAS_2024_I2768921);
 20
 21
 22    // Book histograms and initialise projections before the run
 23    void init() {
 24
 25      // Set running mode (default is both channels)
 26      _mode = 3;
 27      if ( getOption("LMODE") == "SINGLE" )   _mode = 1;
 28      if ( getOption("LMODE") == "DILEPTON" ) _mode = 2;
 29      if ( getOption("LMODE") == "ALL" )      _mode = 3;
 30
 31      // Cuts
 32      Cut eta_full = Cuts::abseta < 5.0;
 33      Cut lepCuts  = Cuts::abseta < 2.5 && Cuts::pT > 25.0*GeV;
 34      Cut lepCutsVeto = Cuts::abseta < 2.5 && Cuts::pT > 7.0*GeV;
 35
 36      // Charged particles for signal photon isolation
 37      declare(ChargedFinalState(), "CFS");
 38
 39      // Signal photons
 40      PromptFinalState photons(Cuts::abspid == PID::PHOTON && Cuts::pT > 15*GeV && Cuts::abseta < 2.5);
 41      declare(photons, "Photons");
 42
 43      // Invisibles
 44      InvisibleFinalState prompt_invis(OnlyPrompt::YES, TauDecaysAs::PROMPT);
 45
 46      // Leptons
 47      PromptFinalState leptons(Cuts::abspid == PID::MUON || Cuts::abspid == PID::ELECTRON, TauDecaysAs::PROMPT);
 48
 49      // Dress the leptons
 50      FinalState dressPhotons(Cuts::abspid == PID::PHOTON);
 51      LeptonFinder dressedLeptons(leptons, dressPhotons, 0.1, lepCuts);
 52      declare(dressedLeptons, "DressedLeptons");
 53      // Full dressed leptons
 54      LeptonFinder fulldressedleptons(leptons, dressPhotons, 0.1, eta_full);
 55
 56      // Add new lepton collection for veto of 2nd lepton pt > 7 GeV
 57      LeptonFinder vetodressedLeptons(leptons, dressPhotons, 0.1, lepCutsVeto);
 58      declare(vetodressedLeptons, "VetoDressedLeptons");
 59
 60      // Jet alg input
 61      VetoedFinalState vfs;
 62      vfs.addVetoOnThisFinalState(fulldressedleptons);
 63      vfs.addVetoOnThisFinalState(prompt_invis);
 64
 65      // Jet clustering
 66      FastJets jets(vfs, JetAlg::ANTIKT, 0.4, JetMuons::ALL, JetInvisibles::ALL);
 67      declare(jets, "Jets");
 68
 69      dualbook("sldl_ph_pt",     18);
 70      dualbook("sldl_ph_abseta", 20);
 71
 72      // Book histograms
 73      if (_mode == 1 || _mode == 3) { // Single lepton channel
 74        dualbook("sl_ph_bjet_dR",  6);
 75        dualbook("sl_ph_l_dR",     8);
 76        dualbook("sl_l_jet_dR",   10);
 77        dualbook("sl_ph_pt",      44);
 78        dualbook("sl_ph_abseta",  46);
 79        dualbook("sl_jet1_pt",    48);
 80      }
 81
 82      if (_mode == 2 || _mode == 3) { // Dilepton channel
 83        dualbook("dl_ph_bjet_dR", 12);
 84        dualbook("dl_ph_l_dR",    14);
 85        dualbook("dl_l_jet_dR",   16);
 86        dualbook("dl_ph_pt",      50);
 87        dualbook("dl_ph_abseta",  52);
 88        dualbook("dl_jet1_pt",    54);
 89      }
 90    }
 91
 92
 93    void analyze(const Event& event) {
 94
 95      FourMomentum lepton(0.,0.,0.,0.);
 96
 97      // Fetch objects
 98      Particles photons          = apply<PromptFinalState>(event, "Photons").particlesByPt();
 99      DressedLeptons leptons     = apply<LeptonFinder>(event, "DressedLeptons").dressedLeptons();
100      DressedLeptons vetoleptons = apply<LeptonFinder>(event, "VetoDressedLeptons").dressedLeptons();
101      Jets jets                  = apply<FastJets>(event, "Jets").jetsByPt(Cuts::abseta < 2.5 && Cuts::pT > 25*GeV);
102
103
104      // Immediate veto on events without photon
105      if (photons.empty())  vetoEvent;
106
107      // Overlap removal of leptons near jets
108      idiscardIfAnyDeltaRLess(leptons, jets, 0.4);
109      idiscardIfAnyDeltaRLess(vetoleptons, jets, 0.4);
110
111      // Photon isolation with overlap removal of jets near isolated photon
112      Particles iso_photons;
113      const ChargedFinalState& chargedFS = apply<ChargedFinalState>(event, "CFS");
114      for (const Particle& ph : photons) {
115        const double conePt = sum(chargedFS.particles(deltaRLess(ph, 0.2)), Kin::pT, 0.0);
116        if (conePt/ph.pT() < 0.05) {
117          iso_photons += ph;
118          idiscard(jets, deltaRLess(ph, 0.4));
119        }
120      }
121
122      // Veto on events without one good photon
123      if (iso_photons.size() != 1) vetoEvent;
124
125      // Require at least one bjet
126      Jets bjets;
127      for (const Jet& jet : jets) {
128       if (jet.bTagged(Cuts::pT > 5*GeV)) bjets += jet;
129      }
130      if (bjets.empty()) vetoEvent;
131
132
133      const Particle& photon = iso_photons[0];
134      // Apply cuts to isolated photon, after requiring exactly 1 photon
135      if (photon.pT() <= 20*GeV || photon.abseta() >= 2.37) vetoEvent;
136
137      // Veto event if photon too close to a lepton
138      if ( any(leptons, deltaRLess(photon, 0.4)) )  vetoEvent;
139
140
141      // Fill histos in single lepton channel, last bin contains overflow events
142      if (_mode == 1 || _mode == 3) {
143        if (jets.size() >= 4 && leptons.size() == 1 && vetoleptons.size() < 2) {
144
145          double mindR_ph_l = 1e3;
146          for (const DressedLepton& lep : leptons) {
147            double dr = deltaR(photon, lep);
148            if (dr <  mindR_ph_l)  mindR_ph_l = dr;
149          }
150          if (mindR_ph_l >= 3.4)  mindR_ph_l = 3.;
151
152          // dRmin(photon,bjet)
153          double mindR_ph_bjet = 1e3;
154          for (const Jet& bjet : bjets) {
155            double dr = deltaR(bjet, photon);
156            if (dr < mindR_ph_bjet) mindR_ph_bjet = dr;
157          }
158          if (mindR_ph_bjet >= 3.4) mindR_ph_bjet = 3.;
159
160          //dRmin(lepton,jet)
161          double mindR_l_jet = 1e3;
162          for (const DressedLepton& lep : leptons) {
163            for (const Jet& jet : jets) {
164              double dr = deltaR(lep, jet);
165              if (dr < mindR_l_jet) mindR_l_jet = dr;
166            }
167          }
168          if (mindR_l_jet >= 3.4) mindR_l_jet = 3.;
169
170          double ph_pt = photon.pT();
171          if (ph_pt >= 500) ph_pt = 450.;
172
173          double jet1_pt = jets[0].pT();
174          if (jet1_pt >= 450.) jet1_pt = 405;
175
176          dualfill("sl_ph_pt",       ph_pt/GeV);
177          dualfill("sl_ph_abseta",   photon.abseta());
178          dualfill("sl_ph_l_dR",     mindR_ph_l);
179          dualfill("sl_ph_bjet_dR",  mindR_ph_bjet);
180          dualfill("sl_l_jet_dR",    mindR_l_jet);
181          dualfill("sl_jet1_pt",     jet1_pt/GeV);
182          dualfill("sldl_ph_pt",     ph_pt/GeV);
183          dualfill("sldl_ph_abseta", photon.abseta());
184
185        }
186      }
187
188
189      // Fill histos in dilepton channel, last bin contains overflow events
190      if ( _mode == 2 || _mode == 3 ) {
191        if ( jets.size() >= 2 && leptons.size() == 2 && vetoleptons.size() < 3) {
192
193          // dRmin(photon, lepton)
194          double mindR_ph_l = 1e3;
195          for (const DressedLepton& lep : leptons) {
196            double dr = deltaR(lep, photon);
197            if (dr < mindR_ph_l)  mindR_ph_l = dr;
198          }
199          if (mindR_ph_l >= 3.4)  mindR_ph_l = 3.;
200
201          // dRmin(photon,bjet)
202          double mindR_ph_bjet = 1e3;
203          for (const Jet& bjet : bjets) {
204            double dr = deltaR(bjet, photon);
205            if (dr < mindR_ph_bjet ) mindR_ph_bjet = dr;
206          }
207          if (mindR_ph_bjet >= 3.4) mindR_ph_bjet = 3.;
208
209          //dRmin(lepton,jet)
210          double mindR_l_jet = 999.0;
211          for (const DressedLepton& lep :leptons) {
212            for (const Jet& jet : jets) {
213             double dr = deltaR(lep, jet);
214             if (dr < mindR_l_jet)  mindR_l_jet = dr;
215            }
216          }
217          if (mindR_l_jet >= 3.4) mindR_l_jet = 3.;
218
219          double ph_pt = photon.pT();
220          if (ph_pt >= 500.) ph_pt = 450.;
221
222          double jet1_pt = jets[0].pT();
223          if (jet1_pt >= 450.) jet1_pt = 405.;
224
225          dualfill("dl_ph_pt",       ph_pt/GeV);
226          dualfill("dl_ph_abseta",   photon.abseta());
227          dualfill("dl_ph_l_dR",     mindR_ph_l);
228          dualfill("dl_ph_bjet_dR",  mindR_ph_bjet);
229          dualfill("dl_l_jet_dR",    mindR_l_jet);
230          dualfill("dl_jet1_pt",     jet1_pt/GeV);
231          dualfill("sldl_ph_pt",     ph_pt/GeV);
232          dualfill("sldl_ph_abseta", photon.abseta());
233
234        }
235      }
236    }
237
238
239
240    void finalize() {
241
242      // Normalise histogram to absolute fiducial cross section or to unity
243      scale(_h, crossSection() / femtobarn / sumOfWeights());
244      normalize(_n, 1.0, true);
245
246    }
247
248
249    void dualbook(const string& name, unsigned int id) {
250      book(_h["abs_"+name],    id, 1, 1);
251      book(_n["norm_"+name], 1+id, 1, 1);
252    }
253
254    void dualfill(const string& name, const double value) {
255      _h["abs_"+name]->fill(value);
256      _n["norm_"+name]->fill(value);
257    }
258
259  private:
260
261    /// @name Histograms
262    map<string, Histo1DPtr> _h, _n;
263
264    size_t _mode;
265
266  };
267
268  // The hook for the plugin system
269  RIVET_DECLARE_PLUGIN(ATLAS_2024_I2768921);
270}