rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

CMS_2011_I930319

Forward energy flow in MB and dijet events at 0.9 and 7 TeV
Experiment: CMS (LHC)
Inspire ID: 930319
Status: VALIDATED
Authors:
  • S. Dooling
  • A. Knutsson
References: Beams: p+ p+
Beam energies: (450.0, 450.0); (3500.0, 3500.0) GeV
Run details:
  • $pp$ MB and QCD interactions at 0.9 and 7 TeV. No pT-cuts. Beam energy must be specified as analysis option "ENERGY" when rivet-merge'ing samples.

Forward energy flow measured by CMS at $\sqrt{s} = 0.9$ and 7 TeV in MB and dijet events. Beam energy must be specified (in GeV) as analysis option "ENERGY" when rivet-merging samples.'

Source code: CMS_2011_I930319.cc
  1// -*- C++ -*-
  2#include "Rivet/Analysis.hh"
  3#include "Rivet/Projections/FinalState.hh"
  4#include "Rivet/Projections/ChargedFinalState.hh"
  5#include "Rivet/Projections/FastJets.hh"
  6#include "Rivet/Projections/VetoedFinalState.hh"
  7
  8namespace Rivet {
  9
 10
 11  /// Forward energy flow in MB and dijet events at 0.9 and 7 TeV
 12  class CMS_2011_I930319 : public Analysis {
 13  public:
 14
 15    /// Constructor
 16    RIVET_DEFAULT_ANALYSIS_CTOR(CMS_2011_I930319);
 17
 18
 19    void init() {
 20      const FinalState fs((Cuts::etaIn(-6.0, 6.0)));
 21      declare(fs, "FS");
 22      declare(FastJets(fs, JetAlg::ANTIKT, 0.5), "Jets");
 23
 24      VetoedFinalState fsv(fs);
 25      fsv.vetoNeutrinos();
 26      fsv.addVetoPairDetail(PID::MUON, 0.0*GeV, 99999.9*GeV);
 27      declare(fsv, "fsv");
 28
 29      // For the MB ND selection
 30      const ChargedFinalState fschrgd(Cuts::abseta < 6.0);
 31      declare(fschrgd, "fschrgd");
 32      VetoedFinalState fschrgdv(fschrgd);
 33      fschrgdv.vetoNeutrinos();
 34      declare(fschrgdv, "fschrgdv");
 35
 36      if (isCompatibleWithSqrtS(900*GeV)) {
 37        book(_hist_mb      ,1, 1, 1); // energy flow in MB, 0.9 TeV
 38        book(_hist_dijet ,2, 1, 1); // energy flow in dijet events, 0.9 TeV
 39      } else if (isCompatibleWithSqrtS(7000*GeV)) {
 40        book(_hist_mb      ,3, 1, 1); // energy flow in MB, 7 TeV
 41        book(_hist_dijet ,4, 1, 1); // energy flow in dijet events, 7 TeV
 42      }
 43
 44      book(_weightMB, "/tmp/weightMB");
 45      book(_weightDiJet, "/tmp/weightDijet");
 46    }
 47
 48
 49    void analyze(const Event& event) {
 50      // Skip if the event is empty
 51      const FinalState& fsv = apply<FinalState>(event, "fsv");
 52      if (fsv.empty()) vetoEvent;
 53
 54      // Veto diffractive topologies according to defined hadron level
 55      double count_chrg_forward = 0;
 56      double count_chrg_backward = 0;
 57      const FinalState& fschrgdv = apply<FinalState>(event, "fschrgdv");
 58      for (const Particle& p : fschrgdv.particles()) {
 59        if (3.9 < p.eta() && p.eta() < 4.4) count_chrg_forward++;
 60        if (-4.4 < p.eta() && p.eta() < -3.9) count_chrg_backward++;
 61      }
 62      if (count_chrg_forward == 0 || count_chrg_backward == 0) vetoEvent;
 63      /// @todo "Diffractive" veto should really also veto dijet events?
 64
 65
 66      // MINIMUM BIAS EVENTS
 67      _weightMB->fill();
 68      for (const Particle& p: fsv.particles()) {
 69        _hist_mb->fill(p.abseta(), p.E()/GeV);
 70      }
 71
 72
 73      // DIJET EVENTS
 74      double PTCUT = -1.0;
 75      if (isCompatibleWithSqrtS(900*GeV)) PTCUT = 8.0*GeV;
 76      else if (isCompatibleWithSqrtS(7000*GeV)) PTCUT = 20.0*GeV;
 77      const FastJets& jetpro = apply<FastJets>(event, "Jets");
 78      const Jets jets = jetpro.jetsByPt(Cuts::pT > PTCUT);
 79      if (jets.size() >= 2) {
 80        // eta cut for the central jets
 81        if (fabs(jets[0].eta()) < 2.5 && fabs(jets[1].eta()) < 2.5) {
 82          // Back to back condition of the jets
 83          const double diffphi = deltaPhi(jets[1].phi(), jets[0].phi());
 84          if (diffphi-PI < 1.0) {
 85            _weightDiJet->fill();
 86            for (const Particle& p: fsv.particles()) {
 87              _hist_dijet->fill(p.abseta(), p.E()/GeV);
 88            }
 89          }
 90        }
 91      }
 92
 93    }
 94
 95
 96    void finalize() {
 97      scale(_hist_mb   , 0.5/_weightMB->sumW());
 98      scale(_hist_dijet, 0.5/_weightDiJet->sumW());
 99    }
100
101
102  private:
103
104    /// @{
105    Histo1DPtr _hist_mb, _hist_dijet;
106    CounterPtr _weightMB, _weightDiJet;
107    /// @}
108
109  };
110
111
112
113  RIVET_DECLARE_ALIASED_PLUGIN(CMS_2011_I930319, CMS_2011_S9215166);
114
115}