rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

ATLAS_2016_I1448301

$Z \gamma (\gamma)$ cross sections at 8 TeV
Experiment: ATLAS (LHC)
Inspire ID: 1448301
Status: VALIDATED
Authors:
  • Hulin Wang
  • Evgeny Soldatov
  • Christian Gutschow
References: Beams: p+ p+
Beam energies: (4000.0, 4000.0) GeV
Run details:
  • p + p -> e+ e- gamma (gamma) at 8 TeV

The production of $Z$ bosons with one or two isolated high-energy photons is studied using pp collisions at $\sqrt{s}=8$ TeV. The analyses use a data sample with an integrated luminosity of 20.3 fb${}^{-1}$ collected by the ATLAS detector during the 2012 LHC data taking. The $Z\gamma$ and $Z\gamma\gamma$ production cross sections are measured with leptonic ($e^+e^-$, $\mu^+\mu^-$, $\nu\bar{\nu}$) decays of the Z boson, in extended fiducial regions defined in terms of the lepton and photon acceptance. They are then compared to cross-section predictions from the Standard Model, where the sources of the photons are radiation off initial-state quarks and radiative Z-boson decay to charged leptons, and from fragmentation of final-state quarks and gluons into photons. The yields of events with photon transverse energy $E_\text{T}>250$ GeV from $\ell^+\ell^-\gamma$ events and with $E_\text{T}>400$ GeV from $\nu\bar{\nu}\gamma$ events are used to search for anomalous triple gauge-boson couplings $ZZ\gamma$ and $Z\gamma\gamma$. The yields of events with diphoton invariant mass $m_{\gamma\gamma}>200$ GeV from $\ell^+\ell^-\gamma\gamma$ events and with $m_{\gamma\gamma}>300$ GeV from $\nu\bar{\nu}\gamma\gamma$ events are used to search for anomalous quartic gauge-boson couplings $ZZ\gamma\gamma$ and $Z\gamma\gamma\gamma$. No deviations from Standard Model predictions are observed and limits are placed on parameters used to describe anomalous triple and quartic gauge-boson couplings. The lepton channel can be specified by using the relevant LMODE (EL, MU, NU, LL). The default plugin will fill all of them with any events passing the cuts. NU will do only the MET channels, EL & MU will do only the electron or muon and will populate the combined plots assuming charged lepton universality. LL with do on the electron and muon but will populate the combined plots with the average of the two.

