rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

CMS_2020_I1776758

Ratios of cross sections in the associated production of a Z boson with at least one charm or bottom quark jet are measured in proton-proton collisions at $\sqrt{s}=13$ TeV
Experiment: CMS (LHC)
Inspire ID: 1776758
Status: VALIDATED
Authors:
  • Hannes Jung
References: Beams: p+ p+
Beam energies: (6500.0, 6500.0) GeV
Run details:
  • Z to leptons

Ratios of cross sections, $\sigma(Z+c jets)/\sigma(Z+jets)$, $\sigma(Z+b jets)/\sigma(Z+jets)$ and $\sigma(Z+c jets)/\sigma(Z+b jets)$ in the associated production of a Z boson with at least one charm or bottom quark jet are measured in proton-proton collisions at $\sqrt{s}=13$ TeV. The data sample, collected by the CMS experiment at the CERN LHC, corresponds to an integrated luminosity of 35.9 $fv^{-1}$, with a fiducial volume of $p_T > 30 $ GeV and $|\eta|<2.4$ for the jets. The Z boson candidates come from leptonic decays into electrons or muons with $p_T > 25 $ GeV and $|\eta|<2.4$ and the dilepton mass satisfies $ 71 < M_Z < 111 $ GeV.

Source code: CMS_2020_I1776758.cc
  1// -*- C++ -*-
  2#include "Rivet/Analysis.hh"
  3#include "Rivet/Tools/Cuts.hh"
  4#include "Rivet/Projections/FinalState.hh"
  5#include "Rivet/Projections/VetoedFinalState.hh"
  6#include "Rivet/Projections/LeptonFinder.hh"
  7#include "Rivet/Projections/FastJets.hh"
  8#include "Rivet/Projections/DileptonFinder.hh"
  9
 10namespace Rivet {
 11
 12
 13  /// @brief Ratios of cross-sections in Z + >= 1 c- or b-jet at 13 TeV
 14  class CMS_2020_I1776758 : public Analysis {
 15  public:
 16
 17    /// Constructor
 18    RIVET_DEFAULT_ANALYSIS_CTOR(CMS_2020_I1776758);
 19
 20
 21    // based on SMP_19_004
 22    
 23    /// @name Analysis methods
 24    ///@{
 25
 26    /// Book histograms and initialise projections before the run
 27    void init() {
 28
 29      // Initialise and register projections
 30      DileptonFinder zeeFinder(91.2*GeV, 0.1, Cuts::abseta < 2.4 && Cuts::pT > 25*GeV &&
 31                               Cuts::abspid == PID::ELECTRON, Cuts::massIn(71*GeV, 111*GeV));
 32      declare(zeeFinder, "ZeeFinder");
 33
 34      DileptonFinder zmumuFinder(91.2*GeV, 0.1, Cuts::abseta < 2.4 && Cuts::pT > 25*GeV &&
 35                                 Cuts::abspid == PID::MUON, Cuts::massIn(71*GeV, 111*GeV));
 36      declare(zmumuFinder, "ZmumuFinder");
 37
 38      VisibleFinalState visfs;
 39      VetoedFinalState jetConstits(visfs);
 40      jetConstits.addVetoOnThisFinalState(zeeFinder);
 41      jetConstits.addVetoOnThisFinalState(zmumuFinder);
 42
 43      FastJets akt04Jets(jetConstits, JetAlg::ANTIKT, 0.4);
 44      declare(akt04Jets, "AntiKt04Jets");
 45      
 46      book(_h_jet_pt_combined,"_TMP/jet_pt_combined", refData(1,1,1));
 47      book(_h_jet_pt_cjet_combined,"_TMP/jet_pt_cjet_combined", refData(1,1,1));
 48      book(_h_jet_pt_bjet_combined,"_TMP/jet_pt_bjet_combined", refData(1,1,1));
 49      
 50      book(_h_Z_pt_combined,"_TMP/Z_pt_combined", refData(2, 1, 1));
 51      book(_h_Z_pt_cjet_combined,"_TMP/Z_pt_cjet_combined", refData(2, 1, 1));
 52      book(_h_Z_pt_bjet_combined,"_TMP/Z_pt_bjet_combined", refData(2, 1, 1));
 53      
 54      // book ratio histos      
 55
 56      book(_h_R_jet_pt_cjet_combined, 1, 1, 1);
 57      book(_h_R_jet_pt_bjet_combined, 3, 1, 1);
 58      book(_h_R_jet_pt_cb_combined, 5, 1, 1);
 59
 60      book(_h_R_Z_pt_cjet_combined, 2, 1, 1);
 61      book(_h_R_Z_pt_bjet_combined, 4, 1, 1);
 62      book(_h_R_Z_pt_cb_combined, 6, 1, 1);
 63
 64
 65    }
 66
 67
 68    /// Perform the per-event analysis
 69    void analyze(const Event& event) {
 70
 71      const DileptonFinder& zeeFS = apply<DileptonFinder>(event, "ZeeFinder");
 72      const DileptonFinder& zmumuFS = apply<DileptonFinder>(event, "ZmumuFinder");
 73
 74      const Particles& zees = zeeFS.bosons();
 75      const Particles& zmumus = zmumuFS.bosons();
 76
 77      // We did not find exactly one Z. No good.
 78      if (zees.size() + zmumus.size() != 1) {
 79        MSG_DEBUG("Did not find exactly one good Z candidate");
 80        vetoEvent;
 81      }
 82
 83      //event identification depending on mass window
 84      bool ee_event=false;
 85      bool mm_event=false;
 86            
 87      if (zees.size() == 1) { ee_event = true; }
 88      if (zmumus.size() == 1) { mm_event = true; }
 89      if (!(ee_event || mm_event)) vetoEvent;
 90
 91      // Cluster jets
 92      // NB. Veto has already been applied on leptons and photons used for dressing
 93      const FastJets& fj = apply<FastJets>(event, "AntiKt04Jets");
 94      Jets jets = fj.jetsByPt(Cuts::abseta < 2.4 && Cuts::pT > 30*GeV);
 95      idiscardIfAnyDeltaRLess(jets, ee_event ? zees : zmumus, 0.4);
 96
 97      // We don't care about events with no isolated jets
 98      if (jets.empty()) {
 99        MSG_DEBUG("No jets in event");
100        vetoEvent;
101      }
102
103      Jets jc_final;
104      Jets jb_final;
105            
106      //identification of bjets
107      int n_ctag = 0; 
108      int n_btag = 0; 
109       for (const Jet& j : jets) {
110        if ( j.cTagged() ) { n_ctag = n_ctag + 1; }
111        if ( j.bTagged() ) { n_btag = n_btag + 1; }
112      }
113            
114      for (const Jet& j : jets) {
115        if ( j.cTagged() && n_btag == 0 ) { jc_final.push_back(j); }
116        if ( j.bTagged() ) { jb_final.push_back(j); }
117      }
118      //histogram filling
119
120      if ((ee_event || mm_event) && jets.size() > 0) {
121        
122        FourMomentum j1(jets[0].momentum());
123
124        if ( ee_event ) {
125           _h_Z_pt_combined->fill(zees[0].pt()/GeV); 
126           _h_jet_pt_combined -> fill(j1.pt()/GeV) ;
127           }
128        if ( mm_event ) {
129           _h_Z_pt_combined->fill(zmumus[0].pt()/GeV); 
130           _h_jet_pt_combined -> fill(j1.pt()/GeV) ;
131           }
132           
133     
134        if ( jc_final.size() > 0 ) { 
135          FourMomentum c1(jc_final[0].momentum());
136          _h_jet_pt_cjet_combined -> fill(c1.pt()/GeV);
137          if ( ee_event ) { 
138             _h_Z_pt_cjet_combined->fill(zees[0].pt()/GeV);
139             }
140          if ( mm_event ) { 
141             _h_Z_pt_cjet_combined->fill(zmumus[0].pt()/GeV);
142             }
143          
144        }
145        if ( jb_final.size() > 0 ) { 
146          FourMomentum b1(jb_final[0].momentum());
147          _h_jet_pt_bjet_combined -> fill(b1.pt()/GeV);
148          if ( ee_event ) { 
149             _h_Z_pt_bjet_combined->fill(zees[0].pt()/GeV);
150             }
151          if ( mm_event ) { 
152             _h_Z_pt_bjet_combined->fill(zmumus[0].pt()/GeV);
153             }
154          
155        }
156      }
157
158    }
159
160
161    /// Normalise histograms etc., after the run
162    void finalize() {
163
164      divide( _h_jet_pt_cjet_combined , _h_jet_pt_combined , _h_R_jet_pt_cjet_combined );
165      divide( _h_jet_pt_bjet_combined , _h_jet_pt_combined , _h_R_jet_pt_bjet_combined );
166      divide( _h_jet_pt_cjet_combined , _h_jet_pt_bjet_combined , _h_R_jet_pt_cb_combined );
167      
168      divide ( _h_Z_pt_cjet_combined , _h_Z_pt_combined, _h_R_Z_pt_cjet_combined);
169      divide ( _h_Z_pt_bjet_combined , _h_Z_pt_combined, _h_R_Z_pt_bjet_combined);
170      divide ( _h_Z_pt_cjet_combined , _h_Z_pt_bjet_combined, _h_R_Z_pt_cb_combined);
171
172    }
173
174    ///@}
175
176
177    /// @name Histograms
178     
179     Histo1DPtr _h_jet_pt_combined;
180     Histo1DPtr _h_jet_pt_cjet_combined;
181     Histo1DPtr _h_jet_pt_bjet_combined;
182     Histo1DPtr _h_Z_pt_combined, _h_Z_pt_cjet_combined, _h_Z_pt_bjet_combined;
183
184
185     Estimate1DPtr _h_R_jet_pt_cjet_combined, _h_R_jet_pt_bjet_combined, _h_R_jet_pt_cb_combined;
186
187 
188     Estimate1DPtr _h_R_Z_pt_cjet_combined, _h_R_Z_pt_bjet_combined, _h_R_Z_pt_cb_combined ;
189
190  };
191
192
193  RIVET_DECLARE_PLUGIN(CMS_2020_I1776758);
194
195}