rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

ATLAS_2020_I1788444

Differential cross-sections for Z + b-jets
Experiment: ATLAS (LHC)
Inspire ID: 1788444
Status: VALIDATED
Authors:
  • Federico Sforza
  • Deepak Kar
References: Beams: p+ p+
Beam energies: (6500.0, 6500.0) GeV
Run details:
  • Z+jets at 13 TeV

This paper presents a measurement of the production cross-section of a $Z$ boson in association with $b$-jets, 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 35.6 fb$^{-1}$. Inclusive and differential cross-sections are measured for events containing a $Z$ boson decaying into electrons or muons and produced in association with at least one or at least two $b$-jets with transverse momentum $p_T> 20$ GeV and rapidity $|y|<2.5$. Predictions from several Monte Carlo generators based on leading-order (LO) or next-to-leading-order (NLO) matrix elements interfaced with a parton-shower simulation and testing different flavour schemes for the choice of initial-state partons are compared with measured cross-sections. The 5-flavour number scheme predictions at NLO accuracy agree better with data than 4-flavour number scheme ones. The 4-flavour number scheme predictions underestimate data in events with at least one $b$-jet. N.B.: Data correspond to average of electron and muon channel. Use LMODE=EL,MU to run on individual channels.

Source code: ATLAS_2020_I1788444.cc
  1// -*- C++ -*-
  2#include "Rivet/Analysis.hh"
  3#include "Rivet/Projections/DileptonFinder.hh"
  4#include "Rivet/Projections/FastJets.hh"
  5#include "Rivet/Projections/VetoedFinalState.hh"
  6#include "Rivet/Projections/HeavyHadrons.hh"
  7#include "Rivet/Projections/PromptFinalState.hh"
  8
  9namespace Rivet {
 10
 11
 12  /// @brief Z + b(b) in pp at 13 TeV
 13  class ATLAS_2020_I1788444 : public Analysis {
 14  public:
 15
 16    /// Constructor
 17    RIVET_DEFAULT_ANALYSIS_CTOR(ATLAS_2020_I1788444);
 18
 19    /// Book histograms and initialise projections before the run
 20    void init() {
 21
 22      _mode = 0;
 23      if ( getOption("LMODE") == "EL" ) _mode = 1;
 24      if ( getOption("LMODE") == "MU" ) _mode = 2;
 25
 26      // Define fiducial cuts for the leptons in the DileptonFinder
 27      Cut lepcuts = (Cuts::pT > 27*GeV) & (Cuts::abseta < 2.5);
 28      DileptonFinder zfinderE(91.2*GeV, 0.1, lepcuts && Cuts::abspid == PID::ELECTRON, Cuts::massIn(76*GeV, 106*GeV));
 29      DileptonFinder zfinderM(91.2*GeV, 0.1, lepcuts && Cuts::abspid == PID::MUON, Cuts::massIn(76*GeV, 106*GeV));
 30      declare(zfinderE, "zfinderE");
 31      declare(zfinderM, "zfinderM");
 32      declare(HeavyHadrons(), "HFHadrons");
 33
 34      // Photons
 35      FinalState photons(Cuts::abspid == PID::PHOTON);
 36      // Muons
 37      PromptFinalState bare_mu(Cuts::abspid == PID::MUON, TauDecaysAs::PROMPT);
 38      LeptonFinder all_dressed_mu(bare_mu, photons, 0.1, Cuts::abseta < 2.5);
 39      // Electrons
 40      PromptFinalState bare_el(Cuts::abspid == PID::ELECTRON, TauDecaysAs::PROMPT);
 41      LeptonFinder all_dressed_el(bare_el, photons, 0.1, Cuts::abseta < 2.5);
 42
 43      //Jet forming
 44      VetoedFinalState vfs(FinalState(Cuts::abseta < 4.5));
 45      vfs.addVetoOnThisFinalState(all_dressed_el);
 46      vfs.addVetoOnThisFinalState(all_dressed_mu);
 47
 48      FastJets jets(vfs, JetAlg::ANTIKT, 0.4, JetMuons::ALL, JetInvisibles::DECAY);
 49      declare(jets, "jets");
 50
 51      // Book histos - binning taken from data.yoda
 52      book(_h["i1b_ZpT"],2,1,1);
 53      book(_h["i1b_ZY"],4,1,1);
 54      book(_h["i1b_dPhiZb"],6,1,1);
 55      book(_h["i1b_dRZb"],8,1,1);
 56      book(_h["i1b_dYZb"],7,1,1);
 57      book(_h["i1b_bpT"],3,1,1);
 58      book(_h["i1b_bY"],5,1,1);
 59
 60      book(_h["i2b_ZpT"],13,1,1);
 61      book(_h["i2b_dPhibb"],9,1,1);
 62      book(_h["i2b_dRbb"],11,1,1);
 63      book(_h["i2b_dYbb"],10,1,1);
 64      book(_h["i2b_Mbb"],12,1,1);
 65      book(_h["i2b_pTbb"],14,1,1);
 66      book(_h["i2b_pTOnMbb"],15,1,1);
 67
 68      book(_h["ib_nBJets"],1,1,1);
 69    }
 70
 71
 72    /// Perform the per-event analysis
 73    void analyze(const Event& event) {
 74
 75      const DileptonFinder& zfinderE = apply<DileptonFinder>(event, "zfinderE");
 76      const Particles& els = zfinderE.constituents();
 77      const DileptonFinder& zfinderM = apply<DileptonFinder>(event, "zfinderM");
 78      const Particles& mus = zfinderM.constituents();
 79
 80      // default is to run average of Z->ee and Z->mm
 81      // use LMODE option to pick one channel
 82      if ( (els.size() + mus.size()) != 2 )  vetoEvent;
 83
 84      if (_mode == 0 && !(els.size()==2 || mus.size()==2) )  vetoEvent;
 85      else if ( _mode == 1 && !(els.size() == 2 && mus.empty()) )  vetoEvent;
 86      else if ( _mode == 2 && !(els.empty() && mus.size() == 2) )  vetoEvent;
 87
 88      double Vpt = 0, Vy = 0, Veta = 0, Vphi = 0;
 89
 90      if ( els.size()==2 ) {
 91        Vpt = zfinderE.boson().pt()/GeV;
 92        Vphi = zfinderE.boson().phi();
 93        Vy = zfinderE.boson().rapidity();
 94        Veta = zfinderE.boson().eta();
 95      } else {
 96        Vpt = zfinderM.boson().pt()/GeV;
 97        Vphi = zfinderM.boson().phi();
 98        Vy = zfinderM.boson().rapidity();
 99        Veta = zfinderM.boson().eta();
100      }
101
102      Jets jets = apply<JetFinder>(event, "jets").jetsByPt(Cuts::pT>20*GeV && Cuts::absrap < 2.5);
103      idiscardIfAnyDeltaRLess(jets, els, 0.4);
104      idiscardIfAnyDeltaRLess(jets, mus, 0.4);
105
106      Jets btagged;
107      const Particles allBs = apply<HeavyHadrons>(event, "HFHadrons").bHadrons(Cuts::pT > 5.0*GeV);
108      Particles matchedBs;
109
110      for (const Jet& j : jets) {
111        Jet closest_j;
112        Particle closest_b;
113        double minDR_j_b = 10;
114
115        for (const Particle& bHad : allBs) {
116          bool alreadyMatched = false;
117          for (const Particle& bMatched : matchedBs) {
118            alreadyMatched |= bMatched.isSame(bHad);
119          }
120          if (alreadyMatched)  continue;
121
122          double DR_j_b = deltaR(j, bHad);
123          if ( DR_j_b <= 0.3 && DR_j_b < minDR_j_b) {
124            minDR_j_b = DR_j_b;
125            closest_j = j;
126            closest_b = bHad;
127          }
128        }
129
130        if (minDR_j_b < 0.3) {
131          btagged += closest_j;
132          matchedBs += closest_b;
133        }
134      }
135      //size_t njets = jets.size();
136      size_t ntags = btagged.size();
137      if (ntags < 1) vetoEvent;
138
139      _h["ib_nBJets"]->fill(1); //inclusive 1-b
140
141      double dYVb = fabs(Vy - btagged[0].rap());
142      double dEtaVb = fabs(Veta - btagged[0].eta());
143      double dPhiVb = deltaPhi(Vphi, btagged[0]);
144      double dRVb = sqrt(dEtaVb*dEtaVb + dPhiVb*dPhiVb);
145
146      _h["i1b_ZpT"]   ->fill(Vpt/GeV);
147      _h["i1b_ZY"]    ->fill(fabs(Vy));
148      _h["i1b_dPhiZb"]->fill(dPhiVb);
149      _h["i1b_dRZb"]->fill(dRVb);
150      _h["i1b_dYZb"]->fill(dYVb);
151      _h["i1b_bpT"]->fill(btagged[0].pt()/GeV);
152      _h["i1b_bY"]->fill(btagged[0].absrap());
153
154      if ( ntags>1 ) {
155        _h["ib_nBJets"]->fill(2); //inclusive 2-b
156
157        double dYbb   = fabs(btagged[0].rap() - btagged[1].rap());
158        double dPhibb = deltaPhi(btagged[0], btagged[1]);
159        double dRbb   = deltaR(btagged[0], btagged[1]);
160        double Mbb    = (btagged[0].mom() + btagged[1].mom()).mass()/GeV;
161        double Ptbb   = (btagged[0].mom() + btagged[1].mom()).pt()/GeV;
162
163        _h["i2b_ZpT"]->fill(Vpt);
164        _h["i2b_dPhibb"]->fill(dPhibb);
165        _h["i2b_dRbb"]->fill(dRbb);
166        _h["i2b_dYbb"]->fill(dYbb);
167        _h["i2b_Mbb"]->fill(Mbb);
168        _h["i2b_pTbb"]->fill(Ptbb);
169        _h["i2b_pTOnMbb"]->fill(Ptbb/Mbb);
170
171      }
172    }
173
174
175    void finalize() {
176      // routine accepts both Z->ee and Z->mm
177      // data corresponds to average
178      const double sf = _mode? 1.0 : 0.5;
179      scale(_h, sf * crossSectionPerEvent());
180    }
181
182
183  private:
184
185    size_t _mode;
186
187
188    map<string, Histo1DPtr> _h;
189
190  };
191
192
193  RIVET_DECLARE_PLUGIN(ATLAS_2020_I1788444);
194
195}