rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

CMS_2017_I1499471

Measurements of the associated production of a Z boson and b jets in pp collisions at $\sqrt{s} = 8$ TeV
Experiment: CMS (LHC)
Inspire ID: 1499471
Status: VALIDATED
Authors:
  • Fabio Cossutti
References: Beams: p+ p+
Beam energies: (4000.0, 4000.0) GeV
Run details:
  • Run MC generators with Z decaying leptonically into both electrons and muons at 8 TeV CoM energy. For jet-flavour inclusive samples, order of 20 million unweighted events can give a reasonable global comparison, but precision in the high jet multiplicity region/high jet pt may require substantially larger samples or statistical enhancement of high jet multiplicities.

Measurements of the associated production of a Z boson with at least one jet originating from a b quark in proton-proton collisions at $\sqrt{s} = 8$ TeV are presented. Differential cross sections are measured with data collected by the CMS experiment corresponding to an integrated luminosity of 19.8 $\text{fb}^{-1}$. Z bosons are reconstructed through their decays to electrons and muons. Cross sections are measured as a function of observables characterizing the kinematics of the b jet and the Z boson. Ratios of differential cross sections for the associated production with at least one b jet to the associated production with any jet are also presented. The production of a Z boson with at least two b jets is investigated, and differential cross sections are measured for the dijet system. Results are compared to theoretical predictions, testing two different flavour schemes for the choice of initial-state partons. Note of the authors: the Z boson selection and b-tagging algorithms used in the code are based on standard Rivet code and do not exactly coincide with those used in the paper. They anyway reproduce the paper results at the percent level or better for most of the bins.

