Rivet Analyses Reference

ATLAS_2010_CONF_2010_049

Cross-section and fragmentation function in anti-$k_t$ track jets
Experiment: ATLAS (LHC 7000GeV)
Spires ID: None
Status: OBSOLETE
Authors:
  • Hendrik Hoeth
References:
  • ATLAS-CONF-2010-049
Beams: p+ p+
Beam energies: (3500.0, 3500.0) GeV
Run details:
  • pp QCD interactions at 7000 GeV including diffractive events.

Jets are identified and their properties studied using tracks measured by the ATLAS Inner Detector. Events are selected using a minimum-bias trigger, allowing the emergence of jets at low transverse momentum to be observed and for jets to be studied independently of the calorimeter. Jets are reconstructed using the anti-kt algorithm applied to tracks with two parameter choices, 0.4 and 0.6. An inclusive jet transverse momentum cross section measurement from 4 GeV to 80 GeV is shown, integrated over $|\eta| < 0.57$ and corrected to charged particle-level truth jets. The probability that a particular particle carries a fixed fraction of the jet momentum (fragmentation function) is also measured. All data is corrected to the particle level. ATTENTION - Data read from plots!

Source code: ATLAS_2010_CONF_2010_049.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
// -*- C++ -*-
#include "Rivet/Analysis.hh"
#include "Rivet/Projections/ChargedFinalState.hh"
#include "Rivet/Projections/FastJets.hh"

namespace Rivet {


  class ATLAS_2010_CONF_2010_049 : public Analysis {
  public:

    ATLAS_2010_CONF_2010_049()
      : Analysis("ATLAS_2010_CONF_2010_049")
    {    }


    void init() {
      ChargedFinalState cfs((Cuts::etaIn(-1.5, 1.5) && Cuts::pT >=  0.5*GeV));
      declare(cfs, "CFS");

      FastJets jetsproj6(cfs, FastJets::ANTIKT, 0.6);
      declare(jetsproj6, "Jets6");

      FastJets jetsproj4(cfs, FastJets::ANTIKT, 0.4);
      declare(jetsproj4, "Jets4");

      // @todo tmp YOs
      for (size_t i=0 ; i<2 ; i++) {
        book(_h_xsec[i]       ,1+i, 1, 1);
        book(_h_frag_04_06[i] ,3+i, 1, 1);
        book(_h_frag_06_10[i] ,3+i, 2, 1);
        book(_h_frag_10_15[i] ,3+i, 3, 1);
        book(_h_frag_15_24[i] ,3+i, 4, 1);
        book(_njets_04_06[i], "njets_04_06_"+to_string(i));
        book(_njets_06_10[i], "njets_06_10_"+to_string(i));
        book(_njets_10_15[i], "njets_10_15_"+to_string(i));
        book(_njets_15_24[i], "njets_15_24_"+to_string(i));
      }
    }


    void analyze(const Event& event) {
      const FastJets & jetsproj6 = apply<FastJets>(event, "Jets6");
      const FastJets & jetsproj4 = apply<FastJets>(event, "Jets4");
      Jets alljets[2];
      alljets[0] = jetsproj6.jetsByPt(4.0*GeV);
      alljets[1] = jetsproj4.jetsByPt(4.0*GeV);

      for (size_t i=0 ; i<2 ; i++) {
        Jets jets;

        // First we want to make sure that we only use jets within |eta|<0.57
        for (const Jet& jet : alljets[i]) {
          if (jet.abseta()<0.57) {
            jets.push_back(jet);
          }
        }
        for (const Jet& jet : jets) {
          const double pTjet = jet.pT();
          const double pjet = jet.p3().mod();
          _h_xsec[i]->fill(pTjet);
          if (pTjet > 24*GeV) continue;
          for (const Particle& p : jet.particles()) {
            double z = p.p3().mod()/pjet;
            if (z >= 1) z = 0.9999; // Make sure that z=1 doesn't go into overflow
            if (pTjet > 15*GeV) {
              _h_frag_15_24[i]->fill(z);
            }
            else if (pTjet > 10*GeV) {
              _h_frag_10_15[i]->fill(z);
            }
            else if (pTjet > 6*GeV) {
              _h_frag_06_10[i]->fill(z);
            }
            else {
              _h_frag_04_06[i]->fill(z);
            }
          }
          if (pTjet > 15*GeV) {
            _njets_15_24[i]->fill();
          }
          else if (pTjet > 10*GeV) {
            _njets_10_15[i]->fill();
          }
          else if (pTjet > 6*GeV) {
            _njets_06_10[i]->fill();
          }
          else {
            _njets_04_06[i]->fill();
          }
        }
      }
    }

    void finalize() {
      for (size_t i=0 ; i<2 ; i++) {
        // deta = 2*0.57
        scale(_h_xsec[i], crossSection()/microbarn/sumOfWeights()/(2*0.57));
        scale(_h_frag_04_06[i], 1./_njets_04_06[i]->val());
        scale(_h_frag_06_10[i], 1./_njets_06_10[i]->val());
        scale(_h_frag_10_15[i], 1./_njets_10_15[i]->val());
        scale(_h_frag_15_24[i], 1./_njets_15_24[i]->val());
      }
    }


  private:

    Histo1DPtr _h_xsec[2];
    Histo1DPtr _h_frag_04_06[2];
    Histo1DPtr _h_frag_06_10[2];
    Histo1DPtr _h_frag_10_15[2];
    Histo1DPtr _h_frag_15_24[2];
    CounterPtr _njets_04_06[2];
    CounterPtr _njets_06_10[2];
    CounterPtr _njets_10_15[2];
    CounterPtr _njets_15_24[2];
  };



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

}