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
  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
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
#include "Rivet/Analysis.hh"
#include "Rivet/Projections/FinalState.hh"
#include "Rivet/Projections/VetoedFinalState.hh"
#include "Rivet/Projections/InvisibleFinalState.hh"
#include "Rivet/Projections/PromptFinalState.hh"
#include "Rivet/Projections/FastJets.hh"

namespace Rivet {

  /// @brief jet fragmentation at 13 TeV
  class ATLAS_2019_I1740909: public Analysis {
  public:

    RIVET_DEFAULT_ANALYSIS_CTOR(ATLAS_2019_I1740909);
    
    /// Book cuts and projections
    void init() {

      const FinalState bare_MU(Cuts::abspid == PID::MUON);

      VetoedFinalState jetinput;
      jetinput.addVetoOnThisFinalState(bare_MU);
      jetinput.addVetoOnThisFinalState(InvisibleFinalState());

      FastJets jetpro(jetinput, FastJets::ANTIKT, 0.4);
      declare(jetpro, "Jets");

      book(_p["nch_jetpt_F"] , 1, 1, 1);
      book(_p["nch_jetpt_C"] , 2, 1, 1);
      book(_p["nch_jetpt_B"] , 9, 1, 1);
    
      for (size_t i_bin = 0; i_bin < 14; ++i_bin) {
        book(_h["nch_B"+to_str(i_bin)] , 13 + i_bin, 1, 1); 
        book(_hr["r_B"+to_str(i_bin)] , 27 + i_bin, 1, 1); 
        book(_h["zeta_B"+to_str(i_bin)] , 41 + i_bin, 1, 1); 
        book(_h["pTrel_B"+to_str(i_bin)] , 55 + i_bin, 1, 1); 
        book(_h["nch_F"+to_str(i_bin)] , 69 + i_bin, 1, 1); 
        book(_hr["r_F"+to_str(i_bin)] , 83 + i_bin, 1, 1); 
        book(_h["zeta_F"+to_str(i_bin)] , 97 + i_bin, 1, 1); 
        book(_h["pTrel_F"+to_str(i_bin)] , 111 + i_bin, 1, 1);  
        book(_h["nch_C"+to_str(i_bin)] , 125 + i_bin, 1, 1); 
        book(_hr["r_C"+to_str(i_bin)] , 139 + i_bin, 1, 1); 
        book(_h["zeta_C"+to_str(i_bin)] , 153 + i_bin, 1, 1); 
        book(_h["pTrel_C"+to_str(i_bin)] , 167 + i_bin, 1, 1);          
        }
    }

