rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

CMS_2017_I1634835

Measurement of associated Z + charm production in proton-proton collisions at sqrts = 8 TeV
Experiment: CMS (LHC)
Inspire ID: 1634835
Status: VALIDATED
Authors:
  • Hannes Jung
  • Juan Pablo Fernandez
References: Beams: p+ p+
Beam energies: (4000.0, 4000.0) GeV
Run details:
  • Z+c production

A study of the associated production of a Z boson and a charm quark jet (Z+c) and a comparison to production with a b quark jet (Z+b), in pp collisions at a centre-of-mass energy of 8 TeV are presented. The analysis uses a data sample cor- responding to an integrated luminosity of 19.7 fb$^{-1}$, collected with the CMS detector at the CERN LHC. The Z boson candidates are identified through their decays into pairs of electrons or muons. Jets originating from heavy flavour quarks are iden- tified using semileptonic decays of c or b flavoured hadrons and hadronic decays of charm hadrons. The measurements are performed in the kinematic region with two leptons with $p_T^l > 20$ GeV, $|\eta^l| < 2.1$, $71 < m_{ll} < 111$ GeV, and heavy flavour jets with $p_T^{jet} > 25$ GeV and $|\eta^{jet}| <2.5$. The Z + c production cross section and the cross section ratio are also measured as a function of the transverse momentum of the Z boson and of the heavy flavour jet.

