rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

ATLAS_2019_I1734263

WW production at 13 TeV
Experiment: ATLAS (LHC)
Inspire ID: 1734263
Status: VALIDATED
Authors:
  • Kristin Lohwasser
  • Christian Gutschow
  • Hannes Mildner
References: Beams: p+ p+
Beam energies: (6500.0, 6500.0) GeV
Run details:
  • pp -> WW production at 13 TeV

A measurement of fiducial and differential cross-sections for $W^+W^-$ production in proton--proton collisions at $\sqrt{s}=$13 TeV with the ATLAS experiment at the Large Hadron Collider using data corresponding to an integrated luminosity of $36.1$ fb$^{-1}$ is presented. Events with one electron and one muon are selected, corresponding to the decay of the diboson system as $WW\rightarrow e^{\pm}\nu\mu^{\mp}\nu$. To suppress top-quark background, events containing jets with a transverse momentum exceeding 35 GeV are not included in the measurement phase space. The fiducial cross-section, six differential distributions and the cross-section as a function of the jet-veto transverse momentum threshold are measured and compared with several theoretical predictions. Constraints on anomalous electroweak gauge boson self-interactions are also presented in the framework of a dimension-six effective field theory.

Source code: ATLAS_2019_I1734263.cc
  1// -*- C++ -*-
  2#include "Rivet/Analysis.hh"
  3#include "Rivet/Projections/FastJets.hh"
  4#include "Rivet/Projections/FinalState.hh"
  5#include "Rivet/Projections/MissingMomentum.hh"
  6#include "Rivet/Projections/PromptFinalState.hh"
  7#include "Rivet/Projections/LeptonFinder.hh"
  8#include "Rivet/Projections/VetoedFinalState.hh"
  9
 10namespace Rivet {
 11
 12
 13  /// @brief WW production at 13 TeV
 14  class ATLAS_2019_I1734263 : public Analysis {
 15  public:
 16
 17    /// Constructor
 18    RIVET_DEFAULT_ANALYSIS_CTOR(ATLAS_2019_I1734263);
 19
 20
 21    /// @name Analysis methods
 22    /// @{
 23
 24    /// Book histograms and initialise projections before the run
 25    void init() {
 26
 27      const FinalState fs(Cuts::abseta < 5.);
 28
 29      // Project photons for dressing
 30      FinalState photon_id(Cuts::abspid == PID::PHOTON);
 31
 32      // Cuts for leptons
 33      Cut lepton_cuts = (Cuts::abseta < 2.5) && (Cuts::pT > 27*GeV);
 34      // Lepton for simplified phase space (e.g. for comparison with CMS)
 35      Cut lepton_cuts_simpl = (Cuts::abseta < 2.5) && (Cuts::pT > 25*GeV);
 36
 37      // Project dressed leptons (e/mu not from tau) with pT > 27 GeV and |eta| < 2.5
 38      // Both for normal and simplified phase space
 39      PromptFinalState lep_bare(Cuts::abspid == PID::MUON || Cuts::abspid == PID::ELECTRON);
 40      LeptonFinder lep_dressed(lep_bare, photon_id, 0.1, lepton_cuts);
 41      declare(lep_dressed,"lep_dressed");
 42      LeptonFinder lep_dressed_simpl(lep_bare, photon_id, 0.1, lepton_cuts_simpl);
 43      declare(lep_dressed_simpl,"lep_dressed_simpl");
 44
 45      // Get MET
 46      MissingMomentum mm(fs);
 47      declare(mm, "met");
 48
 49      // Define hadrons as everything but dressed leptons (for jet clustering)
 50      VetoedFinalState hadrons(fs);
 51      hadrons.addVetoOnThisFinalState(lep_dressed);
 52      declare(hadrons, "hadrons");
 53      VetoedFinalState hadrons_simpl(fs);
 54      hadrons_simpl.addVetoOnThisFinalState(lep_dressed_simpl);
 55      declare(hadrons_simpl, "hadrons_simpl");
 56
 57      // Project jets
 58      FastJets jets(hadrons, JetAlg::ANTIKT, 0.4, JetMuons::ALL, JetInvisibles::NONE);
 59      declare(jets, "jets");
 60      FastJets jets_simpl(hadrons_simpl, JetAlg::ANTIKT, 0.4, JetMuons::ALL, JetInvisibles::NONE);
 61      declare(jets_simpl, "jets_simpl");
 62
 63
 64      // Book histograms
 65      // fiducial differential cross section as a function of the jet-veto pt cut
 66      book(_d, 1, 1, 1);
 67
 68      // fiducial differential cross section (leading lepton pt)
 69      book(_h["ptlead"], 4, 1, 1);
 70      book(_h["ptlead_norm"], 22, 1, 1);
 71      book(_h["ptlead_simpl"], 41, 1, 1);
 72
 73      // fiducial differential cross section (dilepton-system mll)
 74      book(_h["mll"], 7, 1, 1);
 75      book(_h["mll_norm"], 25, 1, 1);
 76      book(_h["mll_simpl"], 42, 1, 1);
 77
 78      // fiducial differential cross section (dilepton-system ptll)
 79      book(_h["ptll"], 10, 1, 1);
 80      book(_h["ptll_norm"], 28, 1, 1);
 81      book(_h["ptll_simpl"], 43, 1, 1);
 82
 83      // fiducial differential cross section (absolute rapidity of dilepton-system y_ll)
 84      book(_h["yll"], 13, 1, 1);
 85      book(_h["yll_norm"], 31, 1, 1);
 86
 87      // fiducial differential cross section (dilepton-system delta_phi_ll)
 88      book(_h["dphill"], 16, 1, 1);
 89      book(_h["dphill_norm"], 34, 1, 1);
 90
 91
 92      // fiducial differential cross section (absolute costheta* of dilepton-system costhetastar_ll)
 93      book(_h["costhetastarll"], 19, 1, 1);
 94      book(_h["costhetastarll_norm"], 37, 1, 1);
 95
 96    }
 97
 98
 99    /// Perform the per-event analysis
100    void analyze(const Event& event) {
101
102      // Get met and find leptons
103      const MissingMomentum& met = apply<MissingMomentum>(event, "met");
104      const DressedLeptons &leptons       = apply<LeptonFinder>(event, "lep_dressed").dressedLeptons();
105      const DressedLeptons &leptons_simpl = apply<LeptonFinder>(event, "lep_dressed_simpl").dressedLeptons();
106
107      // Find jets and jets for simplified phase space (for the latter slightly different leptons are excluded from clustering)
108      const Jets& jets30     = apply<FastJets>(event, "jets").jetsByPt(Cuts::abseta < 4.5 && Cuts::pT > 30*GeV);
109      const Jets& jets_simpl = apply<FastJets>(event, "jets_simpl").jetsByPt(Cuts::abseta < 4.5 && Cuts::pT > 30*GeV);
110
111      // Define observables
112      const FourMomentum dilep  = leptons.size()>1 ? leptons[0].momentum() + leptons[1].momentum() : FourMomentum(0,0,0,0);
113      const double ptll         = leptons.size()>1 ? dilep.pT()/GeV : -1;
114      const double Mll          = leptons.size()>1 ? dilep.mass()/GeV : -1;
115      const double Yll          = leptons.size()>1 ? dilep.absrap() : -5;
116      const double DPhill       = leptons.size()>1 ? fabs(deltaPhi(leptons[0], leptons[1])) : -1.;
117      const double costhetastar = leptons.size()>1 ? fabs(tanh((leptons[0].eta() - leptons[1].eta()) / 2)) : -0.2;
118
119      // Define observables for simplified PS
120      const FourMomentum dilep_simpl = leptons_simpl.size()>1 ? leptons_simpl[0].momentum() + leptons_simpl[1].momentum() : FourMomentum(0,0,0,0);
121      const double ptll_simpl        = leptons_simpl.size()>1 ? dilep_simpl.pT()/GeV : -1;
122      const double Mll_simpl         = leptons_simpl.size()>1 ? dilep_simpl.mass()/GeV : -1;
123
124      // Cuts for simplified PS
125      bool veto_simpl = false;
126      // Remove events that do not contain 2 good leptons (either muons or electrons)
127      if ( leptons_simpl.size() != 2 ) veto_simpl = true;
128      // Veto same-flavour events
129      else if ( leptons_simpl[0].abspid() == leptons_simpl[1].abspid()) veto_simpl = true;
130      // Veto same-charge events
131      else if ( leptons_simpl[0].pid()*leptons_simpl[1].pid()>0) veto_simpl = true;
132      // MET (pt-MET) cut
133      else if (met.missingPt() <= 20*GeV)  veto_simpl = true;
134      // jetveto
135      else if ( !jets_simpl.empty() ) veto_simpl = true;
136
137      // Fill histos for simplified phase space
138      if ( !veto_simpl ){
139        _h["ptlead_simpl"]->fill(leptons_simpl[0].pT()/GeV);
140        _h["ptll_simpl"]->fill(ptll_simpl);
141        _h["mll_simpl"]->fill(Mll_simpl);
142      }
143
144
145      // Event selection for proper fiducial phase space
146      // Remove events that do not contain 2 good leptons (either muons or electrons)
147      if ( leptons.size() != 2)  vetoEvent;
148      // Veto same-flavour events
149      if ( leptons[0].abspid() == leptons[1].abspid())  vetoEvent;
150      // Veto same-charge events
151      if ( leptons[0].pid()*leptons[1].pid()>0)  vetoEvent;
152      // MET (pt-MET) cut
153      if (met.missingPt() <= 20*GeV)  vetoEvent;
154      // m_ll cut
155      if (dilep.mass() <= 55*GeV)  vetoEvent;
156      // pt_ll cut
157      if (dilep.pT() <= 30*GeV)  vetoEvent;
158
159      // Fill cross section as function of veto-jet pt before applying jet veto
160      if (jets30.empty() || jets30[0].pT()/GeV < 30.) _d->fill(30);
161      if (jets30.empty() || jets30[0].pT()/GeV < 35.) _d->fill(35);
162      if (jets30.empty() || jets30[0].pT()/GeV < 40.) _d->fill(40);
163      if (jets30.empty() || jets30[0].pT()/GeV < 45.) _d->fill(45);
164      if (jets30.empty() || jets30[0].pT()/GeV < 50.) _d->fill(50);
165      if (jets30.empty() || jets30[0].pT()/GeV < 55.) _d->fill(55);
166      if (jets30.empty() || jets30[0].pT()/GeV < 60.) _d->fill(60);
167      // Jet veto at 35 GeV is the default
168      if (!jets30.empty() && jets30[0].pT()/GeV > 35.)  vetoEvent;
169
170      // fill histograms
171      _h["ptlead"]->fill(leptons[0].pT()/GeV);
172      _h["ptlead_norm"]->fill(leptons[0].pT()/GeV);
173
174      _h["ptll"]->fill(ptll);
175      _h["ptll_norm"]->fill(ptll);
176
177      _h["mll"]->fill(Mll);
178      _h["mll_norm"]->fill(Mll);
179
180      _h["yll"]->fill(Yll);
181      _h["yll_norm"]->fill(Yll);
182
183      _h["dphill"]->fill(DPhill);
184      _h["dphill_norm"]->fill(DPhill);
185
186      _h["costhetastarll"]->fill(costhetastar);
187      _h["costhetastarll_norm"]->fill(costhetastar);
188
189    }
190
191
192    /// Normalise histograms etc., after the run
193    void finalize() {
194      const double sf(crossSection()/femtobarn/sumOfWeights());
195      // scale histogram by binwidth, as bin content is actually a integrated fiducial cross section
196      // scale to cross section
197      scale(_d, sf);
198      for (auto& hist : _h) {
199        scale(hist.second, sf);
200        if (hist.first.find("norm") != string::npos) normalize(hist.second);
201      }
202    }
203
204    /// @}
205
206  private:
207
208    /// @name Histograms
209    /// @{
210    map<string, Histo1DPtr> _h;
211    BinnedHistoPtr<int> _d;
212    /// @}
213
214  };
215
216
217  RIVET_DECLARE_PLUGIN(ATLAS_2019_I1734263);
218
219}