rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

ATLAS_2018_I1635273

W + jets production at 8 TeV
Experiment: ATLAS (LHC)
Inspire ID: 1635273
Status: VALIDATED
Authors:
  • Valerie Lang
  • Christian Gutschow
References: Beams: p+ p+
Beam energies: (4000.0, 4000.0) GeV
Run details:
  • inclusive W -> ev production at 8 TeV

This paper presents a measurement of the $W$ boson production cross section and the $W^{+}/W^{-}$ cross-section ratio, both in association with jets, in proton-proton collisions at $\sqrt{s}=8$ TeV with the ATLAS experiment at the Large Hadron Collider. The measurement is performed in final states containing one electron and missing transverse momentum using data corresponding to an integrated luminosity of 20.2 fb$^{-1}$. Differential cross sections for events with at least one or two jets are presented for a range of observables, including jet transverse momenta and rapidities, the scalar sum of transverse momenta of the visible particles and the missing transverse momentum in the event, and the transverse momentum of the $W$ boson. For a subset of the observables, the differential cross sections of positively and negatively charged $W$ bosons are measured separately. In the cross-section ratio of $W^{+}/W^{-}$ the dominant systematic uncertainties cancel out, improving the measurement precision by up to a factor of nine. The observables and ratios selected for this paper provide valuable input for the up quark, down quark, and gluon parton distribution functions of the proton.