Source code: CMS_2017_I1634835.cc
  1// -*- C++ -*-
  2#include "Rivet/Analysis.hh"
  3#include "Rivet/Projections/FastJets.hh"
  4#include "Rivet/Projections/FinalState.hh"
  5#include "Rivet/Projections/VetoedFinalState.hh"
  6#include "Rivet/Projections/DileptonFinder.hh"
  7#include "Rivet/Tools/Cuts.hh"
  8
  9namespace Rivet {
 10
 11
 12  /// @brief Z+charm at 8 TeV
 13  class CMS_2017_I1634835 : public Analysis {
 14  public:
 15
 16      /// Constructor
 17      RIVET_DEFAULT_ANALYSIS_CTOR(CMS_2017_I1634835);
 18
 19
 20      /// @name Analysis methods
 21      ///@{
 22
 23      /// Book histograms and initialise projections before the run
 24      void init() {
 25          // Initialise and register projections
 26          DileptonFinder zeeFinder(91.2*GeV, 0.1, Cuts::abseta < 2.1 && Cuts::pT > 20*GeV &&
 27                                   Cuts::abspid == PID::ELECTRON, Cuts::massIn(71*GeV, 111*GeV));
 28          declare(zeeFinder, "ZeeFinder");
 29
 30          DileptonFinder zmumuFinder(91.2*GeV, 0.1, Cuts::abseta < 2.1 && Cuts::pT > 20*GeV &&
 31                                     Cuts::abspid == PID::MUON, Cuts::massIn(71*GeV, 111*GeV));
 32          declare(zmumuFinder, "ZmumuFinder");
 33
 34          VetoedFinalState jetConstits;
 35          jetConstits.addVetoOnThisFinalState(zeeFinder);
 36          jetConstits.addVetoOnThisFinalState(zmumuFinder);
 37
 38          FastJets akt04Jets(jetConstits, JetAlg::ANTIKT, 0.4);
 39          declare(akt04Jets, "AntiKt04Jets");
 40
 41          book(_h_Z_pt_cjet, 4, 1, 1);
 42          book(_h_pt_cjet, 5, 1, 1);
 43          book(_h_Z_pt_bjet, "TMP/Z_pt_bjet", refData(4, 1, 2));
 44          book(_h_pt_bjet, "TMP/pt_bjet", refData(5, 1, 2));
 45
 46          // book ratio histos
 47          book(_h_R_Z_pt_cjet, 4, 1, 2);
 48          book(_h_R_jet_pt_cjet, 5, 1, 2);
 49
 50          book(counter_ee, "TMP/counter_ee");
 51          book(counter_mm, "TMP/counter_mm");
 52      }
 53
 54      /// Perform the per-event analysis
 55      void analyze(const Event& event) {
 56          const DileptonFinder& zeeFS = apply<DileptonFinder>(event, "ZeeFinder");
 57          const DileptonFinder& zmumuFS = apply<DileptonFinder>(event, "ZmumuFinder");
 58
 59          const Particles& zees = zeeFS.bosons();
 60          const Particles& zmumus = zmumuFS.bosons();
 61
 62          // We did not find exactly one Z. No good.
 63          if (zees.size() + zmumus.size() != 1) {
 64              MSG_DEBUG("Did not find exactly one good Z candidate");
 65              vetoEvent;
 66          }
 67
 68          Particles leptons;
 69          Particle zcand;
 70          if (zees.size() == 1) {
 71              leptons = zeeFS.leptons();
 72              zcand = zees[0];
 73              counter_ee->fill();
 74          }
 75          if (zmumus.size() == 1) {
 76              leptons = zmumuFS.leptons();
 77              zcand = zmumus[0];
 78              counter_mm->fill();
 79          }
 80
 81          // Cluster jets
 82          // NB. Veto has already been applied on leptons and photons used for dressing
 83          const FastJets& fj = apply<FastJets>(event, "AntiKt04Jets");
 84          Jets goodjets = fj.jetsByPt(Cuts::abseta < 2.5 && Cuts::pT > 25*GeV);
 85          idiscardIfAnyDeltaRLess(goodjets, leptons, 0.5);
 86
 87          // We don't care about events with no isolated jets
 88          if (goodjets.empty()) {
 89              MSG_DEBUG("No jets in event");
 90              vetoEvent;
 91          }
 92
 93          Jets jc_final;
 94          Jets jb_final;
 95
 96          // identification of bjets
 97          int n_btag = count(goodjets, hasBTag());
 98
 99          for (const Jet& j : goodjets) {
100              if (j.cTagged() && n_btag == 0) {
101                  jc_final.push_back(j);
102              }
103              if (j.bTagged()) {
104                  jb_final.push_back(j);
105              }
106          }
107          // histogram filling
108
109          if (goodjets.size() > 0) {
110              if (jc_final.size() > 0) {
111                  _h_pt_cjet->fill(jc_final[0].pt()/GeV);
112                  _h_Z_pt_cjet->fill(zcand.pt()/GeV);
113              }
114              if (jb_final.size() > 0) {
115                  _h_pt_bjet->fill(jb_final[0].pt()/GeV);
116                  _h_Z_pt_bjet->fill(zcand.pt()/GeV);
117              }
118          }
119      }
120
121      /// Normalise histograms etc., after the run
122      void finalize() {
123          double norm = (sumOfWeights() != 0) ? crossSection() / picobarn / sumOfWeights() : 1.0;
124          // account if we have electrons and muons in the sample
125          if ((counter_ee->val() > 1.) && (counter_mm->val() > 1.)) {
126              norm = norm / 2.;
127          }
128
129          divide(_h_pt_cjet, _h_pt_bjet, _h_R_jet_pt_cjet);
130          divide(_h_Z_pt_cjet, _h_Z_pt_bjet, _h_R_Z_pt_cjet);
131
132          scale(_h_Z_pt_cjet, norm);
133          scale(_h_pt_cjet, norm);
134      }
135
136      ///@}
137
138      /// @name Histograms
139
140      Histo1DPtr _h_pt_cjet;
141      Histo1DPtr _h_pt_bjet;
142      Histo1DPtr _h_Z_pt_cjet, _h_Z_pt_bjet;
143
144      Estimate1DPtr _h_R_jet_pt_cjet;
145      Estimate1DPtr _h_R_Z_pt_cjet;
146
147      CounterPtr counter_ee, counter_mm;
148  };
149
150  RIVET_DECLARE_PLUGIN(CMS_2017_I1634835);
151
152}  // namespace Rivet