rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

CMS_2012_I1111014

Shape, Transverse Size, and Charged Hadron Multiplicity of Jets in pp Collisions at $\sqrt{s}=7$ TeV
Experiment: CMS (LHC)
Inspire ID: 1111014
Status: VALIDATED
Authors:
  • Andreas Hinzmann
References:
  • arXiv: 1204.3170
  • inspirehep 1111014
  • http://cms-results.web.cern.ch/cms-results/public-results/publications/QCD-10-029/
Beams: p+ p+
Beam energies: (3500.0, 3500.0) GeV
Run details:
  • pp QCD interactions at $\sqrt{s} = 7$ TeV. Data collected by CMS during the year 2010.

Measurements of jet characteristics from inclusive jet production in proton-proton collisions at a centre-of-mass energy of 7 TeV are presented. The data sample was collected with the CMS detector at the LHC during 2010 and corresponds to an integrated luminosity of 36 inverse picobarns. The mean charged hadron multiplicity, the differential and integral jet shape distributions, and two independent moments of the shape distributions are measured as functions of the jet transverse momentum for jets reconstructed with the anti-kT algorithm. The measured observables are corrected to the particle level and compared with predictions from various QCD Monte Carlo generators.

Source code: CMS_2012_I1111014.cc
  1// -*- C++ -*-
  2#include "Rivet/Analysis.hh"
  3#include "Rivet/Projections/FastJets.hh"
  4#include "Rivet/Projections/VetoedFinalState.hh"
  5#include "Rivet/Projections/VisibleFinalState.hh"
  6#include "Rivet/Projections/JetShape.hh"
  7
  8namespace Rivet {
  9
 10
 11  /// @brief CMS jet shape analysis
 12  /// @author Andreas Hinzmann
 13  class CMS_2012_I1111014 : public Analysis {
 14  public:
 15
 16    /// Constructor
 17    CMS_2012_I1111014()
 18      : Analysis("CMS_2012_I1111014")
 19    {    }
 20
 21
 22    /// @name Analysis methods
 23    /// @{
 24
 25    void init() {
 26      // Set up projections
 27      const FinalState fs((Cuts::etaIn(-5.0, 5.0)));
 28      declare(fs, "FS");
 29      FastJets fj5(fs, JetAlg::ANTIKT, 0.5);
 30      declare(fj5, "Jets5");
 31      FastJets fj7(fs, JetAlg::ANTIKT, 0.7);
 32      declare(fj7, "Jets7");
 33
 34      // Specify pT bins
 35      _ptedges = {{ 20.,25.,30.,40.,50.,60.,70.,80.,90.,100.,110.,125.,140.,160.,180.,200.,225.,250.,300.,400.,500.,600.,1000. }};
 36      _yedges  = {{ 0.0, 0.5, 1.0, 1.5, 2.0, 2.5, 3.0 }};
 37
 38      // Register a jet shape projection and histogram for each pT bin
 39      unsigned int histo_counter=1;
 40      for (size_t j = 0; j < 6; ++j) {
 41        for (size_t i = 0; i < 22; ++i) {
 42          if (i > 20 && j == 3) continue;
 43          if (i > 18 && j >= 4) continue;
 44
 45          // Set up projections for each (pT,y) bin
 46          _jsnames_pT[i][j] = "JetShape" + to_str(i) + "_" + to_str(j);
 47          const JetShape jsp(fj7, 0.0, 0.7, 7, _ptedges[i], _ptedges[i+1], _yedges[j], _yedges[j+1], RAPIDITY);
 48          declare(jsp, _jsnames_pT[i][j]);
 49
 50          // Book profile histograms for each (pT,y) bin
 51          book(_profhistRho_pT[i][j], histo_counter, 1, 1);
 52          histo_counter+=1;
 53        }
 54      }
 55      book(_profhistNch[0], 126, 1, 1);
 56      book(_profhistNch[1], 126, 1, 2);
 57      book(_profhistDr[0], 127, 1, 1);
 58      book(_profhistDr[1], 127, 1, 2);
 59      book(_profhistDeta, "TMP/Deta", refData(127,1,1));
 60      book(_profhistDphi, "TMP/Dphi", refData(127,1,1));
 61      book(_profhistAsym, "d128-x01-y01");
 62
 63    }
 64
 65
 66
 67    /// Do the analysis
 68    void analyze(const Event& evt) {
 69
 70      // Get jets and require at least one to pass pT and y cuts
 71      Jets jets7 = apply<FastJets>(evt, "Jets7")
 72        .jetsByPt(Cuts::ptIn(_ptedges.front()*GeV, _ptedges.back()*GeV) && Cuts::absrap < 3.0);
 73      if(jets7.size()>2) jets7.resize(2); // Use only the first two jets
 74      MSG_DEBUG("Jet (R=0.7) multiplicity before cuts = " << jets7.size());
 75      if (jets7.size() == 0) {
 76        MSG_DEBUG("No jets (R=0.7) found in required pT and rapidity range");
 77        vetoEvent;
 78      }
 79      // Calculate and histogram jet shapes
 80      for (size_t jy = 0; jy < 6; ++jy) {
 81        for (size_t ipt = 0; ipt < 22; ++ipt) {
 82          if (ipt > 20 && jy == 3) continue;
 83          if (ipt > 18 && jy >= 4) continue;
 84          JetShape jsipt = apply<JetShape>(evt, _jsnames_pT[ipt][jy]);
 85          jsipt.calc(jets7);
 86          for (size_t ijet = 0; ijet < jsipt.numJets(); ++ijet) {
 87            for (size_t rbin = 0; rbin < jsipt.numBins(); ++rbin) {
 88              const double r_rho = jsipt.rBinMid(rbin);
 89              _profhistRho_pT[ipt][jy]->fill(r_rho, (1./0.1)*jsipt.diffJetShape(ijet, rbin));
 90            }
 91          }
 92        }
 93      }
 94
 95      // Get jets and require at least one to pass pT and y cuts
 96      Jets jets5 = apply<FastJets>(evt, "Jets5")
 97        .jetsByPt(Cuts::ptIn(50*GeV, 1000*GeV) && Cuts::absrap < 2.0);
 98      // Calculate and histogram charged jet shapes
 99      for (const Jet& jet : jets5) {
100        double ncharge = 0;
101        double eta=0;
102        double phi=0;
103        double sumpt=0;
104        for (const Particle& p : jet.particles()) {
105          if ((p.pT() < 0.5) || (p.charge3()==0) || (abs(p.pid())==11) || (abs(p.pid())==13)) continue;
106          ncharge++;
107          sumpt+=p.pT();
108          eta+=p.pT()*p.eta();
109          phi+=p.pT()*mapAngleMPiToPi(p.phi()-jet.phi());
110        }
111        if (jet.absrap()<1.0) {
112          _profhistNch[0]->fill(jet.pT(), ncharge );
113        } else if (jet.absrap()<2.0) {
114          _profhistNch[1]->fill(jet.pT(), ncharge );
115        }
116        if (sumpt==0) continue;
117        eta/=sumpt;
118        phi/=sumpt;
119        double deta=0;
120        double dphi=0;
121        for (const Particle& p : jet.particles()) {
122          if ((p.pT() < 0.5) || (p.charge3()==0) || (abs(p.pid())==11) || (abs(p.pid())==13)) continue;
123          deta+=p.pT()*pow(p.eta()-eta,2);
124          dphi+=p.pT()*pow(mapAngleMPiToPi(p.phi()-phi-jet.phi()),2);
125        }
126        deta/=sumpt;
127        dphi/=sumpt;
128        if ((dphi==0)||(deta==0)) continue;
129        if (jet.absrap()<1.0) {
130          _profhistDr[0]->fill(jet.pT(), deta+dphi );
131          _profhistDeta->fill(jet.pT(), deta );
132          _profhistDphi->fill(jet.pT(), dphi );
133        } else if (jet.absrap()<2.0) {
134          _profhistDr[1]->fill(jet.pT(), deta+dphi );
135        }
136      }
137    }
138
139
140    // Finalize
141    void finalize() {
142      for (size_t i = 1; i < _profhistAsym->numBins()+1; ++i) {
143        if((_profhistDeta->bin(i).effNumEntries()<2)||(_profhistDphi->bin(i).effNumEntries()<2)) continue;
144        if((_profhistDeta->bin(i).yMean()==0)||(_profhistDphi->bin(i).yMean()==0)) continue;
145        double mean_ratio=_profhistDeta->bin(i).yMean() / _profhistDphi->bin(i).yMean();
146        double mean_error=mean_ratio*sqrt(sqr(_profhistDeta->bin(i).yStdErr()/_profhistDeta->bin(i).yMean())
147                                         +sqr(_profhistDphi->bin(i).yStdErr()/_profhistDphi->bin(i).yMean()));
148        _profhistAsym->bin(i).set(mean_ratio,mean_error);
149      }
150    }
151
152    /// @}
153
154
155  private:
156
157    /// @name Analysis data
158    /// @{
159
160    /// Jet \f$ p_\perp\f$ bins.
161    vector<double> _ptedges; // This can't be a raw array if we want to initialise it non-painfully
162    vector<double> _yedges;
163
164    /// JetShape projection name for each \f$p_\perp\f$ bin.
165    string _jsnames_pT[22][6];
166
167    /// @}
168
169    /// @name Histograms
170    /// @{
171    Profile1DPtr _profhistRho_pT[22][6];
172    Profile1DPtr _profhistNch[2];
173    Profile1DPtr _profhistDr[2];
174    Profile1DPtr _profhistDeta;
175    Profile1DPtr _profhistDphi;
176    Estimate1DPtr _profhistAsym;
177    /// @}
178
179  };
180
181
182
183  RIVET_DECLARE_PLUGIN(CMS_2012_I1111014);
184
185}