Source code: ATLAS_2016_I1448301.cc
  1// -*- C++ -*-
  2#include "Rivet/Analysis.hh"
  3#include "Rivet/Projections/FinalState.hh"
  4#include "Rivet/Projections/FastJets.hh"
  5#include "Rivet/Projections/VetoedFinalState.hh"
  6#include "Rivet/Projections/IdentifiedFinalState.hh"
  7#include "Rivet/Projections/LeptonFinder.hh"
  8#include "Rivet/Projections/PromptFinalState.hh"
  9#include "Rivet/Projections/VisibleFinalState.hh"
 10
 11namespace Rivet {
 12
 13
 14  /// @brief Z/gamma cross section measurement at 8 TeV
 15  class ATLAS_2016_I1448301 : public Analysis {
 16  public:
 17
 18    /// Constructor
 19    RIVET_DEFAULT_ANALYSIS_CTOR(ATLAS_2016_I1448301);
 20
 21
 22    /// @name Analysis methods
 23    /// @{
 24
 25    /// Book histograms and initialise projections before the run
 26    void init() {
 27
 28      // Get options from the new option system
 29      // Default tries to fill everything
 30      // NU fills only the MET plots
 31      // EL fills only the Electron plots and combined plots assuming lepton univerality
 32      // MU fills only the Muon plots and combined plots assuming lepton univerality
 33      // LL fills electron and Muon plots and combined plots from correct average
 34      _mode = 0;
 35      if ( getOption("LMODE") == "NU" ) _mode = 1;
 36      if ( getOption("LMODE") == "EL" ) _mode = 2;
 37      if ( getOption("LMODE") == "MU" ) _mode = 3;
 38      if ( getOption("LMODE") == "LL" ) _mode = 4;
 39
 40      // Prompt photons
 41      const Cut photoncut = Cuts::abspid == PID::PHOTON && Cuts::pT > 15*GeV && Cuts::abseta < 2.37;
 42      PromptFinalState photon_fs(photoncut);
 43      declare(photon_fs, "Photons");
 44
 45      // Prompt leptons
 46      const PromptFinalState bareelectron_fs = Cuts::abspid == PID::ELECTRON;
 47      const PromptFinalState baremuon_fs = Cuts::abspid == PID::MUON;
 48
 49      // Dressed leptons
 50      const IdentifiedFinalState allphoton_fs(PID::PHOTON); // photons used for lepton dressing
 51      const Cut leptoncut = Cuts::pT > 25*GeV && Cuts::abseta < 2.47;
 52      const LeptonFinder dressedelectron_fs(bareelectron_fs, allphoton_fs, 0.1, leptoncut);
 53      const LeptonFinder dressedmuon_fs(baremuon_fs, allphoton_fs, 0.1, leptoncut);
 54
 55      declare(dressedelectron_fs, "Electrons");
 56      declare(dressedmuon_fs, "Muons");
 57
 58      // MET (prompt neutrinos)
 59      VetoedFinalState ivfs;
 60      ivfs.addVetoOnThisFinalState(VisibleFinalState());
 61      declare(PromptFinalState(ivfs), "MET");
 62
 63      // Jets
 64      VetoedFinalState jet_fs;
 65      jet_fs.vetoNeutrinos();
 66      jet_fs.addVetoPairId(PID::MUON);
 67      const FastJets fastjets(jet_fs, JetAlg::ANTIKT, 0.4);
 68      declare(fastjets, "Jets");
 69
 70      // Histograms
 71
 72      // MET
 73      if (_mode == 0 || _mode == 1){
 74        book(_d["vvg"],     2, 1, 1);
 75        book(_d["vvgg"],    4, 1, 1);
 76        book(_h["pT"],      7, 1, 1);
 77        book(_h["pT_0jet"], 8, 1, 1);
 78      }
 79
 80      // always book e and mu in charged lepton modes; there are sometimes 4 leptons.
 81      if (_mode != 1){
 82        // electron
 83        book(_d["eeg"],  1, 1, 1);
 84        book(_d["eegg"], 3, 1, 1);
 85        // muon
 86        book(_d["mmg"],  1, 1, 2);
 87        book(_d["mmgg"], 3, 1, 2);
 88
 89        // combined
 90        book(_d["llgg"], 3, 1, 3);
 91        book(_d["llg"], 1, 1, 3);
 92        book(_h["pT"], 5, 1, 1);
 93        book(_h["pT_0jet"], 6, 1, 1);
 94        book(_h["M"], 9, 1, 1);
 95        book(_h["M_0jet"], 10, 1, 1);
 96        book(_d["Njets"], 11, 1, 1);
 97      }
 98    }
 99
100
101    /// Perform the per-event analysis
102    void analyze(const Event& event) {
103      // Get objects
104      DressedLeptons electrons = apply<LeptonFinder>(event, "Electrons").dressedLeptons();
105      DressedLeptons muons = apply<LeptonFinder>(event, "Muons").dressedLeptons();
106
107      const Particles& photons = apply<PromptFinalState>(event, "Photons").particlesByPt();
108      const Jets jets = apply<FastJets>(event, "Jets").jetsByPt();
109
110      const FinalState& metfs = apply<PromptFinalState>(event, "MET");
111      Vector3 met_vec;
112      for (const Particle& p : metfs.particles()) met_vec += p.mom().perpVec();
113
114      if (met_vec.mod() >= 100*GeV && !photons.empty() && _mode < 2) {
115
116        if (photons.size() > 1) { // nu nu y y
117
118          bool yy_veto = false;
119          yy_veto |= photons[0].pT() < 22*GeV;
120          yy_veto |= photons[1].pT() < 22*GeV;
121          yy_veto |= met_vec.mod() < 110*GeV;
122          const double yyPhi = (photons[0].momentum() + photons[1].momentum()).phi();
123          yy_veto |= fabs(yyPhi - met_vec.phi()) < 2.62 || fabs(yyPhi - met_vec.phi()) > 3.66;
124          yy_veto |= deltaR(photons[0], photons[1]) < 0.4;
125
126          // Photon isolation calculated by jets, count jets
127          Jet ph0_jet, ph1_jet;
128          double min_dR_ph0_jet = 999., min_dR_ph1_jet = 999.;
129          size_t njets = 0;
130          for (const Jet& j : jets) {
131            if (j.pT() > 30*GeV && j.abseta() < 4.5) {
132              if (deltaR(j, photons[0]) > 0.3 && deltaR(j, photons[1]) > 0.3)  ++njets;
133            }
134            if (deltaR(j, photons[0]) < min_dR_ph0_jet) {
135              min_dR_ph0_jet = deltaR(j, photons[0]);
136              ph0_jet = j;
137            }
138            if (deltaR(j, photons[1]) < min_dR_ph1_jet) {
139              min_dR_ph1_jet = deltaR(j, photons[1]);
140              ph1_jet = j;
141            }
142          }
143          double photon0iso = 0., photon1iso = 0.;
144          if (min_dR_ph0_jet < 0.4)  photon0iso = ph0_jet.pT() - photons[0].pT();
145          if (min_dR_ph1_jet < 0.4)  photon1iso = ph1_jet.pT() - photons[1].pT();
146          yy_veto |= photon0iso/photons[0].pT() > 0.5;
147          yy_veto |= photon1iso/photons[1].pT() > 0.5;
148
149          if (!yy_veto) {
150            _d["vvgg"]->fill(">= 0.0"s);
151            if (!njets)  _d["vvgg"]->fill("= 0"s);
152          }
153        } // end of nu nu y y section
154
155
156        if ((photons[0].pT() >= 130*GeV)  &&
157            (fabs(fabs(deltaPhi(photons[0], met_vec)) - 3.14) <= 1.57)) {
158
159          // Photon isolation calculated by jets, count jets
160          Jet ph_jet;
161          double min_dR_ph_jet = 999.;
162          size_t njets = 0;
163          for (const Jet& j : jets) {
164            if (j.pT() > 30*GeV && j.abseta() < 4.5) {
165              if (deltaR(j, photons[0]) > 0.3)  ++njets;
166            }
167            if (deltaR(j, photons[0]) < min_dR_ph_jet) {
168              min_dR_ph_jet = deltaR(j, photons[0]);
169              ph_jet = j;
170            }
171          }
172          double photoniso = 0;
173          if (min_dR_ph_jet < 0.4)  photoniso = ph_jet.pT() - photons[0].pT();
174          if (photoniso/photons[0].pT() > 0.5)  vetoEvent;
175
176          const double pTgamma = photons[0].pT()/GeV;
177          _h["pT"]->fill(pTgamma);
178          _d["vvg"]->fill(">= 0.0"s);
179          if (!njets) {
180            _d["vvg"]->fill("= 0"s);
181            _h["pT_0jet"]->fill(pTgamma);
182          }
183
184        }
185      }  // end of nu nu y (y) section
186
187      // Dilepton candidate
188      bool el = false;
189
190      if (_mode != 1) {
191	  if (_sedges.empty())  _sedges = _d["Njets"]->xEdges();
192      }
193
194      if ( (_mode != 1) &&
195           (( electrons.size() >= 2 && _mode != 3 ) ||
196            ( muons.size()     >= 2 && _mode != 2 ) )) {
197
198        DressedLeptons lep_p, lep_m;
199
200        // Sort the dressed leptons by pt
201        if (electrons.size() >= 2) {
202          el = true;
203          isortByPt(electrons);
204          for (const DressedLepton& lep : electrons) {
205            if (lep.charge() > 0.)  lep_p.push_back(lep);
206            if (lep.charge() < 0.)  lep_m.push_back(lep);
207          }
208        } else {
209          isortByPt(muons);
210          for (const DressedLepton& lep : muons) {
211            if (lep.charge() > 0.)  lep_p.push_back(lep);
212            if (lep.charge() < 0.)  lep_m.push_back(lep);
213          }
214        }
215
216
217        if (!lep_p.empty() && !lep_m.empty() &&
218            (lep_p[0].abspid() == lep_m[0].abspid()) &&
219            ((lep_p[0].momentum() + lep_m[0].momentum()).mass() >= 40*GeV)){
220
221          // Photon lepton overlap removal
222          if (photons.empty())  vetoEvent;
223
224          if (photons.size() > 1) {
225
226            bool veto = false;
227            veto |= deltaR(photons[0], lep_p[0]) < 0.4;
228            veto |= deltaR(photons[0], lep_m[0]) < 0.4;
229            veto |= deltaR(photons[1], lep_p[0]) < 0.4;
230            veto |= deltaR(photons[1], lep_m[0]) < 0.4;
231            veto |= deltaR(photons[0], photons[1]) < 0.4;
232
233            Jet ph0_jet, ph1_jet;
234            double min_dR_ph0_jet = 999., min_dR_ph1_jet=999.;
235            int njets = 0;
236            for (const Jet& j : jets){
237              if (j.pT() > 30*GeV && j.abseta() < 4.5) {
238                if (deltaR(j, lep_p[0]) > 0.3 && deltaR(j, lep_m[0]) > 0.3) {
239                  if (deltaR(j, photons[0]) > 0.3 && deltaR(j, photons[1]) > 0.3 )  ++njets;
240                }
241              }
242              if (deltaR(j, photons[0]) < min_dR_ph0_jet) {
243                min_dR_ph0_jet = deltaR(j, photons[0]);
244                ph0_jet = j;
245              }
246              if (deltaR(j, photons[1]) < min_dR_ph1_jet) {
247                min_dR_ph1_jet = deltaR(j, photons[1]);
248                ph1_jet = j;
249              }
250            }
251            double photon0iso = 0, photon1iso = 0;
252            if (min_dR_ph0_jet < 0.4) photon0iso = ph0_jet.pT() - photons[0].pT();
253            if (min_dR_ph1_jet < 0.4) photon1iso = ph1_jet.pT() - photons[1].pT();
254            veto |= photon0iso/photons[0].pT() > 0.5;
255            veto |= photon1iso/photons[1].pT() > 0.5;
256
257            // Fill plots
258            // ee and mm need doing.
259            if (!veto) {
260              _d["llgg"]->fill(">= 0.0"s);
261              if (el) {
262                _d["eegg"]->fill(">= 0.0"s);
263              } else {
264                _d["mmgg"]->fill(">= 0.0"s);
265              }
266
267              if (!njets) {
268                _d["llgg"]->fill("= 0"s);
269                if (el) {
270                  _d["eegg"]->fill("= 0"s);
271                } else {
272                  _d["mmgg"]->fill("= 0"s);
273                }
274              }
275            }
276          }
277
278          if (deltaR(photons[0], lep_p[0]) < 0.7)  vetoEvent;
279          if (deltaR(photons[0], lep_m[0]) < 0.7)  vetoEvent;
280
281          // Photon isolation calculated by jets, count jets
282          Jet ph_jet;
283          double min_dR_ph_jet = 999.;
284          size_t njets = 0;
285          for (const Jet& j : jets) {
286            if (j.pT() > 30*GeV && j.abseta() < 4.5) {
287              if (deltaR(j, lep_p[0]) > 0.3 && deltaR(j, lep_m[0]) > 0.3 && deltaR(j, photons[0]) > 0.3)  ++njets;
288            }
289            if (deltaR(j, photons[0]) < min_dR_ph_jet) {
290              min_dR_ph_jet = deltaR(j, photons[0]);
291              ph_jet = j;
292            }
293          }
294
295          double photoniso = 0;
296          if (min_dR_ph_jet < 0.4)  photoniso = ph_jet.pT() - photons[0].pT();
297          if (photoniso/photons[0].pT() > 0.5)  vetoEvent;
298
299
300          // Fill plots
301          const double pTgamma = photons[0].pT()/GeV;
302          const double mllgamma = (lep_p[0].momentum() + lep_m[0].momentum() + photons[0].momentum()).mass()/GeV;
303
304          _h["pT"]->fill(pTgamma);
305          _h["M"]->fill(mllgamma);
306          _d["Njets"]->fill(_sedges[min(3, njets)]);
307
308          _d["llg"]->fill(">= 0.0"s);
309          if (el) {
310            _d["eeg"]->fill(">= 0.0"s);
311          } else {
312            _d["mmg"]->fill(">= 0.0"s);
313          }
314
315          if (!njets) {
316            _h["pT_0jet"]->fill(pTgamma);
317            _h["M_0jet"]->fill(mllgamma);
318            _d["llg"]->fill("= 0"s);
319            if (el) {
320              _d["eeg"]->fill("= 0"s);
321            } else {
322              _d["mmg"]->fill("= 0"s);
323            }
324          }
325        }
326      }
327    } // end of analysis
328
329
330    /// Normalise histograms etc., after the run
331    void finalize() {
332      const double sf = crossSection()/femtobarn/sumOfWeights();
333      scale(_h, sf);
334      scale(_d, sf);
335      // if we are running both e and mu, the combined lepton histos
336      // need to be divided by two to get the average
337      if (_mode == 0 || _mode == 4){
338        scale(_d["llgg"], 0.5);
339        scale(_d["llg"], 0.5);
340        scale(_h["pT"], 0.5);
341        scale(_h["pT_0jet"], 0.5);
342        scale(_h["M"], 0.5);
343        scale(_h["M_0jet"], 0.5);
344        scale(_d["Njets"], 0.5);
345      }
346    }
347
348    /// @}
349
350
351  protected:
352
353    // Data members like post-cuts event weight counters go here
354    size_t _mode;
355
356  private:
357
358    /// Histograms
359    map<string, Histo1DPtr> _h;
360    map<string, BinnedHistoPtr<string>> _d;
361    vector<string> _sedges;
362  };
363
364
365  RIVET_DECLARE_PLUGIN(ATLAS_2016_I1448301);
366
367}