rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

ATLAS_2019_I1740909

Jet fragmentation using charged particles
Experiment: ATLAS (LHC)
Inspire ID: 1740909
Status: VALIDATED
Authors:
  • Ben Nachman
  • Deepak Kar
References: Beams: p+ p+
Beam energies: (6500.0, 6500.0) GeV
Run details:
  • p + p -> jj

This paper presents a measurement of quantities related to the formation of jets from high-energy quarks and gluons (fragmentation). Jets with transverse momentum 100 GeV $< p_T <$ 2.5 TeV and pseudorapidity $|\eta|<2.1$ from an integrated luminosity of 33 fb$^{-1}$ of $\sqrt{s}=$13 TeV proton-proton collisions are reconstructed with the ATLAS detector at the Large Hadron Collider. Charged-particle tracks with $p_T >$ 500 MeV and $|\eta|<2.5$ are used to probe the detailed structure of the jet. The fragmentation properties of the more forward and the more central of the two leading jets from each event are studied. The data are unfolded to correct for detector resolution and acceptance effects. Comparisons with parton shower Monte Carlo generators indicate that existing models provide a reasonable description of the data across a wide range of phase space, but there are also significant differences. Furthermore, the data are interpreted in the context of quark- and gluon-initiated jets by exploiting the rapidity dependence of the jet flavor fraction. A first measurement of the charged-particle multiplicity using model-independent jet labels (topic modeling) provides a promising alternative to traditional quark and gluon extractions using input from simulation. The simulations provide a reasonable description of the quark-like data across the jet $p_T$ range presented in this measurement, but the gluon-like data have systematically fewer charged particles than the simulations.

