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
  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
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
#include "Rivet/Analysis.hh"
#include "Rivet/Projections/FinalState.hh"
#include "Rivet/Projections/InvisibleFinalState.hh"
#include "Rivet/Projections/FastJets.hh"
#include "Rivet/Projections/PromptFinalState.hh"
#include "Rivet/Projections/DressedLeptons.hh"

namespace Rivet {

  /// @brief W + jets production at 8 TeV
  class ATLAS_2018_I1635273 : public Analysis {
  public:

    RIVET_DEFAULT_ANALYSIS_CTOR(ATLAS_2018_I1635273);

	  // Book histograms and initialise projections before the run
    void init() {

      // Get options from the new option system
      // Default uses electrons
      // EL looks for W->ev candidate
      // MU looks for W->mv candidate
      _mode = 0;
      if ( getOption("LMODE") == "EL" ) _mode = 0;
      if ( getOption("LMODE") == "MU" ) _mode = 1;

      FinalState fs;
      Cut cuts = Cuts::pT > 25*GeV && Cuts::abseta < 2.5;

      // Get photons to dress leptons
      // (Paper says "radiated photons", but there was
      // no promptness requirement in the analysis code)
      FinalState photons(Cuts::abspid == PID::PHOTON);

      // Get dressed leptons
      PromptFinalState leptons(Cuts::abspid == (_mode? PID::MUON : PID::ELECTRON), false);
      DressedLeptons dressedleptons(photons, leptons, 0.1, cuts, true);
      declare(dressedleptons, "DressedLeptons");

      // Get neutrinos for MET calculation
      declare(InvisibleFinalState(true), "InvFS"); // true = only allow prompt invisibles

      // jets
      FastJets jets(fs, FastJets::ANTIKT, 0.4, JetAlg::Muons::NONE);
      declare(jets, "Jets");

      // book histograms
      book(_h["N_incl_pb"],            1, 1, 1);
      book(_h["HT_1j_fb"],             6, 1, 1);
      book(_h["W_pt_1j_fb"],          11, 1, 1);
      book(_h["jet_pt1_1j_fb"],       16, 1, 1); 
      book(_h["jet_y1_1j_fb"],        21, 1, 1);
      book(_h["jet_pt2_2j_fb"],       26, 1, 1);
      book(_h["jet_y2_2j_fb"],        28, 1, 1);
      book(_h["DeltaRj12_2j_fb"],     30, 1, 1);
      book(_h["jet_mass12_2j_fb"],    32, 1, 1);
      book(_h["N_pb"],                34, 1, 1);
      book(_h["HT_2j_fb"],            36, 1, 1);
      book(_h["W_pt_2j_fb"],          41, 1, 1);
      book(_h["jet_pt1_2j_fb"],       46, 1, 1); 
      book(_h["el_eta_0j_pb"],        51, 1, 1);  
      book(_h["el_eta_1j_pb"],        56, 1, 1);   

      book(_h["Wplus_N_incl_pb"],      3, 1, 1);
      book(_h["Wplus_HT_1j_fb"],       8, 1, 1);
      book(_h["Wplus_W_pt_1j_fb"],    13, 1, 1);
      book(_h["Wplus_jet_pt1_1j_fb"], 18, 1, 1); 
      book(_h["Wplus_jet_y1_1j_fb"],  23, 1, 1);
      book(_h["Wplus_HT_2j_fb"],      38, 1, 1);
      book(_h["Wplus_W_pt_2j_fb"],    43, 1, 1);
      book(_h["Wplus_jet_pt1_2j_fb"], 48, 1, 1); 
      book(_h["Wplus_el_eta_0j_pb"],  53, 1, 1);
      book(_h["Wplus_el_eta_1j_pb"],  58, 1, 1);

      book(_h["Wminus_N_incl_pb"],     3, 1, 2);
      book(_h["Wminus_HT_1j_fb"],      8, 1, 2);
      book(_h["Wminus_W_pt_1j_fb"],   13, 1, 2);
      book(_h["Wminus_jet_pt1_1j_fb"],18, 1, 2);  
      book(_h["Wminus_jet_y1_1j_fb"], 23, 1, 2);
      book(_h["Wminus_HT_2j_fb"],     38, 1, 2);
      book(_h["Wminus_W_pt_2j_fb"],   43, 1, 2);
      book(_h["Wminus_jet_pt1_2j_fb"],48, 1, 2); 
      book(_h["Wminus_el_eta_0j_pb"], 53, 1, 2);
      book(_h["Wminus_el_eta_1j_pb"], 58, 1, 2);

      // and ratios
      book(_r["WplusOverWminus_N_incl_pb"],      3, 1, 3);
      book(_r["WplusOverWminus_HT_1j_fb"],       8, 1, 3);
      book(_r["WplusOverWminus_W_pt_1j_fb"],    13, 1, 3);
      book(_r["WplusOverWminus_jet_pt1_1j_fb"], 18, 1, 3);
      book(_r["WplusOverWminus_jet_y1_1j_fb"],  23, 1, 3);
      book(_r["WplusOverWminus_HT_2j_fb"],      38, 1, 3);
      book(_r["WplusOverWminus_W_pt_2j_fb"],    43, 1, 3);
      book(_r["WplusOverWminus_jet_pt1_2j_fb"], 48, 1, 3);
      book(_r["WplusOverWminus_el_eta_0j_pb"],  53, 1, 3);
      book(_r["WplusOverWminus_el_eta_1j_pb"],  58, 1, 3);

    }


