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      for (size_t i_bin = 0; i_bin < 14; ++i_bin) {
 33        book(_h["nch_B"+to_str(i_bin)] , 13 + i_bin, 1, 1);
 34        book(_hr["r_B"+to_str(i_bin)] , 27 + i_bin, 1, 1);
 35        book(_h["zeta_B"+to_str(i_bin)] , 41 + i_bin, 1, 1);
 36        book(_h["pTrel_B"+to_str(i_bin)] , 55 + i_bin, 1, 1);
 37        book(_h["nch_F"+to_str(i_bin)] , 69 + i_bin, 1, 1);
 38        book(_hr["r_F"+to_str(i_bin)] , 83 + i_bin, 1, 1);
 39        book(_h["zeta_F"+to_str(i_bin)] , 97 + i_bin, 1, 1);
 40        book(_h["pTrel_F"+to_str(i_bin)] , 111 + i_bin, 1, 1);
 41        book(_h["nch_C"+to_str(i_bin)] , 125 + i_bin, 1, 1);
 42        book(_hr["r_C"+to_str(i_bin)] , 139 + i_bin, 1, 1);
 43        book(_h["zeta_C"+to_str(i_bin)] , 153 + i_bin, 1, 1);
 44        book(_h["pTrel_C"+to_str(i_bin)] , 167 + i_bin, 1, 1);
 45        }
 46    }
 47
 48    void analyze(const Event& event) {
 49
 50      //Init
 51      double fnch=0;
 52      double cnch=0;
 53      double fzval=0;
 54      double czval=0;
 55      double frval=0;
 56      double crval=0;
 57      double ftval=0;
 58      double ctval=0;
 59
 60      // Event selection
 61      Jets m_goodJets = apply<JetFinder>(event, "Jets").jetsByPt(Cuts::pT > 100*GeV && Cuts::abseta < 2.1);
 62      if (m_goodJets.size() < 2) vetoEvent;
 63      if (fabs(1.0 - m_goodJets[0].pT()/m_goodJets[1].pT()) > 0.5)  vetoEvent;
 64      // Decide forward or central
 65      bool check = m_goodJets[0].abseta() < m_goodJets[1].abseta();
 66      int pos_f = int(check);
 67      int pos_c = int(!check);
 68
 69      // Calculate obs, separately for central and fwd, also get bin
 70      double fpt = m_goodJets[pos_f].pT();
 71      size_t f_bin = GetJetBin(fpt);
 72      double cpt = m_goodJets[pos_c].pT();
 73      size_t c_bin = GetJetBin(cpt);
 74
 75      for (const Particle& p : m_goodJets[pos_f].particles()) {
 76        if (p.pT() < 0.5*GeV)  continue;
 77        if (p.charge() != 0) {
 78          ++fnch;
 79
 80          fzval = p.pT() / m_goodJets[pos_f].pt();
 81          ftval = p.pT()*sin(p.phi()-m_goodJets[pos_f].phi());
 82          frval = deltaR(m_goodJets[pos_f], p);
 83
 84          _hr["r_F"+to_str(f_bin)]->fill(frval);
 85          _hr["r_B"+to_str(f_bin)]->fill(frval);
 86          _h["zeta_F"+to_str(f_bin)]->fill(fzval);
 87          _h["zeta_B"+to_str(f_bin)]->fill(fzval);
 88          _h["pTrel_F"+to_str(f_bin)]->fill(ftval);
 89          _h["pTrel_B"+to_str(f_bin)]->fill(ftval);
 90        }
 91      }
 92
 93
 94	    for (const Particle& p : m_goodJets[pos_c].particles()) {
 95        if (p.pT() < 0.5*GeV)  continue;
 96        if (p.charge() != 0) {
 97          ++cnch;
 98          czval = p.pT() / m_goodJets[pos_c].pt();
 99          ctval = p.pT()*sin(p.phi()-m_goodJets[pos_c].phi());
100          crval = deltaR(m_goodJets[pos_c], p);
101
102          _hr["r_C"+to_str(c_bin)]->fill(crval);
103          _hr["r_B"+to_str(c_bin)]->fill(crval);
104          _h["zeta_C"+to_str(c_bin)]->fill(czval);
105          _h["zeta_B"+to_str(c_bin)]->fill(czval);
106          _h["pTrel_C"+to_str(c_bin)]->fill(ctval);
107          _h["pTrel_B"+to_str(c_bin)]->fill(ctval);
108        }
109      }
110
111      if (fnch > 63)  fnch = 63;
112      if (cnch > 63)  cnch = 63;
113
114       //Fill nchg histo
115
116      _p["nch_jetpt_F"]->fill(fpt,fnch);
117      _p["nch_jetpt_C"]->fill(cpt,cnch);
118      _p["nch_jetpt_B"]->fill(fpt,fnch);
119      _p["nch_jetpt_B"]->fill(cpt,cnch);
120
121      _h["nch_F"+to_str(f_bin)]->fill(fnch);
122      _h["nch_C"+to_str(c_bin)]->fill(cnch);
123      _h["nch_B"+to_str(f_bin)]->fill(fnch);
124      _h["nch_B"+to_str(c_bin)]->fill(cnch);
125    }
126
127
128    void finalize() {
129
130      for (size_t i_bin = 0; i_bin < 14; ++i_bin) {
131
132        const double sfB =  _h["nch_B"+to_str(i_bin)]->sumW();
133        if (sfB) {
134          scale(_h["zeta_B"+to_str(i_bin)],2.0/sfB);
135          scale(_h["pTrel_B"+to_str(i_bin)],2.0/sfB);
136          scale(_hr["r_B"+to_str(i_bin)],2.0/sfB);
137          scale(_h["nch_B"+to_str(i_bin)], 2.0/sfB);
138        }
139
140        const double sfF =  _h["nch_F"+to_str(i_bin)]->sumW();
141        if (sfF) {
142          scale(_h["zeta_F"+to_str(i_bin)],1.0/sfF);
143          scale(_h["pTrel_F"+to_str(i_bin)],1.0/sfF);
144          scale(_hr["r_F"+to_str(i_bin)],1.0/sfF);
145          scale(_h["nch_F"+to_str(i_bin)], 1.0/sfF);
146        }
147
148        const double sfC =  _h["nch_C"+to_str(i_bin)]->sumW();
149        if (sfC) {
150          scale(_h["zeta_C"+to_str(i_bin)],1.0/sfC);
151          scale(_h["pTrel_C"+to_str(i_bin)],1.0/sfC);
152          scale(_hr["r_C"+to_str(i_bin)],1.0/sfC);
153          scale(_h["nch_C"+to_str(i_bin)], 1.0/sfC);
154        }
155      }
156
157      // For r only
158      /// @todo Replace with barchart()
159      for (auto& hist : _hr) {
160        for(size_t i=1; i < hist.second->numBins()+1; ++i) {
161          const double x = hist.second->bin(i).xMid();
162          const double bW = hist.second->bin(i).xWidth();
163          hist.second->bin(i).scaleW(bW/(2.0*M_PI*x));
164        }
165      }
166
167      // The rest
168      /// @todo Replace with barchart()
169      for (auto& hist : _h) {
170        for (size_t i=1; i < hist.second->numBins()+1; ++i) {
171          const double bW = hist.second->bin(i).xWidth();
172          hist.second->bin(i).scaleW(bW);
173        }
174      }
175
176    }
177
178  private:
179
180    size_t GetJetBin(const double jetpt){
181      size_t i_bin = 0;
182      if (inRange(jetpt,100,200)) i_bin=0;
183      if (inRange(jetpt,200,300)) i_bin=1;
184      if (inRange(jetpt,300,400)) i_bin=2;
185      if (inRange(jetpt,400,500)) i_bin=3;
186      if (inRange(jetpt,500,600)) i_bin=4;
187      if (inRange(jetpt,600,700)) i_bin=5;
188      if (inRange(jetpt,700,800)) i_bin=6;
189      if (inRange(jetpt,800,900)) i_bin=7;
190      if (inRange(jetpt,900,1000)) i_bin=8;
191      if (inRange(jetpt,1000,1200)) i_bin=9;
192      if (inRange(jetpt,1200,1400)) i_bin=10;
193      if (inRange(jetpt,1400,1600)) i_bin=11;
194      if (inRange(jetpt,1600,2000)) i_bin=12;
195      if (inRange(jetpt,2000,2500)) i_bin=13;
196      if(jetpt < 100) i_bin=0;
197      if(jetpt > 2500) i_bin=13;
198      return i_bin;
199    }
200
201    map<string, Histo1DPtr> _h;
202    map<string, Histo1DPtr> _hr;
203    map<string, Profile1DPtr> _p;
204  };
205
206  RIVET_DECLARE_PLUGIN(ATLAS_2019_I1740909);
207
208}