rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

ATLAS_2014_I1319490

W + jets
Experiment: ATLAS (LHC)
Inspire ID: 1319490
Status: VALIDATED
Authors:
  • Matthew Mondragon
References: Beams: p+ p+
Beam energies: (3500.0, 3500.0) GeV
Run details:
  • inclusive W production in the electron channel

Measurements of cross sections for the production of a $W$ boson in association with jets in proton-proton collisions at $\sqrt{s} = 7$ TeV with the ATLAS experiment at the Large Hadron Collider. With an integrated luminosity of 4.6 $\text{fb}^{-1}$, this data set allows for an exploration of a large kinematic range, including jet production up to a transverse momentum of 1 TeV and multiplicities up to seven associated jets. The production cross sections for W bosons are measured in both the electron and muon decay channels. Differential cross sections for many observables are also presented including measurements of the jet observables such as the rapidities and the transverse momenta as well as measurements of event observables such as the scalar sums of the transverse momenta of the jets. The default routine assumes both muon and electron decay channel of the $W$ boson are being generate and the average will be constructed. Individual lepton channels can be specified using options LMODE=EL and LMODE:MU respectively.

Source code: ATLAS_2014_I1319490.cc
  1#include "Rivet/Analysis.hh"
  2#include "Rivet/Projections/FinalState.hh"
  3#include "Rivet/Projections/PromptFinalState.hh"
  4#include "Rivet/Projections/FastJets.hh"
  5#include "Rivet/Projections/MissingMomentum.hh"
  6#include "Rivet/Projections/LeptonFinder.hh"
  7#include "Rivet/Projections/VetoedFinalState.hh"
  8
  9namespace Rivet {
 10
 11
 12  /// @brief Electroweak Wjj production at 8 TeV
 13  class ATLAS_2014_I1319490 : public Analysis {
 14  public:
 15
 16    RIVET_DEFAULT_ANALYSIS_CTOR(ATLAS_2014_I1319490);
 17
 18
 19    /// Book histograms and initialise projections before the run
 20    void init() {
 21
 22     // Get options from the new option system
 23      _mode = 0;
 24      if ( getOption("LMODE") == "EL" ) _mode = 1;
 25      if ( getOption("LMODE") == "MU" ) _mode = 2;
 26
 27      Cut cuts;
 28      if (_mode == 2) { // muon channel
 29        cuts = Cuts::pT > 25*GeV && Cuts::abseta < 2.4;
 30      } else if (_mode) { // electron channel
 31        cuts = Cuts::pT > 25*GeV && ( Cuts::abseta < 1.37 || Cuts::absetaIn(1.52, 2.47) );
 32      } else { // combined data extrapolated to common phase space
 33        cuts = Cuts::pT > 25*GeV && Cuts::abseta < 2.5;
 34      }
 35
 36      // Bosons
 37      LeptonFinder ef(0.1, cuts && Cuts::abspid == PID::ELECTRON);
 38      declare(ef, "Elecs");
 39      LeptonFinder mf(0.1, cuts && Cuts::abspid == PID::MUON);
 40      declare(mf, "Muons");
 41      declare(MissingMomentum(), "MET");
 42
 43      // Jets
 44      VetoedFinalState jet_fs;
 45      jet_fs.addVetoOnThisFinalState(ef);
 46      jet_fs.addVetoOnThisFinalState(mf);
 47      FastJets jets(jet_fs, JetAlg::ANTIKT, 0.4, JetMuons::ALL, JetInvisibles::DECAY);
 48      declare(jets, "Jets");
 49
 50      // Book histograms
 51      book(histos["h_N_incl"]            ,1,1,_mode+1);
 52      book(histos["h_N"]                 ,4,1,_mode+1);
 53      book(histos["h_pt_jet1_1jet"]      ,5,1,_mode+1);
 54      book(histos["h_pt_jet1_1jet_excl"] ,6,1,_mode+1);
 55      book(histos["h_pt_jet1_2jet"]      ,7,1,_mode+1);
 56      book(histos["h_pt_jet1_3jet"]      ,8,1,_mode+1);
 57      book(histos["h_pt_jet2_2jet"]      ,9,1,_mode+1);
 58      book(histos["h_pt_jet3_3jet"]      ,10,1,_mode+1);
 59      book(histos["h_pt_jet4_4jet"]      ,11,1,_mode+1);
 60      book(histos["h_pt_jet5_5jet"]      ,12,1,_mode+1);
 61      book(histos["h_y_jet1_1jet"]       ,13,1,_mode+1);
 62      book(histos["h_y_jet2_2jet"]       ,14,1,_mode+1);
 63      book(histos["h_HT_1jet"]           ,15,1,_mode+1);
 64      book(histos["h_HT_1jet_excl"]      ,16,1,_mode+1);
 65      book(histos["h_HT_2jet"]           ,17,1,_mode+1);
 66      book(histos["h_HT_2jet_excl"]      ,18,1,_mode+1);
 67      book(histos["h_HT_3jet"]           ,19,1,_mode+1);
 68      book(histos["h_HT_3jet_excl"]      ,20,1,_mode+1);
 69      book(histos["h_HT_4jet"]           ,21,1,_mode+1);
 70      book(histos["h_HT_5jet"]           ,22,1,_mode+1);
 71      book(histos["h_deltaPhi_jet12"]    ,23,1,_mode+1);
 72      book(histos["h_deltaRap_jet12"]    ,24,1,_mode+1);
 73      book(histos["h_deltaR_jet12"]      ,25,1,_mode+1);
 74      book(histos["h_M_Jet12_2jet"]      ,26,1,_mode+1);
 75      book(histos["h_y_jet3_3jet"]       ,27,1,_mode+1);
 76      book(histos["h_y_jet4_4jet"]       ,28,1,_mode+1);
 77      book(histos["h_y_jet5_5jet"]       ,29,1,_mode+1);
 78      book(histos["h_ST_1jet"]           ,30,1,_mode+1);
 79      book(histos["h_ST_2jet"]           ,31,1,_mode+1);
 80      book(histos["h_ST_2jet_excl"]      ,32,1,_mode+1);
 81      book(histos["h_ST_3jet"]           ,33,1,_mode+1);
 82      book(histos["h_ST_3jet_excl"]      ,34,1,_mode+1);
 83      book(histos["h_ST_4jet"]           ,35,1,_mode+1);
 84      book(histos["h_ST_5jet"]           ,36,1,_mode+1);
 85    }
 86
 87
 88    /// A convenience method to do all the plot-filling work
 89    void fillPlots(const Particle& lepton, const double& missET, Jets& all_jets) {
 90      // do jet-lepton overlap removal
 91      Jets jets;
 92      double ST = 0.0; // scalar pT sum of all selected jets
 93      for (const Jet &j : all_jets) {
 94        if (deltaR(j, lepton) > 0.5) {
 95          jets += j;
 96          ST += j.pT() / GeV;
 97        }
 98      }
 99
100      const size_t njets = jets.size();
101
102      const double HT = ST + lepton.pT() / GeV + missET;
103
104      histos["h_N"]->fill(njets + 0.5);
105      for (size_t i = 0; i <= njets; ++i) {
106        histos["h_N_incl"]->fill(i + 0.5);
107      }
108
109      if (njets) {
110        const double pT1  = jets[0].pT() / GeV;
111        const double rap1 = jets[0].absrap();
112        histos["h_pt_jet1_1jet" ]->fill(pT1);
113        histos["h_y_jet1_1jet"]->fill(rap1);
114        histos["h_HT_1jet"]->fill(HT);
115        histos["h_ST_1jet"]->fill(ST);
116        if (njets == 1) {
117          histos["h_pt_jet1_1jet_excl"]->fill(pT1);
118          histos["h_HT_1jet_excl"]->fill(HT);
119        } else {
120          const double pT2  = jets[1].pT() / GeV;
121          const double rap2 = jets[1].absrap();
122          const double dR   = deltaR(jets[0], jets[1]);
123          const double dRap = deltaRap(jets[0], jets[1]);
124          const double dPhi = deltaPhi(jets[0], jets[1]);
125          const double mjj  = (jets[0].momentum() + jets[1].momentum()).mass() / GeV;
126          histos["h_pt_jet1_2jet"]->fill(pT1);
127          histos["h_pt_jet2_2jet"]->fill(pT2);
128          histos["h_y_jet2_2jet"]->fill(rap2);
129          histos["h_M_Jet12_2jet"]->fill(mjj);
130          histos["h_HT_2jet"]->fill(HT);
131          histos["h_ST_2jet"]->fill(ST);
132          histos["h_deltaPhi_jet12"]->fill(dPhi);
133          histos["h_deltaRap_jet12"]->fill(dRap);
134          histos["h_deltaR_jet12"]->fill(dR);
135          if (njets == 2) {
136            histos["h_ST_2jet_excl"]->fill(ST);
137            histos["h_HT_2jet_excl"]->fill(HT);
138          } else {
139            const double pT3  = jets[2].pT() / GeV;
140            const double rap3 = jets[2].absrap();
141            histos["h_pt_jet1_3jet"]->fill(pT1);
142            histos["h_pt_jet3_3jet"]->fill(pT3);
143            histos["h_y_jet3_3jet"]->fill(rap3);
144            histos["h_HT_3jet"]->fill(HT);
145            histos["h_ST_3jet"]->fill(ST);
146            if(njets == 3) {
147              histos["h_ST_3jet_excl"]->fill(ST);
148              histos["h_HT_3jet_excl"]->fill(HT);
149            } else {
150              const double pT4  = jets[3].pT() / GeV;
151              const double rap4 = jets[3].absrap();
152              histos["h_pt_jet4_4jet"]->fill(pT4);
153              histos["h_y_jet4_4jet"]->fill(rap4);
154              histos["h_HT_4jet"]->fill(HT);
155              histos["h_ST_4jet"]->fill(ST);
156              if (njets > 4) {
157                const double pT5  = jets[4].pT() / GeV;
158                const double rap5 = jets[4].absrap();
159                histos["h_pt_jet5_5jet"]->fill(pT5);
160                histos["h_y_jet5_5jet"]->fill(rap5);
161                histos["h_HT_5jet"]->fill(HT);
162                histos["h_ST_5jet"]->fill(ST);
163              }
164            }
165          }
166        }
167      }
168    }
169
170
171    /// Perform the per-event analysis
172    void analyze(const Event& event) {
173      // MET cut
174      const P4& pmiss = apply<MissingMom>(event, "MET").missingMom();
175      if (pmiss.pT() < 25*GeV) vetoEvent;
176
177      // Retrieve boson candidate
178      const Particles& es = apply<LeptonFinder>(event, "Elecs").particles();
179      const Particles es_mtfilt = select(es, [&](const Particle& e){ return mT(e, pmiss) > 40*GeV; });
180      const int iefound = closestMatchIndex(es_mtfilt, pmiss, Kin::mass, 80.4*GeV);
181      const Particles& mus = apply<LeptonFinder>(event, "Muons").particles();
182      const Particles mus_mtfilt = select(mus, [&](const Particle& m){ return mT(m, pmiss) > 40*GeV; });
183      const int imfound = closestMatchIndex(mus_mtfilt, pmiss, Kin::mass, 80.4*GeV);
184
185      // Restrict allowed W-candidate combinations
186      size_t nWmu = (imfound >= 0);
187      size_t nWel = (iefound >= 0);
188      if (_mode == 0 && !((nWmu == 1 && !nWel) || (!nWmu && nWel == 1)))  vetoEvent; // one W->munu OR W->elnu candidate, otherwise veto
189      if (_mode == 1 && !(!nWmu && nWel == 1))  vetoEvent; // one W->elnu candidate, otherwise veto
190      if (_mode == 2 && !(nWmu == 1 && !nWel))  vetoEvent; // one W->munu candidate, otherwise veto
191
192      // Retrieve jets
193      const JetFinder& jetfs = apply<JetFinder>(event, "Jets");
194      Jets all_jets = jetfs.jetsByPt(Cuts::pT > 30*GeV && Cuts::absrap < 4.4);
195
196      // Fill histograms
197      const Particle& lepton = nWmu ? mus_mtfilt[imfound] : es_mtfilt[iefound];
198      fillPlots(lepton, pmiss.pT()/GeV, all_jets);
199    }
200
201
202    /// Finalize data objects after the run
203    void finalize() {
204      const double sf = _mode? 1.0 : 0.5;
205      const double scalefactor = sf * crossSection()/picobarn / sumOfWeights();
206      scale(histos, scalefactor);
207    }
208
209
210  protected:
211
212    size_t _mode;
213
214
215  private:
216
217    map<string, Histo1DPtr> histos;
218
219  };
220
221
222  RIVET_DECLARE_PLUGIN(ATLAS_2014_I1319490);
223
224}