    // Perform the per-event analysis
    void analyze(const Event& event) {
    
      // retrieve the dressed electrons
      const Particles& signal_leptons = apply<DressedLeptons>(event, "DressedLeptons").particlesByPt();
      if (signal_leptons.size() != 1 ) vetoEvent;
      const Particle& lepton = signal_leptons[0];

      // calulate MET and mT
      const Particles& invisibles = apply<InvisibleFinalState>(event, "InvFS").particles();
      FourMomentum pMET = sum(invisibles, Kin::mom, FourMomentum()).setZ(0);
      const double MET = pMET.pT() / GeV;
      const double mT = sqrt( 2 * lepton.Et() / GeV * MET * (1 - cos(deltaPhi(lepton,pMET))));
      if ( MET <= 25. ) vetoEvent;
      if ( mT  <= 40. ) vetoEvent;

      // retrieve jets
      Jets jets = apply<FastJets>(event, "Jets").jetsByPt(Cuts::pT > 30*GeV && Cuts::absrap < 4.4);
   
      idiscardIfAnyDeltaRLess(jets, signal_leptons, 0.2);

      // apply event selection on dR
      for (const Jet& j : jets) {
          if (deltaR(j, lepton) < 0.4) vetoEvent;
      }

      // calculate the observables
      const double w_pt = (lepton.momentum() + pMET).pT() / GeV;   
      const size_t njets = jets.size();
      double ST = sum(jets, Kin::pT, 0.0); // scalar pT sum of all selected jets
      const double HT = ST + lepton.pT() / GeV + MET; //missET;

      // fill W histograms		
      _h["N_pb"]->fill(njets);
      for (size_t i = 0; i <= njets; ++i) {
          _h["N_incl_pb"]->fill(i);
      }
      _h["el_eta_0j_pb"]->fill(lepton.abseta());

      if (njets > 0) {
          _h["HT_1j_fb"]->fill(HT);
          _h["W_pt_1j_fb"]->fill(w_pt);
          _h["jet_pt1_1j_fb"]->fill(jets[0].pT()/GeV);
          _h["jet_y1_1j_fb"]->fill(jets[0].absrap());
          _h["el_eta_1j_pb"]->fill(lepton.abseta());
      }
      if (njets > 1) {
          _h["HT_2j_fb"]->fill(HT);
          _h["W_pt_2j_fb"]->fill(w_pt);
          _h["jet_pt1_2j_fb"]->fill(jets[0].pT()/GeV);
          _h["DeltaRj12_2j_fb"]->fill(deltaR(jets[0],jets[1]));
          _h["jet_pt2_2j_fb"]->fill(jets[1].pT()/GeV);
          _h["jet_y2_2j_fb"]->fill(jets[1].absrap());
          _h["jet_mass12_2j_fb"]->fill( (jets[0].mom()+jets[1].mom()).mass()/GeV);
      }		
        // fill W+ histograms
      if (lepton.charge() > 0) {
          for (size_t i = 0; i <= njets; ++i) {
              _h["Wplus_N_incl_pb"]->fill(i);
          }
          _h["Wplus_el_eta_0j_pb"]->fill(lepton.abseta());
          if (njets > 0) {
              _h["Wplus_HT_1j_fb"]->fill(HT);
              _h["Wplus_W_pt_1j_fb"]->fill(w_pt);
              _h["Wplus_jet_pt1_1j_fb"]->fill(jets[0].pT()/GeV);
              _h["Wplus_jet_y1_1j_fb"]->fill(jets[0].absrap());
              _h["Wplus_el_eta_1j_pb"]->fill(lepton.abseta());
              if (njets > 1) {
                  _h["Wplus_HT_2j_fb"]->fill(HT);
                  _h["Wplus_W_pt_2j_fb"]->fill(w_pt);
                  _h["Wplus_jet_pt1_2j_fb"]->fill(jets[0].pT()/GeV);
              }		
          }
      }
      // fill W- histograms
      if (lepton.charge() < 0) {
        for (size_t i = 0; i <= njets; ++i) {
          _h["Wminus_N_incl_pb"]->fill(i);
        }
        _h["Wminus_el_eta_0j_pb"]->fill(lepton.abseta());
        if (njets > 0) {
          _h["Wminus_HT_1j_fb"]->fill(HT);
          _h["Wminus_W_pt_1j_fb"]->fill(w_pt);
          _h["Wminus_jet_pt1_1j_fb"]->fill(jets[0].pT()/GeV);
          _h["Wminus_jet_y1_1j_fb"]->fill(jets[0].absrap());
          _h["Wminus_el_eta_1j_pb"]->fill(lepton.abseta());
          if (njets > 1) {
            _h["Wminus_HT_2j_fb"]->fill(HT);
            _h["Wminus_W_pt_2j_fb"]->fill(w_pt);
            _h["Wminus_jet_pt1_2j_fb"]->fill(jets[0].pT()/GeV);
          }		
        }
      }
    }
    


    void finalize() {
      const double scalefactor_fb = crossSection() / sumOfWeights() / femtobarn;
      const double scalefactor_pb = crossSection() / sumOfWeights() / picobarn;

      for (auto& hit : _h){
        if (hit.first.find("_fb") != string::npos) scale(hit.second, scalefactor_fb);
        else                                       scale(hit.second, scalefactor_pb);
      }
      for (auto& rit : _r) {
        string ratio_label = "WplusOverWminus";
        string num_name   = rit.first;
        string denom_name = rit.first;
        num_name.replace(rit.first.find(ratio_label),ratio_label.length(),"Wplus");
        denom_name.replace(rit.first.find(ratio_label),ratio_label.length(),"Wminus");                   
        divide(_h[num_name], _h[denom_name], rit.second);
      }
    }

    protected:
      // Data members like post-cuts event weight counters go here
      size_t _mode;

    private:
      map<string, Histo1DPtr> _h;
      map<string, Scatter2DPtr> _r;
  };

  // The hook for the plugin system
  RIVET_DECLARE_PLUGIN(ATLAS_2018_I1635273);
}