Rivet Analyses Reference

ATLAS_2018_I1698006

Z(vv)+gamma measurement at 13 TeV
Experiment: ATLAS (LHC)
Inspire ID: 1698006
Status: VALIDATED
Authors:
  • J.Butterworth@ucl.ac.uk
  • sokratis.michail.19@ucl.ac.uk
  • yoran.yeh.20@ucl.ac.uk
References:
  • JHEP 12 (2018) 010
  • 10.1007/JHEP12(2018)010
  • arXiv: 1810.04995
Beams: p+ p+
Beam energies: (6500.0, 6500.0) GeV
Run details:
  • p + p -> Z0 nu nubar gamma at 13 TeV

The production of Z bosons in association with a high-energy photon (Z$\gamma$ production) is studied in the neutrino decay channel of the Z boson using pp collisions at $\sqrt s$ = 13 TeV. The analysis uses a data sample with an integrated luminosity of 36.1 fb−1 collected by the ATLAS detector at the LHC in 2015 and 2016. Candidate Z$\gamma$ events with invisible decays of the Z boson are selected by requiring significant transverse momentum (pT) of the dineutrino system in conjunction with a single isolated photon with large transverse energy (ET). The rate of Z$\gamma$ production is measured as a function of photon ET, dineutrino system pT and jet multiplicity. Evidence of anomalous triple gauge-boson couplings is sought in Z$\gamma$ production with photo ET greater than 600 GeV. No excess is observed relative to the Standard Model expectation, and upper limits are set on the strength of ZZ$\gamma$ and Z$\gamma\gamma$ couplings.

Source code: ATLAS_2018_I1698006.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
#include "Rivet/Analysis.hh"
#include "Rivet/Projections/FinalState.hh"
#include "Rivet/Projections/FastJets.hh"
#include "Rivet/Projections/PromptFinalState.hh"
#include "Rivet/Projections/InvisibleFinalState.hh"
#include "Rivet/Projections/VetoedFinalState.hh"
#include "Rivet/Projections/DressedLeptons.hh"



namespace Rivet {
  
  /// @brief ATLAS pTmiss+gamma measurement at 13 TeV 
  class ATLAS_2018_I1698006 : public Analysis {
  public:
    
    /// Default constructor
    RIVET_DEFAULT_ANALYSIS_CTOR(ATLAS_2018_I1698006);
    
    /// @name Analysis methods
    //@{
    void init() {

      // Get options
      // Default (OFF) uses the extended phase space of the paper.
      // LVETO ON will add an additonal lepton veto to reject non-Zgamma events
      // (for conservative BSM limits, for example).
      _mode = 0;
      if ( getOption("LVETO") == "ON" ) _mode = 1;

      //prompt photons
      const Cut photoncut = Cuts::abspid == PID::PHOTON && Cuts::Et > 150*GeV && Cuts::abseta < 2.37;
      const PromptFinalState photon_fs(photoncut);
      declare(photon_fs, "Photons");

      //missing energy (prompt neutrinos)
      declare(InvisibleFinalState(true), "MET");
      
      if (_mode==1) {
	FinalState allLeps(Cuts::abspid == PID::ELECTRON || Cuts::abspid == PID::MUON);
	FinalState photons(Cuts::abspid == PID::PHOTON);
	PromptFinalState promptLeps(allLeps);
	Cut dressedLep_cuts = (Cuts::abseta < 2.7) && (Cuts::pT > 7*GeV);
	const DressedLeptons dressedLeps(photons, promptLeps, 0.1, dressedLep_cuts, true);
	declare(dressedLeps, "dressedLeptons");
      }
      
      //jets. run the jet finder on a final state without the prompt photons, and without neutrinos or muons
      VetoedFinalState jet_fs(Cuts::abseta > 4.5);
      jet_fs.addVetoOnThisFinalState(photon_fs);
      FastJets fastjets(jet_fs, FastJets::ANTIKT, 0.4, JetAlg::Muons::NONE, JetAlg::Invisibles::NONE);
      declare(fastjets, "Jets");


      //books histograms
      //fig.4a
      book(_h["Et_inc"],2,1,1);
      //fig.4b
      book(_h["Et_exc"],3,1,1);
      //fig.5a
      book(_h["pT_inc"],4,1,1);
      //fig.5b
      book(_h["pT_exc"],5,1,1);
      //fig.6
      book(_h["Njets"],6,1,1);

    }
    
    void analyze(const Event& event) {


      const Particles& photons = apply<PromptFinalState>(event,"Photons").particlesByPt();
      const Jets& jets = apply<FastJets>(event,"Jets").jetsByPt(Cuts::pT > 50*GeV);
      const FinalState& metfs = apply<InvisibleFinalState>(event,"MET");
      Vector3 met_vec;
      for (const Particle& p : metfs.particles()) met_vec += p.mom().perpVec();


      if (_mode==1) {
	const vector<DressedLepton> &dressedLeptons = apply<DressedLeptons>(event, "dressedLeptons").dressedLeptons();
	if (dressedLeptons.size() > 0) vetoEvent;
      }
      
      
      //Nγ==1 and Emiss > 150 GeV
      if (met_vec.mod() > 150*GeV && photons.size()==1){

	  //inclusive case (Njet>=0)
	  if (jets.size()>=0){
	    bool dR_veto = any(jets, DeltaRLess(photons[0], 0.3));
	    if (not dR_veto) {
	      double Et_photon = photons[0].Et()/GeV;
	      _h["Et_inc"]->fill(Et_photon);
	      //fill in missing energy (neutrino) pT histogram (inclusive)
	      _h["pT_inc"]->fill(met_vec.mod()/GeV);
	    }
	  }

	  //exclusive case (Njet==0)
	  if (jets.size() == 0){
	    double Et_photon = photons[0].Et()/GeV;
	    _h["Et_exc"]->fill(Et_photon);
	     //fill in missing energy (neutrino) pT histogram (exclusive)
	    _h["pT_exc"]->fill(met_vec.mod()/GeV);
	  }

	  _h["Njets"]->fill(jets.size());
	  
      }
      
    }

    void finalize() {
      
      const double sf = crossSection()/femtobarn/sumOfWeights();
      scale(_h, sf);

    }

    //@}

    private:
    map<string, Histo1DPtr> _h;
    size_t _mode;
    
  };

  // Magic required by the plugin system 
  RIVET_DECLARE_PLUGIN(ATLAS_2018_I1698006);
  
}