rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

ATLAS_2012_I1082009

$D^{*\pm}$ production in jets
Experiment: ATLAS (LHC)
Inspire ID: 1082009
Status: VALIDATED
Authors:
  • Peter Richardson
References: Beams: p+ p+
Beam energies: (3500.0, 3500.0) GeV
Run details:
  • All flavours of quark and gluon jet production at 7 TeV

Measurement of $D^{*\pm}$ meson production in jets from proton-proton collisions at a centre-of-mass energy of $\sqrt{s}=7$ TeV at the LHC. The measurement is based on a data sample recorded with the ATLAS detector with an integrated luminosity of $0.30\,\text{pb}^{-1}$ for jets with transverse momentum between 25 and 70 GeV in the pseudorapidity range $|eta| < 2.5$.

Source code: ATLAS_2012_I1082009.cc
  1// -*- C++ -*-
  2#include "Rivet/Analysis.hh"
  3#include "Rivet/Projections/IdentifiedFinalState.hh"
  4#include "Rivet/Projections/VetoedFinalState.hh"
  5#include "Rivet/Projections/MissingMomentum.hh"
  6#include "Rivet/Projections/FastJets.hh"
  7#include "Rivet/Projections/LeadingParticlesFinalState.hh"
  8#include "Rivet/Projections/UnstableParticles.hh"
  9
 10namespace Rivet {
 11
 12
 13  class ATLAS_2012_I1082009 : public Analysis {
 14  public:
 15
 16    /// @name Constructors etc.
 17    /// @{
 18
 19    /// Constructor
 20    ATLAS_2012_I1082009()
 21      : Analysis("ATLAS_2012_I1082009")
 22    {    }
 23
 24    /// @}
 25
 26
 27  public:
 28
 29    /// @name Analysis methods
 30    /// @{
 31
 32    /// Book histograms and initialise projections before the run
 33    void init() {
 34
 35      // Input for the jets: No neutrinos, no muons
 36      VetoedFinalState veto;
 37      veto.addVetoPairId(PID::MUON);
 38      veto.vetoNeutrinos();
 39      FastJets jets(veto, JetAlg::ANTIKT, 0.6);
 40      declare(jets, "jets");
 41      // unstable final-state for D*
 42      declare(UnstableParticles(), "UFS");
 43
 44      book(_weight25_30, "_weight_25_30");
 45      book(_weight30_40, "_weight_30_40");
 46      book(_weight40_50, "_weight_40_50");
 47      book(_weight50_60, "_weight_50_60");
 48      book(_weight60_70, "_weight_60_70");
 49      book(_weight25_70, "_weight_25_70");
 50
 51      book(_h_pt25_30 , 8,1,1);
 52      book(_h_pt30_40 , 9,1,1);
 53      book(_h_pt40_50 ,10,1,1);
 54      book(_h_pt50_60 ,11,1,1);
 55      book(_h_pt60_70 ,12,1,1);
 56      book(_h_pt25_70 ,13,1,1);
 57    }
 58
 59
 60    /// Perform the per-event analysis
 61    void analyze(const Event& event) {
 62      // get the jets
 63      Jets jets = apply<FastJets>(event, "jets").jetsByPt(Cuts::pT > 25*GeV && Cuts::abseta < 2.5);
 64      // get the D* mesons
 65      const UnstableParticles& ufs = apply<UnstableParticles>(event, "UFS");
 66      Particles Dstar;
 67      for (const Particle& p : ufs.particles()) {
 68        const int id = p.abspid();
 69        if(id==413) Dstar.push_back(p);
 70      }
 71
 72      // loop over the jobs
 73      for (const Jet& jet : jets ) {
 74        double perp = jet.perp();
 75        bool found = false;
 76        double z(0.);
 77        if(perp<25.||perp>70.) continue;
 78        for(const Particle & p : Dstar) {
 79          if(p.perp()<7.5) continue;
 80          if(deltaR(p, jet.momentum())<0.6) {
 81            Vector3 axis = jet.p3().unit();
 82            z = axis.dot(p.p3())/jet.E();
 83            if(z<0.3) continue;
 84            found = true;
 85            break;
 86          }
 87        }
 88        _weight25_70->fill();
 89        if(found) _h_pt25_70->fill(z);
 90        if(perp>=25.&&perp<30.) {
 91          _weight25_30->fill();
 92          if(found) _h_pt25_30->fill(z);
 93        }
 94        else if(perp>=30.&&perp<40.) {
 95          _weight30_40->fill();
 96          if(found) _h_pt30_40->fill(z);
 97        }
 98        else if(perp>=40.&&perp<50.) {
 99          _weight40_50->fill();
100          if(found) _h_pt40_50->fill(z);
101        }
102        else if(perp>=50.&&perp<60.) {
103          _weight50_60->fill();
104          if(found) _h_pt50_60->fill(z);
105        }
106        else if(perp>=60.&&perp<70.) {
107          _weight60_70->fill();
108          if(found) _h_pt60_70->fill(z);
109        }
110      }
111    }
112
113
114    /// Normalise histograms etc., after the run
115    void finalize() {
116      scale(_h_pt25_30,1./ *_weight25_30);
117      scale(_h_pt30_40,1./ *_weight30_40);
118      scale(_h_pt40_50,1./ *_weight40_50);
119      scale(_h_pt50_60,1./ *_weight50_60);
120      scale(_h_pt60_70,1./ *_weight60_70);
121      scale(_h_pt25_70,1./ *_weight25_70);
122    }
123
124    /// @}
125
126
127  private:
128
129    /// @name Histograms
130    /// @{
131    CounterPtr _weight25_30,_weight30_40,_weight40_50;
132    CounterPtr _weight50_60,_weight60_70,_weight25_70;
133
134    Histo1DPtr _h_pt25_30;
135    Histo1DPtr _h_pt30_40;
136    Histo1DPtr _h_pt40_50;
137    Histo1DPtr _h_pt50_60;
138    Histo1DPtr _h_pt60_70;
139    Histo1DPtr _h_pt25_70;
140    /// @}
141
142  };
143
144  RIVET_DECLARE_PLUGIN(ATLAS_2012_I1082009);
145
146}