Source code: CMS_2017_I1499471.cc
  1#include "Rivet/Analysis.hh"
  2#include "Rivet/Tools/Cuts.hh"
  3#include "Rivet/Projections/FinalState.hh"
  4#include "Rivet/Projections/VetoedFinalState.hh"
  5#include "Rivet/Projections/DressedLeptons.hh"
  6#include "Rivet/Projections/FastJets.hh"
  7#include "Rivet/Projections/ZFinder.hh"
  8
  9//#define DebugLog 
 10
 11namespace Rivet {
 12
 13  class CMS_2017_I1499471 : public Analysis {
 14  public:
 15    
 16    /// Constructor
 17    RIVET_DEFAULT_ANALYSIS_CTOR(CMS_2017_I1499471);
 18    
 19    /// Book histograms and initialise projections before the run
 20    void init() {
 21
 22#ifdef DebugLog
 23      // set optionally the verbosity for the internal Rivet message system
 24      getLog().setLevel(0);
 25#endif      
 26
 27      FinalState fs; ///< @todo No cuts?
 28      VisibleFinalState visfs(fs);
 29
 30      ZFinder zeeFinder(fs, Cuts::abseta < 2.4 && Cuts::pT > 20*GeV, PID::ELECTRON, 71.0*GeV, 111.0*GeV, 0.1 );
 31      declare(zeeFinder, "ZeeFinder");
 32
 33      ZFinder zmumuFinder(fs, Cuts::abseta < 2.4 && Cuts::pT > 20*GeV, PID::MUON, 71.0*GeV, 111.0*GeV, 0.1 );
 34      declare(zmumuFinder, "ZmumuFinder");
 35
 36      VetoedFinalState jetConstits(visfs);
 37      jetConstits.addVetoOnThisFinalState(zeeFinder);
 38      jetConstits.addVetoOnThisFinalState(zmumuFinder);
 39
 40      FastJets akt05Jets(jetConstits, FastJets::ANTIKT, 0.5);
 41      declare(akt05Jets, "AntiKt05Jets");
 42      
 43      //Histograms booking
 44      
 45      book(_h_first_bjet_pt_b ,1,1,1);
 46      book(_h_first_bjet_abseta_b ,3,1,1);
 47      book(_h_Z_pt_b ,5,1,1);
 48      book(_h_HT_b ,7,1,1);
 49      book(_h_Dphi_Zb_b ,9,1,1);
 50      
 51      book(_h_first_jet_pt_ratio ,2,1,1);
 52      book(_h_first_jet_abseta_ratio ,4,1,1);
 53      book(_h_Z_pt_ratio ,6,1,1);
 54      book(_h_HT_ratio ,8,1,1);
 55      book(_h_Dphi_Zj_ratio ,10,1,1);
 56      
 57      book(_h_first_jet_pt, "first_jet_pt", refData(1,1,1) ); // (*_h_first_bjet_pt_b);
 58      book(_h_first_jet_abseta, "first_jet_abseta", refData(3,1,1) ); // (*_h_first_bjet_abseta_b);
 59      book(_h_Z_pt, "Z_pt", refData(5,1,1) ); // (*_h_Z_pt_b);
 60      book(_h_HT, "HT", refData(7,1,1) ); // (*_h_HT_b);
 61      book(_h_Dphi_Zj, "Dphi_Zj", refData(9,1,1) ); // (*_h_Dphi_Zb_b);
 62
 63      book(_h_first_bjet_pt_bb ,11,1,1);
 64      book(_h_second_bjet_pt_bb ,12,1,1);
 65      book(_h_Z_pt_bb ,13,1,1);
 66      book(_h_bb_mass_bb ,14,1,1);
 67      book(_h_Zbb_mass_bb ,15,1,1);
 68      book(_h_Dphi_bb ,16,1,1);
 69      book(_h_DR_bb ,17,1,1);
 70      book(_h_DR_Zbmin_bb ,18,1,1);
 71      book(_h_A_DR_Zb_bb ,19,1,1);
 72
 73      book(_h_bjet_multiplicity ,20,1,1);
 74
 75    }
 76                
 77        
 78    /// Perform the per-event analysis
 79    void analyze(const Event& event) {
 80            
 81      const ZFinder& zeeFS = applyProjection<ZFinder>(event, "ZeeFinder");
 82      const ZFinder& zmumuFS = applyProjection<ZFinder>(event, "ZmumuFinder");
 83
 84      const Particles& zees = zeeFS.bosons();
 85      const Particles& zmumus = zmumuFS.bosons();
 86
 87      // We did not find exactly one Z. No good.
 88      if (zees.size() + zmumus.size() != 1) {
 89        MSG_DEBUG("Did not find exactly one good Z candidate");
 90        vetoEvent;
 91      }
 92
 93      //event identification depending on mass window
 94      bool ee_event=false;
 95      bool mm_event=false;
 96            
 97      if (zees.size() == 1) { ee_event = true; }
 98      if (zmumus.size() == 1) { mm_event = true; }
 99      const Particles& theLeptons = zees.size() ? zeeFS.constituents() : zmumuFS.constituents();
100
101      // Cluster jets
102      // NB. Veto has already been applied on leptons and photons used for dressing
103      const FastJets& fj = applyProjection<FastJets>(event, "AntiKt05Jets");
104      const Jets& jets = fj.jetsByPt(Cuts::abseta < 2.4 && Cuts::pT > 30*GeV);
105
106      // Perform lepton-jet overlap and HT calculation
107      double Ht = 0;
108      Jets goodjets;
109      for (const Jet& j : jets) {
110        // Decide if this jet is "good", i.e. isolated from the leptons
111        /// @todo Nice use-case for any() and a C++11 lambda
112        bool overlap = false;
113        for (const Particle& l : theLeptons) {
114          if (Rivet::deltaR(j, l) < 0.5) {
115            overlap = true;
116            break;
117          }
118        }
119
120        // Fill HT and good-jets collection
121        if (overlap) continue;
122        goodjets.push_back(j);
123        Ht += j.pT();
124      }
125
126      // We don't care about events with no isolated jets
127      if (goodjets.empty()) {
128        MSG_DEBUG("No jets in event");
129        vetoEvent;
130      }
131
132      Jets jb_final;
133            
134      //identification of bjets
135            
136      for (const Jet& j : goodjets) {
137        if ( j.bTagged() ) { jb_final.push_back(j); }
138      }
139            
140      //Event weight
141      const double w = 0.5;
142            
143      //histogram filling
144
145      if ((ee_event || mm_event) && goodjets.size() > 0) {
146        
147        FourMomentum j1(goodjets[0].momentum());
148
149        _h_first_jet_pt->fill(j1.pt(),w);
150        _h_first_jet_abseta->fill(fabs(j1.eta()),w);
151        if ( ee_event ) _h_Z_pt->fill(zees[0].pt(),w);
152        if ( mm_event ) _h_Z_pt->fill(zmumus[0].pt(),w);
153        _h_HT->fill(Ht,w);
154        if ( ee_event ) _h_Dphi_Zj->fill(deltaPhi(zees[0], j1),w);
155        if ( mm_event ) _h_Dphi_Zj->fill(deltaPhi(zmumus[0], j1),w);
156        
157        if ( jb_final.size() > 0 ) { 
158
159          FourMomentum b1(jb_final[0].momentum());
160
161          _h_bjet_multiplicity->fill(1.,w);
162
163          _h_first_bjet_pt_b->fill(b1.pt(),w);
164          _h_first_bjet_abseta_b->fill(fabs(b1.eta()),w);
165          if ( ee_event ) _h_Z_pt_b->fill(zees[0].pt(),w);
166          if ( mm_event ) _h_Z_pt_b->fill(zmumus[0].pt(),w);
167          _h_HT_b->fill(Ht,w);
168          if ( ee_event ) _h_Dphi_Zb_b->fill(deltaPhi(zees[0], b1.phi()),w);
169          if ( mm_event ) _h_Dphi_Zb_b->fill(deltaPhi(zmumus[0], b1.phi()),w);
170
171          if ( jb_final.size() > 1 ) {
172
173            FourMomentum b2(jb_final[1].momentum());
174
175            _h_bjet_multiplicity->fill(2.,w);
176
177            _h_first_bjet_pt_bb->fill(b1.pt(),w);
178            _h_second_bjet_pt_bb->fill(b2.pt(),w);
179            if ( ee_event ) _h_Z_pt_bb->fill(zees[0].pt(),w);
180            if ( mm_event ) _h_Z_pt_bb->fill(zmumus[0].pt(),w);
181
182            FourMomentum bb = add(b1,b2);
183            FourMomentum Zbb;
184            if (ee_event) Zbb = add(zees[0],bb);
185            if (mm_event) Zbb = add(zmumus[0],bb);
186
187            _h_bb_mass_bb->fill(bb.mass(),w);
188            _h_Zbb_mass_bb->fill(Zbb.mass(),w);
189
190            _h_Dphi_bb->fill(deltaPhi(b1,b2),w);
191            if (deltaR(b1,b2)>0.5) {
192              _h_DR_bb->fill(deltaR(b1,b2),w);
193            }
194
195            double DR_Z_b1(0.), DR_Z_b2(0.);
196            if ( ee_event ) {
197              DR_Z_b1 = deltaR(zees[0],b1);
198              DR_Z_b2 = deltaR(zees[0],b2);
199            }
200            if ( mm_event ) {
201              DR_Z_b1 = deltaR(zmumus[0],b1);
202              DR_Z_b2 = deltaR(zmumus[0],b2);
203            }
204
205            double DR_Zb_min = DR_Z_b1;
206            double DR_Zb_max = DR_Z_b2;
207            if ( DR_Zb_min > DR_Zb_max ) {
208              DR_Zb_min = DR_Z_b2;
209              DR_Zb_max = DR_Z_b1;
210            }
211            double A_Zbb = (DR_Zb_max - DR_Zb_min)/(DR_Zb_max + DR_Zb_min);
212
213            _h_DR_Zbmin_bb->fill(DR_Zb_min,w);
214            _h_A_DR_Zb_bb->fill(A_Zbb,w);
215
216          }
217          
218        }
219                                           
220      }
221
222    }
223   
224        
225    /// Normalise histograms etc., after the run
226    void finalize() {
227
228      const double norm = (sumOfWeights() != 0) ? crossSection()/picobarn/sumOfWeights() : 1.0;
229
230      MSG_INFO("Cross section = " << std::setfill(' ') << std::setw(14) << std::fixed << std::setprecision(3) << crossSection() << " pb");
231      MSG_INFO("# Events      = " << std::setfill(' ') << std::setw(14) << std::fixed << std::setprecision(3) << numEvents() );
232      MSG_INFO("SumW          = " << std::setfill(' ') << std::setw(14) << std::fixed << std::setprecision(3) << sumOfWeights());
233      MSG_INFO("Norm factor   = " << std::setfill(' ') << std::setw(14) << std::fixed << std::setprecision(6) << norm);
234
235      scale( _h_first_bjet_pt_b, 100. );
236      scale( _h_first_bjet_abseta_b, 100. );
237      scale( _h_Z_pt_b, 100. );
238      scale( _h_HT_b, 100. );
239      scale( _h_Dphi_Zb_b, 100. );
240
241      divide( _h_first_bjet_pt_b , _h_first_jet_pt , _h_first_jet_pt_ratio );
242      divide( _h_first_bjet_abseta_b , _h_first_jet_abseta , _h_first_jet_abseta_ratio );
243      divide( _h_Z_pt_b , _h_Z_pt , _h_Z_pt_ratio );
244      divide( _h_HT_b , _h_HT , _h_HT_ratio );
245      divide( _h_Dphi_Zb_b , _h_Dphi_Zj , _h_Dphi_Zj_ratio );
246
247      scale( _h_first_bjet_pt_b, norm/100. );
248      scale( _h_first_bjet_abseta_b, norm/100. );
249      scale( _h_Z_pt_b, norm/100. );
250      scale( _h_HT_b, norm/100. );
251      scale( _h_Dphi_Zb_b, norm/100. );
252
253      scale( _h_first_bjet_pt_bb, norm);
254      scale( _h_second_bjet_pt_bb, norm);
255      scale( _h_Z_pt_bb, norm);
256      scale( _h_bb_mass_bb, norm);
257      scale( _h_Zbb_mass_bb, norm);
258      scale( _h_Dphi_bb, norm);
259      scale( _h_DR_bb, norm);
260      scale( _h_DR_Zbmin_bb, norm);
261      scale( _h_A_DR_Zb_bb, norm);
262
263      scale( _h_bjet_multiplicity, norm );
264
265    }
266
267
268  private:
269
270    /// @name Histograms
271    
272    Histo1DPtr     _h_first_jet_pt, _h_first_bjet_pt_b;
273    Histo1DPtr     _h_first_jet_abseta, _h_first_bjet_abseta_b;
274    Histo1DPtr     _h_Z_pt, _h_Z_pt_b;
275    Histo1DPtr     _h_HT, _h_HT_b;
276    Histo1DPtr     _h_Dphi_Zj, _h_Dphi_Zb_b;
277
278    Scatter2DPtr     _h_first_jet_pt_ratio;
279    Scatter2DPtr     _h_first_jet_abseta_ratio;
280    Scatter2DPtr     _h_Z_pt_ratio;
281    Scatter2DPtr     _h_HT_ratio;
282    Scatter2DPtr     _h_Dphi_Zj_ratio;
283    
284    Histo1DPtr     _h_first_bjet_pt_bb, _h_second_bjet_pt_bb;
285    Histo1DPtr     _h_Z_pt_bb;
286    Histo1DPtr     _h_bb_mass_bb, _h_Zbb_mass_bb;
287    Histo1DPtr     _h_Dphi_bb, _h_DR_bb, _h_DR_Zbmin_bb, _h_A_DR_Zb_bb;
288    
289    Histo1DPtr     _h_bjet_multiplicity;
290
291  };
292  
293  
294  // The hook for the plugin system
295  RIVET_DECLARE_PLUGIN(CMS_2017_I1499471);
296  
297}