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