rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

CMS_2011_S9215166

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.

Forward energy flow measured by CMS at $\sqrt{s} = 0.9$ and 7 TeV in MB and dijet events.

Source code: CMS_2011_S9215166.cc
  1
  2
  3
  4
  5
  6
  7
  8
  9
 10
 11
 12
 13
 14
 15
 16
 17
 18
 19
 20
 21
 22
 23
 24
 25
 26
 27
 28
 29
 30
 31
 32
 33
 34
 35
 36
 37
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
// -*- C++ -*-
#include "Rivet/Analysis.hh"
#include "Rivet/Projections/FinalState.hh"
#include "Rivet/Projections/ChargedFinalState.hh"
#include "Rivet/Projections/FastJets.hh"
#include "Rivet/Projections/VetoedFinalState.hh"

namespace Rivet {


    class CMS_2011_S9215166 : public Analysis {
        public:

            /// Constructor
            CMS_2011_S9215166() : Analysis("CMS_2011_S9215166") {  }


            void init() {
                const FinalState fs((Cuts::etaIn(-6.0, 6.0)));
                declare(fs, "FS");
                declare(FastJets(fs, FastJets::ANTIKT, 0.5), "Jets");

                VetoedFinalState fsv(fs);
                fsv.vetoNeutrinos();
                fsv.addVetoPairDetail(PID::MUON, 0.0*GeV, 99999.9*GeV);
                declare(fsv, "fsv");

                // For the MB ND selection
                const ChargedFinalState fschrgd((Cuts::etaIn(-6.0,6.0)));
                declare(fschrgd, "fschrgd");
                VetoedFinalState fschrgdv(fschrgd);
                fschrgdv.vetoNeutrinos();
                declare(fschrgdv, "fschrgdv");

                if (fuzzyEquals(sqrtS()/GeV, 900, 1E-3)) {
                    book(_hist_mb      ,1, 1, 1); // energy flow in MB, 0.9 TeV
                    book(_hist_dijet ,2, 1, 1); // energy flow in dijet events, 0.9 TeV
                } else if (fuzzyEquals(sqrtS()/GeV, 7000, 1E-3)) {
                    book(_hist_mb      ,3, 1, 1); // energy flow in MB, 7 TeV
                    book(_hist_dijet ,4, 1, 1); // energy flow in dijet events, 7 TeV
                }

                book(_weightMB, "/tmp/weightMB");
                book(_weightDiJet, "/tmp/weightDijet");
            }


            void analyze(const Event& event) {
                // Skip if the event is empty
                const FinalState& fsv = apply<FinalState>(event, "fsv");
                if (fsv.empty()) vetoEvent;

                // Veto diffractive topologies according to defined hadron level
                double count_chrg_forward = 0;
                double count_chrg_backward = 0;
                const FinalState& fschrgdv = apply<FinalState>(event, "fschrgdv");
                for (const Particle& p : fschrgdv.particles()) {
                    if (3.9 < p.eta() && p.eta() < 4.4) count_chrg_forward++;
                    if (-4.4 < p.eta() && p.eta() < -3.9) count_chrg_backward++;
                }
                if (count_chrg_forward == 0 || count_chrg_backward == 0) vetoEvent;
                /// @todo "Diffractive" veto should really also veto dijet events?


                // MINIMUM BIAS EVENTS
                _weightMB->fill();
                for (const Particle& p: fsv.particles()) {
                    _hist_mb->fill(p.abseta(), p.E()/GeV);
                }


                // DIJET EVENTS
                double PTCUT = -1.0;
                if (fuzzyEquals(sqrtS()/GeV, 900, 1E-3)) PTCUT = 8.0*GeV;
                else if (fuzzyEquals(sqrtS()/GeV, 7000, 1E-3)) PTCUT = 20.0*GeV;
                const FastJets& jetpro = apply<FastJets>(event, "Jets");
                const Jets jets = jetpro.jetsByPt(PTCUT);
                if (jets.size() >= 2) {
                    // eta cut for the central jets
                    if (fabs(jets[0].eta()) < 2.5 && fabs(jets[1].eta()) < 2.5) {
                        // Back to back condition of the jets
                        const double diffphi = deltaPhi(jets[1].phi(), jets[0].phi());
                        if (diffphi-PI < 1.0) {
                            _weightDiJet->fill();
                            for (const Particle& p: fsv.particles()) {
                                _hist_dijet->fill(p.abseta(), p.E()/GeV);
                            }
                        }
                    }
                }

            }


            void finalize() {
                scale(_hist_mb   , 0.5/_weightMB->sumW());
                scale(_hist_dijet, 0.5/_weightDiJet->sumW());
            }


        private:

            Histo1DPtr _hist_mb, _hist_dijet;
            CounterPtr _weightMB, _weightDiJet;

    };


    // Hook for the plugin system
    DECLARE_RIVET_PLUGIN(CMS_2011_S9215166);

}