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
  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
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
// -*- C++ -*-
#include "Rivet/Analysis.hh"
#include "Rivet/Projections/IdentifiedFinalState.hh"
#include "Rivet/Projections/VetoedFinalState.hh"
#include "Rivet/Projections/MissingMomentum.hh"
#include "Rivet/Projections/FastJets.hh"
#include "Rivet/Projections/LeadingParticlesFinalState.hh"
#include "Rivet/Projections/UnstableParticles.hh"

namespace Rivet {


  class ATLAS_2012_I1082009 : public Analysis {
  public:

    /// @name Constructors etc.
    //@{

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

    //@}


  public:

    /// @name Analysis methods
    //@{

    /// Book histograms and initialise projections before the run
    void init() {

      // Input for the jets: No neutrinos, no muons
      VetoedFinalState veto;
      veto.addVetoPairId(PID::MUON);
      veto.vetoNeutrinos();
      FastJets jets(veto, FastJets::ANTIKT, 0.6);
      declare(jets, "jets");
      // unstable final-state for D*
      declare(UnstableParticles(), "UFS");

      book(_weight25_30, "_weight_25_30");
      book(_weight30_40, "_weight_30_40");
      book(_weight40_50, "_weight_40_50");
      book(_weight50_60, "_weight_50_60");
      book(_weight60_70, "_weight_60_70");
      book(_weight25_70, "_weight_25_70");

      book(_h_pt25_30 , 8,1,1);
      book(_h_pt30_40 , 9,1,1);
      book(_h_pt40_50 ,10,1,1);
      book(_h_pt50_60 ,11,1,1);
      book(_h_pt60_70 ,12,1,1);
      book(_h_pt25_70 ,13,1,1);
    }


    /// Perform the per-event analysis
    void analyze(const Event& event) {
      // get the jets
      Jets jets;
      for (const Jet& jet : apply<FastJets>(event, "jets").jetsByPt(25.0*GeV)) {
        if ( jet.abseta() < 2.5 ) jets.push_back(jet);
      }
      // get the D* mesons
      const UnstableParticles& ufs = apply<UnstableParticles>(event, "UFS");
      Particles Dstar;
      for (const Particle& p : ufs.particles()) {
        const int id = p.abspid();
        if(id==413) Dstar.push_back(p);
      }

      // loop over the jobs
      for (const Jet& jet : jets ) {
        double perp = jet.perp();
        bool found = false;
        double z(0.);
        if(perp<25.||perp>70.) continue;
        for(const Particle & p : Dstar) {
          if(p.perp()<7.5) continue;
          if(deltaR(p, jet.momentum())<0.6) {
            Vector3 axis = jet.p3().unit();
            z = axis.dot(p.p3())/jet.E();
            if(z<0.3) continue;
            found = true;
            break;
          }
        }
        _weight25_70->fill();
        if(found) _h_pt25_70->fill(z);
        if(perp>=25.&&perp<30.) {
          _weight25_30->fill();
          if(found) _h_pt25_30->fill(z);
        }
        else if(perp>=30.&&perp<40.) {
          _weight30_40->fill();
          if(found) _h_pt30_40->fill(z);
        }
        else if(perp>=40.&&perp<50.) {
          _weight40_50->fill();
          if(found) _h_pt40_50->fill(z);
        }
        else if(perp>=50.&&perp<60.) {
          _weight50_60->fill();
          if(found) _h_pt50_60->fill(z);
        }
        else if(perp>=60.&&perp<70.) {
          _weight60_70->fill();
          if(found) _h_pt60_70->fill(z);
        }
      }
    }


    /// Normalise histograms etc., after the run
    void finalize() {
      scale(_h_pt25_30,1./ *_weight25_30);
      scale(_h_pt30_40,1./ *_weight30_40);
      scale(_h_pt40_50,1./ *_weight40_50);
      scale(_h_pt50_60,1./ *_weight50_60);
      scale(_h_pt60_70,1./ *_weight60_70);
      scale(_h_pt25_70,1./ *_weight25_70);
    }

    //@}


  private:

    /// @name Histograms
    //@{
    CounterPtr _weight25_30,_weight30_40,_weight40_50;
    CounterPtr _weight50_60,_weight60_70,_weight25_70;

    Histo1DPtr _h_pt25_30;
    Histo1DPtr _h_pt30_40;
    Histo1DPtr _h_pt40_50;
    Histo1DPtr _h_pt50_60;
    Histo1DPtr _h_pt60_70;
    Histo1DPtr _h_pt25_70;
    //@}

  };

  // The hook for the plugin system
  RIVET_DECLARE_PLUGIN(ATLAS_2012_I1082009);

}