Source code: ATLAS_2019_I1740909.cc
  1#include "Rivet/Analysis.hh"
  2#include "Rivet/Projections/FinalState.hh"
  3#include "Rivet/Projections/VetoedFinalState.hh"
  4#include "Rivet/Projections/InvisibleFinalState.hh"
  5#include "Rivet/Projections/PromptFinalState.hh"
  6#include "Rivet/Projections/FastJets.hh"
  7
  8namespace Rivet {
  9
 10  /// @brief jet fragmentation at 13 TeV
 11  class ATLAS_2019_I1740909: public Analysis {
 12  public:
 13
 14    RIVET_DEFAULT_ANALYSIS_CTOR(ATLAS_2019_I1740909);
 15
 16    /// Book cuts and projections
 17    void init() {
 18
 19      const FinalState bare_MU(Cuts::abspid == PID::MUON);
 20
 21      VetoedFinalState jetinput;
 22      jetinput.addVetoOnThisFinalState(bare_MU);
 23      jetinput.addVetoOnThisFinalState(InvisibleFinalState());
 24
 25      FastJets jetpro(jetinput, JetAlg::ANTIKT, 0.4);
 26      declare(jetpro, "Jets");
 27
 28      book(_p["nch_jetpt_F"], 1, 1, 1);
 29      book(_p["nch_jetpt_C"], 2, 1, 1);
 30      book(_p["nch_jetpt_B"], 9, 1, 1);
 31
 32      vector<double> ptbins{100., 200., 300., 400., 500., 600., 700., 800.,
 33                            900., 1000., 1200., 1400., 1600., 2000., 2500.};
 34      for (const string& suff : vector<string>{"B", "F", "C"}) {
 35        book(_h["nch_"+suff], ptbins);
 36        book(_h["r_"+suff], ptbins);
 37        book(_h["zeta_"+suff], ptbins);
 38        book(_h["pTrel_"+suff], ptbins);
 39      }
 40      for (size_t idx = 0; idx < ptbins.size()-1; ++idx) {
 41        dualbook(_h["nch_B"]->bin(idx+1),    13 + idx, 1, 1);
 42        dualbook(_h["r_B"]->bin(idx+1),      27 + idx, 1, 1);
 43        dualbook(_h["zeta_B"]->bin(idx+1),   41 + idx, 1, 1);
 44        dualbook(_h["pTrel_B"]->bin(idx+1),  55 + idx, 1, 1);
 45        dualbook(_h["nch_F"]->bin(idx+1),    69 + idx, 1, 1);
 46        dualbook(_h["r_F"]->bin(idx+1),      83 + idx, 1, 1);
 47        dualbook(_h["zeta_F"]->bin(idx+1),   97 + idx, 1, 1);
 48        dualbook(_h["pTrel_F"]->bin(idx+1), 111 + idx, 1, 1);
 49        dualbook(_h["nch_C"]->bin(idx+1),   125 + idx, 1, 1);
 50        dualbook(_h["r_C"]->bin(idx+1),     139 + idx, 1, 1);
 51        dualbook(_h["zeta_C"]->bin(idx+1),  153 + idx, 1, 1);
 52        dualbook(_h["pTrel_C"]->bin(idx+1), 167 + idx, 1, 1);
 53      }
 54    }
 55
 56    void dualbook(Histo1DPtr& hist, unsigned int d, unsigned int x, unsigned int y) {
 57      const string hname = mkAxisCode(d, x, y);
 58      book(hist, "_aux_"+hname, refData(hname));
 59      book(_s[hist.get()->basePath()], hname);
 60    }
 61
 62    void analyze(const Event& event) {
 63
 64      //Init
 65      double fnch=0;
 66      double cnch=0;
 67      double fzval=0;
 68      double czval=0;
 69      double frval=0;
 70      double crval=0;
 71      double ftval=0;
 72      double ctval=0;
 73
 74      // Event selection
 75      Jets m_goodJets = apply<JetFinder>(event, "Jets").jetsByPt(Cuts::pT > 100*GeV && Cuts::abseta < 2.1);
 76      if (m_goodJets.size() < 2) vetoEvent;
 77      if (fabs(1.0 - m_goodJets[0].pT()/m_goodJets[1].pT()) > 0.5)  vetoEvent;
 78      // Decide forward or central
 79      bool check = m_goodJets[0].abseta() < m_goodJets[1].abseta();
 80      int pos_f = int(check);
 81      int pos_c = int(!check);
 82
 83      // Calculate obs, separately for central and fwd, also get bin
 84      double fpt = m_goodJets[pos_f].pT();
 85      double cpt = m_goodJets[pos_c].pT();
 86
 87      for (const Particle& p : m_goodJets[pos_f].particles()) {
 88        if (p.pT() < 0.5*GeV)  continue;
 89        if (p.charge() != 0) {
 90          ++fnch;
 91
 92          fzval = p.pT() / m_goodJets[pos_f].pt();
 93          ftval = p.pT()*sin(p.phi()-m_goodJets[pos_f].phi());
 94          frval = deltaR(m_goodJets[pos_f], p);
 95
 96          _h["r_F"]->fill(fpt, frval);
 97          _h["r_B"]->fill(fpt, frval);
 98          _h["zeta_F"]->fill(fpt, fzval);
 99          _h["zeta_B"]->fill(fpt, fzval);
100          _h["pTrel_F"]->fill(fpt, ftval);
101          _h["pTrel_B"]->fill(fpt, ftval);
102        }
103      }
104
105
106	    for (const Particle& p : m_goodJets[pos_c].particles()) {
107        if (p.pT() < 0.5*GeV)  continue;
108        if (p.charge() != 0) {
109          ++cnch;
110          czval = p.pT() / m_goodJets[pos_c].pt();
111          ctval = p.pT()*sin(p.phi()-m_goodJets[pos_c].phi());
112          crval = deltaR(m_goodJets[pos_c], p);
113
114          _h["r_C"]->fill(cpt, crval);
115          _h["r_B"]->fill(cpt, crval);
116          _h["zeta_C"]->fill(cpt, czval);
117          _h["zeta_B"]->fill(cpt, czval);
118          _h["pTrel_C"]->fill(cpt, ctval);
119          _h["pTrel_B"]->fill(cpt, ctval);
120        }
121      }
122
123      if (fnch > 63)  fnch = 63;
124      if (cnch > 63)  cnch = 63;
125
126       //Fill nchg histo
127
128      _p["nch_jetpt_F"]->fill(fpt,fnch);
129      _p["nch_jetpt_C"]->fill(cpt,cnch);
130      _p["nch_jetpt_B"]->fill(fpt,fnch);
131      _p["nch_jetpt_B"]->fill(cpt,cnch);
132
133      _h["nch_F"]->fill(fpt, fnch);
134      _h["nch_C"]->fill(cpt, cnch);
135      _h["nch_B"]->fill(fpt, fnch);
136      _h["nch_B"]->fill(cpt, cnch);
137    }
138
139
140    void finalize() {
141
142      for (const string& suff : vector<string>{"B", "F", "C"}) {
143        const double num = suff=="B"? 2.0 : 1.0;
144        vector<double> sf = _h["nch_"+suff]->sumWGroup();
145        for (double& f : sf) {
146          f = safediv(num, f, 0.0);
147        }
148        scale({_h["nch_"+suff], _h["zeta_"+suff], _h["pTrel_"+suff], _h["r_"+suff]}, sf);
149        for (auto& r : _h["r_"+suff]->bins()) {
150          for (auto& b : r->bins()) {
151            b.scaleW(1.0/(2*M_PI*b.xMid()));
152          }
153        }
154      }
155
156      for (auto& hist : _h) {
157        for (auto& b : hist.second->bins()) {
158          barchart(b, _s[b.get()->basePath()]);
159        }
160      }
161
162    }
163
164  private:
165
166    map<string, Histo1DGroupPtr> _h;
167    map<string, Estimate1DPtr> _s;
168    map<string, Profile1DPtr> _p;
169  };
170
171  RIVET_DECLARE_PLUGIN(ATLAS_2019_I1740909);
172
173}