Source code: ATLAS_2018_I1635273.cc
  1#include "Rivet/Analysis.hh"
  2#include "Rivet/Projections/FinalState.hh"
  3#include "Rivet/Projections/InvisibleFinalState.hh"
  4#include "Rivet/Projections/FastJets.hh"
  5#include "Rivet/Projections/PromptFinalState.hh"
  6#include "Rivet/Projections/LeptonFinder.hh"
  7
  8namespace Rivet {
  9
 10
 11  /// @brief W + jets production at 8 TeV
 12  class ATLAS_2018_I1635273 : public Analysis {
 13  public:
 14
 15    /// Constructor
 16    RIVET_DEFAULT_ANALYSIS_CTOR(ATLAS_2018_I1635273);
 17
 18
 19    // Book histograms and initialise projections before the run
 20    void init() {
 21
 22      // Get options from the new option system
 23      // Default uses electrons
 24      // EL looks for W->ev candidate
 25      // MU looks for W->mv candidate
 26      _mode = 0;
 27      if ( getOption("LMODE") == "EL" ) _mode = 0;
 28      if ( getOption("LMODE") == "MU" ) _mode = 1;
 29
 30      FinalState fs;
 31      Cut cuts = Cuts::pT > 25*GeV && Cuts::abseta < 2.5;
 32
 33      // Get photons to dress leptons
 34      // (Paper says "radiated photons", but there was
 35      // no promptness requirement in the analysis code)
 36      FinalState photons(Cuts::abspid == PID::PHOTON);
 37
 38      // Get dressed leptons
 39      PromptFinalState leptons(Cuts::abspid == (_mode? PID::MUON : PID::ELECTRON), TauDecaysAs::NONPROMPT);
 40      LeptonFinder dressedleptons(leptons, photons, 0.1, cuts);
 41      declare(dressedleptons, "LeptonFinder");
 42
 43      // Get neutrinos for MET calculation
 44      declare(InvisibleFinalState(OnlyPrompt::YES), "InvFS");
 45
 46      // jets
 47      FastJets jets(fs, JetAlg::ANTIKT, 0.4, JetMuons::NONE);
 48      declare(jets, "Jets");
 49
 50      // book histograms
 51      book(_d["N_incl_pb"],            1, 1, 1);
 52      book(_h["HT_1j_fb"],             6, 1, 1);
 53      book(_h["W_pt_1j_fb"],          11, 1, 1);
 54      book(_h["jet_pt1_1j_fb"],       16, 1, 1);
 55      book(_h["jet_y1_1j_fb"],        21, 1, 1);
 56      book(_h["jet_pt2_2j_fb"],       26, 1, 1);
 57      book(_h["jet_y2_2j_fb"],        28, 1, 1);
 58      book(_h["DeltaRj12_2j_fb"],     30, 1, 1);
 59      book(_h["jet_mass12_2j_fb"],    32, 1, 1);
 60      book(_d["N_pb"],                34, 1, 1);
 61      book(_h["HT_2j_fb"],            36, 1, 1);
 62      book(_h["W_pt_2j_fb"],          41, 1, 1);
 63      book(_h["jet_pt1_2j_fb"],       46, 1, 1);
 64      book(_h["el_eta_0j_pb"],        51, 1, 1);
 65      book(_h["el_eta_1j_pb"],        56, 1, 1);
 66
 67      book(_d["Wplus_N_incl_pb"],      3, 1, 1);
 68      book(_h["Wplus_HT_1j_fb"],       8, 1, 1);
 69      book(_h["Wplus_W_pt_1j_fb"],    13, 1, 1);
 70      book(_h["Wplus_jet_pt1_1j_fb"], 18, 1, 1);
 71      book(_h["Wplus_jet_y1_1j_fb"],  23, 1, 1);
 72      book(_h["Wplus_HT_2j_fb"],      38, 1, 1);
 73      book(_h["Wplus_W_pt_2j_fb"],    43, 1, 1);
 74      book(_h["Wplus_jet_pt1_2j_fb"], 48, 1, 1);
 75      book(_h["Wplus_el_eta_0j_pb"],  53, 1, 1);
 76      book(_h["Wplus_el_eta_1j_pb"],  58, 1, 1);
 77
 78      book(_d["Wminus_N_incl_pb"],     3, 1, 2);
 79      book(_h["Wminus_HT_1j_fb"],      8, 1, 2);
 80      book(_h["Wminus_W_pt_1j_fb"],   13, 1, 2);
 81      book(_h["Wminus_jet_pt1_1j_fb"],18, 1, 2);
 82      book(_h["Wminus_jet_y1_1j_fb"], 23, 1, 2);
 83      book(_h["Wminus_HT_2j_fb"],     38, 1, 2);
 84      book(_h["Wminus_W_pt_2j_fb"],   43, 1, 2);
 85      book(_h["Wminus_jet_pt1_2j_fb"],48, 1, 2);
 86      book(_h["Wminus_el_eta_0j_pb"], 53, 1, 2);
 87      book(_h["Wminus_el_eta_1j_pb"], 58, 1, 2);
 88
 89      // and ratios
 90      book(_r_disc, 3, 1, 3);
 91      book(_r["WplusOverWminus_HT_1j_fb"],       8, 1, 3);
 92      book(_r["WplusOverWminus_W_pt_1j_fb"],    13, 1, 3);
 93      book(_r["WplusOverWminus_jet_pt1_1j_fb"], 18, 1, 3);
 94      book(_r["WplusOverWminus_jet_y1_1j_fb"],  23, 1, 3);
 95      book(_r["WplusOverWminus_HT_2j_fb"],      38, 1, 3);
 96      book(_r["WplusOverWminus_W_pt_2j_fb"],    43, 1, 3);
 97      book(_r["WplusOverWminus_jet_pt1_2j_fb"], 48, 1, 3);
 98      book(_r["WplusOverWminus_el_eta_0j_pb"],  53, 1, 3);
 99      book(_r["WplusOverWminus_el_eta_1j_pb"],  58, 1, 3);
100
101    }
102
103
104    // Perform the per-event analysis
105    void analyze(const Event& event) {
106
107      // retrieve the dressed electrons
108      const Particles& signal_leptons = apply<LeptonFinder>(event, "LeptonFinder").particlesByPt();
109      if (signal_leptons.size() != 1 ) vetoEvent;
110      const Particle& lepton = signal_leptons[0];
111
112      // calulate MET and mT
113      const Particles& invisibles = apply<InvisibleFinalState>(event, "InvFS").particles();
114      FourMomentum pMET = sum(invisibles, Kin::mom, FourMomentum()).setZ(0);
115      const double MET = pMET.pT() / GeV;
116      const double mT = sqrt( 2 * lepton.Et() / GeV * MET * (1 - cos(deltaPhi(lepton,pMET))));
117      if ( MET <= 25. ) vetoEvent;
118      if ( mT  <= 40. ) vetoEvent;
119
120      // retrieve jets
121      Jets jets = apply<FastJets>(event, "Jets").jetsByPt(Cuts::pT > 30*GeV && Cuts::absrap < 4.4);
122
123      idiscardIfAnyDeltaRLess(jets, signal_leptons, 0.2);
124
125      // apply event selection on dR
126      for (const Jet& j : jets) {
127          if (deltaR(j, lepton) < 0.4) vetoEvent;
128      }
129
130      // calculate the observables
131      const double w_pt = (lepton.momentum() + pMET).pT() / GeV;
132      const size_t njets = jets.size();
133      double ST = sum(jets, Kin::pT, 0.0); // scalar pT sum of all selected jets
134      const double HT = ST + lepton.pT() / GeV + MET; //missET;
135
136      // fill W histograms
137      _d["N_pb"]->fill(njets);
138      for (size_t i = 0; i <= njets; ++i) {
139          _d["N_incl_pb"]->fill(i);
140      }
141      _h["el_eta_0j_pb"]->fill(lepton.abseta());
142
143      if (njets > 0) {
144          _h["HT_1j_fb"]->fill(HT);
145          _h["W_pt_1j_fb"]->fill(w_pt);
146          _h["jet_pt1_1j_fb"]->fill(jets[0].pT()/GeV);
147          _h["jet_y1_1j_fb"]->fill(jets[0].absrap());
148          _h["el_eta_1j_pb"]->fill(lepton.abseta());
149      }
150      if (njets > 1) {
151          _h["HT_2j_fb"]->fill(HT);
152          _h["W_pt_2j_fb"]->fill(w_pt);
153          _h["jet_pt1_2j_fb"]->fill(jets[0].pT()/GeV);
154          _h["DeltaRj12_2j_fb"]->fill(deltaR(jets[0],jets[1]));
155          _h["jet_pt2_2j_fb"]->fill(jets[1].pT()/GeV);
156          _h["jet_y2_2j_fb"]->fill(jets[1].absrap());
157          _h["jet_mass12_2j_fb"]->fill( (jets[0].mom()+jets[1].mom()).mass()/GeV);
158      }
159      // fill W+ histograms
160      if (lepton.charge() > 0) {
161          for (size_t i = 0; i <= njets; ++i) {
162            _d["Wplus_N_incl_pb"]->fill(i);
163          }
164          _h["Wplus_el_eta_0j_pb"]->fill(lepton.abseta());
165          if (njets > 0) {
166              _h["Wplus_HT_1j_fb"]->fill(HT);
167              _h["Wplus_W_pt_1j_fb"]->fill(w_pt);
168              _h["Wplus_jet_pt1_1j_fb"]->fill(jets[0].pT()/GeV);
169              _h["Wplus_jet_y1_1j_fb"]->fill(jets[0].absrap());
170              _h["Wplus_el_eta_1j_pb"]->fill(lepton.abseta());
171              if (njets > 1) {
172                  _h["Wplus_HT_2j_fb"]->fill(HT);
173                  _h["Wplus_W_pt_2j_fb"]->fill(w_pt);
174                  _h["Wplus_jet_pt1_2j_fb"]->fill(jets[0].pT()/GeV);
175              }
176          }
177      }
178      // fill W- histograms
179      if (lepton.charge() < 0) {
180        for (size_t i = 0; i <= njets; ++i) {
181          _d["Wminus_N_incl_pb"]->fill(i);
182        }
183        _h["Wminus_el_eta_0j_pb"]->fill(lepton.abseta());
184        if (njets > 0) {
185          _h["Wminus_HT_1j_fb"]->fill(HT);
186          _h["Wminus_W_pt_1j_fb"]->fill(w_pt);
187          _h["Wminus_jet_pt1_1j_fb"]->fill(jets[0].pT()/GeV);
188          _h["Wminus_jet_y1_1j_fb"]->fill(jets[0].absrap());
189          _h["Wminus_el_eta_1j_pb"]->fill(lepton.abseta());
190          if (njets > 1) {
191            _h["Wminus_HT_2j_fb"]->fill(HT);
192            _h["Wminus_W_pt_2j_fb"]->fill(w_pt);
193            _h["Wminus_jet_pt1_2j_fb"]->fill(jets[0].pT()/GeV);
194          }
195        }
196      }
197    }
198
199
200
201    void finalize() {
202      const double scalefactor_fb = crossSection() / sumOfWeights() / femtobarn;
203      const double scalefactor_pb = crossSection() / sumOfWeights() / picobarn;
204
205      scale(_d, scalefactor_pb);
206      for (auto& hit : _h){
207        if (hit.first.find("_fb") != string::npos) scale(hit.second, scalefactor_fb);
208        else                                       scale(hit.second, scalefactor_pb);
209      }
210      for (auto& rit : _r) {
211        string ratio_label = "WplusOverWminus";
212        string num_name   = rit.first;
213        string denom_name = rit.first;
214        num_name.replace(rit.first.find(ratio_label),ratio_label.length(),"Wplus");
215        denom_name.replace(rit.first.find(ratio_label),ratio_label.length(),"Wminus");
216        divide(_h[num_name], _h[denom_name], rit.second);
217      }
218      divide(_d["Wplus_N_incl_pb"], _d["Wminus_N_incl_pb"], _r_disc);
219    }
220
221    protected:
222      // Data members like post-cuts event weight counters go here
223      size_t _mode;
224
225    private:
226      map<string, BinnedHistoPtr<int>> _d;
227      map<string, Histo1DPtr> _h;
228      map<string, Estimate1DPtr> _r;
229      BinnedEstimatePtr<int> _r_disc;
230  };
231
232  RIVET_DECLARE_PLUGIN(ATLAS_2018_I1635273);
233}