    void analyze(const Event& event) {
    
      //Init
      double fnch=0;
      double cnch=0;
      double fzval=0;
      double czval=0;
      double frval=0;
      double crval=0;
      double ftval=0;
      double ctval=0;

      // Event selection 
      Jets m_goodJets = apply<JetAlg>(event, "Jets").jetsByPt(Cuts::pT > 100*GeV && Cuts::abseta < 2.1);
      if (m_goodJets.size() < 2) vetoEvent;
      if (fabs(1.0 - m_goodJets[0].pT()/m_goodJets[1].pT()) > 0.5)  vetoEvent;
      // Decide forward or central   
      bool check = m_goodJets[0].abseta() < m_goodJets[1].abseta();
      int pos_f = int(check);
      int pos_c = int(!check);
      
      // Calculate obs, separately for central and fwd, also get bin
      double fpt = m_goodJets[pos_f].pT();
      size_t f_bin = GetJetBin(fpt);
      double cpt = m_goodJets[pos_c].pT();
      size_t c_bin = GetJetBin(cpt);

      for (const Particle& p : m_goodJets[pos_f].particles()) {
        if (p.pT() < 0.5*GeV)  continue;
        if (p.charge() != 0) {
          ++fnch;
	                
          fzval = p.pT() / m_goodJets[pos_f].pt();        
          ftval = p.pT()*sin(p.phi()-m_goodJets[pos_f].phi());                
          frval = deltaR(m_goodJets[pos_f], p);
               
          _hr["r_F"+to_str(f_bin)]->fill(frval);
          _hr["r_B"+to_str(f_bin)]->fill(frval);                               
          _h["zeta_F"+to_str(f_bin)]->fill(fzval);
          _h["zeta_B"+to_str(f_bin)]->fill(fzval);     
          _h["pTrel_F"+to_str(f_bin)]->fill(ftval);
          _h["pTrel_B"+to_str(f_bin)]->fill(ftval);
        }    
      }
	               
	               
	    for (const Particle& p : m_goodJets[pos_c].particles()) {            
        if (p.pT() < 0.5*GeV)  continue;
        if (p.charge() != 0) {
          ++cnch;
          czval = p.pT() / m_goodJets[pos_c].pt();
          ctval = p.pT()*sin(p.phi()-m_goodJets[pos_c].phi());
          crval = deltaR(m_goodJets[pos_c], p);
                       
          _hr["r_C"+to_str(c_bin)]->fill(crval);
          _hr["r_B"+to_str(c_bin)]->fill(crval);
          _h["zeta_C"+to_str(c_bin)]->fill(czval);
          _h["zeta_B"+to_str(c_bin)]->fill(czval);
          _h["pTrel_C"+to_str(c_bin)]->fill(ctval);
          _h["pTrel_B"+to_str(c_bin)]->fill(ctval);
        }
      }

      if (fnch > 63)  fnch = 63;
      if (cnch > 63)  cnch = 63;
       
       //Fill nchg histo
       
      _p["nch_jetpt_F"]->fill(fpt,fnch);
      _p["nch_jetpt_C"]->fill(cpt,cnch);
      _p["nch_jetpt_B"]->fill(fpt,fnch);
      _p["nch_jetpt_B"]->fill(cpt,cnch);
            
      _h["nch_F"+to_str(f_bin)]->fill(fnch);
      _h["nch_C"+to_str(c_bin)]->fill(cnch);
      _h["nch_B"+to_str(f_bin)]->fill(fnch);
      _h["nch_B"+to_str(c_bin)]->fill(cnch);
    }
 
    
    void finalize() {
    
      // For r only
      /// @todo Replace with barchart()
      for (auto& hist : _hr) {
        for(size_t i=0; i < hist.second->numBins(); ++i) {
          double x = hist.second->bin(i).xMid();
          double bW = hist.second->bin(i).xWidth();
          hist.second->bin(i).scaleW(bW/(2.0*M_PI*x)); 
        }  
      }

      // The rest
      /// @todo Replace with barchart()
      for (auto& hist : _h) {                
        for (size_t i=0; i < hist.second->numBins(); ++i) {
          double bW = hist.second->bin(i).xWidth();
          hist.second->bin(i).scaleW(bW); 
        }
      }
                          
      for (size_t i_bin = 0; i_bin < 14; ++i_bin) {
      
        double sfB =  _h["nch_B"+to_str(i_bin)]->sumW();                         
        if (sfB) {
          scale(_h["zeta_B"+to_str(i_bin)],2.0/sfB);
          scale(_h["pTrel_B"+to_str(i_bin)],2.0/sfB);
          scale(_hr["r_B"+to_str(i_bin)],2.0/sfB);                      
          scale(_h["nch_B"+to_str(i_bin)], 2.0/sfB);
        }
        
        double sfF =  _h["nch_F"+to_str(i_bin)]->sumW();                         
        if (sfF) {
          scale(_h["zeta_F"+to_str(i_bin)],1.0/sfF);
          scale(_h["pTrel_F"+to_str(i_bin)],1.0/sfF);
          scale(_hr["r_F"+to_str(i_bin)],1.0/sfF);                      
          scale(_h["nch_F"+to_str(i_bin)], 1.0/sfF);
        }
        
        double sfC =  _h["nch_C"+to_str(i_bin)]->sumW();                         
        if (sfC) {
          scale(_h["zeta_C"+to_str(i_bin)],1.0/sfC);
          scale(_h["pTrel_C"+to_str(i_bin)],1.0/sfC);
          scale(_hr["r_C"+to_str(i_bin)],1.0/sfC);                      
          scale(_h["nch_C"+to_str(i_bin)], 1.0/sfC);                 
        }
      }
    }

  private:
    
    size_t GetJetBin(const double jetpt){
      size_t i_bin = 0;
      if (inRange(jetpt,100,200)) i_bin=0;
      if (inRange(jetpt,200,300)) i_bin=1;
      if (inRange(jetpt,300,400)) i_bin=2;
      if (inRange(jetpt,400,500)) i_bin=3;
      if (inRange(jetpt,500,600)) i_bin=4;
      if (inRange(jetpt,600,700)) i_bin=5;
      if (inRange(jetpt,700,800)) i_bin=6;
      if (inRange(jetpt,800,900)) i_bin=7;
      if (inRange(jetpt,900,1000)) i_bin=8;
      if (inRange(jetpt,1000,1200)) i_bin=9;
      if (inRange(jetpt,1200,1400)) i_bin=10;
      if (inRange(jetpt,1400,1600)) i_bin=11;
      if (inRange(jetpt,1600,2000)) i_bin=12;
      if (inRange(jetpt,2000,2500)) i_bin=13;   
      if(jetpt < 100) i_bin=0;
      if(jetpt > 2500) i_bin=13; 
      return i_bin;
    } 
  
    map<string, Histo1DPtr> _h;
    map<string, Histo1DPtr> _hr;
    map<string, Profile1DPtr> _p; 
  };
  
  // The hook for the plugin system
  RIVET_DECLARE_PLUGIN(ATLAS_2019_I1740909);
  
}