Rivet Analyses Reference

SLD_2004_S5693039

Production of $\pi^+$, $\pi^-$, $K^+$, $K^-$, $p$ and $\bar p$ in Light ($uds$), $c$ and $b$ Jets from Z Decays
Experiment: SLD (SLC)
Inspire ID: 630327
Status: VALIDATED
Authors:
  • Peter Richardson
References:Beams: e+ e-
Beam energies: (45.6, 45.6) GeV
Run details:
  • Hadronic Z decay events generated on the Z pole ($\sqrt{s} = 91.2$ GeV)

Measurements of the differential production rates of stable charged particles in hadronic $Z^0$ decays, and of charged pions, kaons and protons identified over a wide momentum range. In addition to flavour-inclusive $Z^0$ decays, measurements are made for $Z^0$ decays into light ($u$, $d$, $s$), $c$ and $b$ primary flavors.

Source code: SLD_2004_S5693039.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
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
// -*- C++ -*-
#include "Rivet/Analysis.hh"
#include "Rivet/Projections/Beam.hh"
#include "Rivet/Projections/FinalState.hh"
#include "Rivet/Projections/ChargedFinalState.hh"
#include "Rivet/Projections/Thrust.hh"

#define I_KNOW_THE_INITIAL_QUARKS_PROJECTION_IS_DODGY_BUT_NEED_TO_USE_IT
#include "Rivet/Projections/InitialQuarks.hh"

namespace Rivet {


  /// @brief SLD flavour-dependent fragmentation paper
  ///
  /// @author Peter Richardson
  class SLD_2004_S5693039 : public Analysis {
  public:

    /// Constructor
    RIVET_DEFAULT_ANALYSIS_CTOR(SLD_2004_S5693039);


    /// @name Analysis methods
    //@{

    void analyze(const Event& e) {
      // First, veto on leptonic events by requiring at least 2 charged FS particles
      const FinalState& fs = apply<FinalState>(e, "FS");
      const size_t numParticles = fs.particles().size();

      // Even if we only generate hadronic events, we still need a cut on numCharged >= 2.
      if (numParticles < 2) {
        MSG_DEBUG("Failed ncharged cut");
        vetoEvent;
      }
      MSG_DEBUG("Passed ncharged cut");

      // Get beams and average beam momentum
      const ParticlePair& beams = apply<Beam>(e, "Beams").beams();
      const double meanBeamMom = ( beams.first.p3().mod() +
                                   beams.second.p3().mod() ) / 2.0;
      MSG_DEBUG("Avg beam momentum = " << meanBeamMom);
      int flavour = 0;
      const InitialQuarks& iqf = apply<InitialQuarks>(e, "IQF");

      // If we only have two quarks (qqbar), just take the flavour.
      // If we have more than two quarks, look for the highest energetic q-qbar pair.
      Particles quarks;
      if (iqf.particles().size() == 2) {
        flavour = iqf.particles().front().abspid();
        quarks = iqf.particles();
      }
      else {
        map<int, Particle > quarkmap;
        for (const Particle& p : iqf.particles()) {
          if (quarkmap.find(p.pid())==quarkmap.end())
            quarkmap[p.pid()] = p;
          else if (quarkmap[p.pid()].E() < p.E())
            quarkmap[p.pid()] = p;
        }
        double maxenergy = 0.;
        for (int i = 1; i <= 5; ++i) {
          double energy(0.);
          if(quarkmap.find( i)!=quarkmap.end())
            energy += quarkmap[ i].E();
          if(quarkmap.find(-i)!=quarkmap.end())
            energy += quarkmap[-i].E();
          if (energy > maxenergy) flavour = i;
        }
        if(quarkmap.find( flavour)!=quarkmap.end())
          quarks.push_back(quarkmap[ flavour]);
        if(quarkmap.find(-flavour)!=quarkmap.end())
          quarks.push_back(quarkmap[-flavour]);
      }

      // total multiplicities
      switch (flavour) {
      case PID::DQUARK:
      case PID::UQUARK:
      case PID::SQUARK:
        _weightLight  ->fill();
        _weightedTotalChargedPartNumLight  ->fill(numParticles);
        break;
      case PID::CQUARK:
        _weightCharm  ->fill();
        _weightedTotalChargedPartNumCharm  ->fill(numParticles);
        break;
      case PID::BQUARK:
        _weightBottom ->fill();
        _weightedTotalChargedPartNumBottom ->fill(numParticles);
        break;
      }
      // thrust axis for projections
      Vector3 axis = apply<Thrust>(e, "Thrust").thrustAxis();
      double dot(0.);
      if(!quarks.empty()) {
        dot = quarks[0].p3().dot(axis);
        if(quarks[0].pid()<0) dot *= -1.;
      }
      // spectra and individual multiplicities
      for (const Particle& p : fs.particles()) {
        double pcm = p.p3().mod();
        const double xp = pcm/meanBeamMom;

        // if in quark or antiquark hemisphere
        bool quark = p.p3().dot(axis)*dot>0.;

        _h_PCharged ->fill(pcm     );
        // all charged
        switch (flavour) {
        case PID::DQUARK:
        case PID::UQUARK:
        case PID::SQUARK:
          _h_XpChargedL->fill(xp);
          break;
        case PID::CQUARK:
          _h_XpChargedC->fill(xp);
          break;
        case PID::BQUARK:
          _h_XpChargedB->fill(xp);
          break;
        }

        int id = p.abspid();
        // charged pions
        if (id == PID::PIPLUS) {
          _h_XpPiPlus->fill(xp);
          _h_XpPiPlusTotal->fill(xp);
          switch (flavour) {
          case PID::DQUARK:
          case PID::UQUARK:
          case PID::SQUARK:
            _h_XpPiPlusL->fill(xp);
            _h_NPiPlusL->fill(sqrtS());
            if( ( quark && p.pid()>0 ) || ( !quark && p.pid()<0 ))
              _h_RPiPlus->fill(xp);
            else
              _h_RPiMinus->fill(xp);
            break;
          case PID::CQUARK:
            _h_XpPiPlusC->fill(xp);
            _h_NPiPlusC->fill(sqrtS());
            break;
          case PID::BQUARK:
            _h_XpPiPlusB->fill(xp);
            _h_NPiPlusB->fill(sqrtS());
            break;
          }
        }
        else if (id == PID::KPLUS) {
          _h_XpKPlus->fill(xp);
          _h_XpKPlusTotal->fill(xp);
          switch (flavour) {
          case PID::DQUARK:
          case PID::UQUARK:
          case PID::SQUARK:
            _h_XpKPlusL->fill(xp);
            _h_NKPlusL->fill(sqrtS());
            if( ( quark && p.pid()>0 ) || ( !quark && p.pid()<0 ))
              _h_RKPlus->fill(xp);
            else
              _h_RKMinus->fill(xp);
            break;
          case PID::CQUARK:
            _h_XpKPlusC->fill(xp);
            _h_NKPlusC->fill(sqrtS());
            break;
          case PID::BQUARK:
            _h_XpKPlusB->fill(xp);
            _h_NKPlusB->fill(sqrtS());
            break;
          }
        }
        else if (id == PID::PROTON) {
          _h_XpProton->fill(xp);
          _h_XpProtonTotal->fill(xp);
          switch (flavour) {
          case PID::DQUARK:
          case PID::UQUARK:
          case PID::SQUARK:
            _h_XpProtonL->fill(xp);
            _h_NProtonL->fill(sqrtS());
            if( ( quark && p.pid()>0 ) || ( !quark && p.pid()<0 ))
              _h_RProton->fill(xp);
            else
              _h_RPBar  ->fill(xp);
            break;
          case PID::CQUARK:
            _h_XpProtonC->fill(xp);
            _h_NProtonC->fill(sqrtS());
            break;
          case PID::BQUARK:
            _h_XpProtonB->fill(xp);
            _h_NProtonB->fill(sqrtS());
            break;
          }
        }
      }
    }


    void init() {
      // Projections
      declare(Beam(), "Beams");
      declare(ChargedFinalState(), "FS");
      declare(InitialQuarks(), "IQF");
      declare(Thrust(FinalState()), "Thrust");

      // Book histograms
      book(_h_PCharged   , 1, 1, 1);
      book(_h_XpPiPlus   , 2, 1, 2);
      book(_h_XpKPlus    , 3, 1, 2);
      book(_h_XpProton   , 4, 1, 2);
      book(_h_XpPiPlusTotal , 2, 2, 2);
      book(_h_XpKPlusTotal  , 3, 2, 2);
      book(_h_XpProtonTotal , 4, 2, 2);
      book(_h_XpPiPlusL  , 5, 1, 1);
      book(_h_XpPiPlusC  , 5, 1, 2);
      book(_h_XpPiPlusB  , 5, 1, 3);
      book(_h_XpKPlusL   , 6, 1, 1);
      book(_h_XpKPlusC   , 6, 1, 2);
      book(_h_XpKPlusB   , 6, 1, 3);
      book(_h_XpProtonL  , 7, 1, 1);
      book(_h_XpProtonC  , 7, 1, 2);
      book(_h_XpProtonB  , 7, 1, 3);
      book(_h_XpChargedL , 8, 1, 1);
      book(_h_XpChargedC , 8, 1, 2);
      book(_h_XpChargedB , 8, 1, 3);

      book(_h_NPiPlusL  , 5, 2, 1);
      book(_h_NPiPlusC  , 5, 2, 2);
      book(_h_NPiPlusB  , 5, 2, 3);
      book(_h_NKPlusL   , 6, 2, 1);
      book(_h_NKPlusC   , 6, 2, 2);
      book(_h_NKPlusB   , 6, 2, 3);
      book(_h_NProtonL  , 7, 2, 1);
      book(_h_NProtonC  , 7, 2, 2);
      book(_h_NProtonB  , 7, 2, 3);

      book(_h_RPiPlus  , 9, 1, 1);
      book(_h_RPiMinus , 9, 1, 2);
      book(_h_RKPlus   ,10, 1, 1);
      book(_h_RKMinus  ,10, 1, 2);
      book(_h_RProton  ,11, 1, 1);
      book(_h_RPBar    ,11, 1, 2);

      // Ratios: used as target of divide() later
      book(_s_PiM_PiP,  9, 1, 3);
      book(_s_KM_KP  , 10, 1, 3);
      book(_s_Pr_PBar, 11, 1, 3);

      book(_weightedTotalChargedPartNumLight, "_weightedTotalChargedPartNumLight");
      book(_weightedTotalChargedPartNumCharm, "_weightedTotalChargedPartNumCharm");
      book(_weightedTotalChargedPartNumBottom, "_weightedTotalChargedPartNumBottom");
      book(_weightLight, "_weightLight");
      book(_weightCharm, "_weightCharm");
      book(_weightBottom, "_weightBottom");

      book(tmp1, 8, 2, 1, true);
      book(tmp2, 8, 2, 2, true);
      book(tmp3, 8, 2, 3, true);
      book(tmp4, 8, 3, 2, true);
      book(tmp5, 8, 3, 3, true);


    }


    /// Finalize
    void finalize() {

      // Multiplicities
      /// @todo Include errors
      const double avgNumPartsLight = _weightedTotalChargedPartNumLight->val() / _weightLight->val();
      const double avgNumPartsCharm = _weightedTotalChargedPartNumCharm->val() / _weightCharm->val();
      const double avgNumPartsBottom = _weightedTotalChargedPartNumBottom->val() / _weightBottom->val();
      tmp1->point(0).setY(avgNumPartsLight);
      tmp2->point(0).setY(avgNumPartsCharm);
      tmp3->point(0).setY(avgNumPartsBottom);
      tmp4->point(0).setY(avgNumPartsCharm - avgNumPartsLight);
      tmp5->point(0).setY(avgNumPartsBottom - avgNumPartsLight);

      // Do divisions
      divide(*_h_RPiMinus - *_h_RPiPlus, *_h_RPiMinus + *_h_RPiPlus, _s_PiM_PiP);
      divide(*_h_RKMinus - *_h_RKPlus, *_h_RKMinus + *_h_RKPlus, _s_KM_KP);
      divide(*_h_RProton - *_h_RPBar, *_h_RProton + *_h_RPBar, _s_Pr_PBar);

      // Scale histograms
      scale(_h_PCharged,      1./sumOfWeights());
      scale(_h_XpPiPlus,      1./sumOfWeights());
      scale(_h_XpKPlus,       1./sumOfWeights());
      scale(_h_XpProton,      1./sumOfWeights());
      scale(_h_XpPiPlusTotal, 1./sumOfWeights());
      scale(_h_XpKPlusTotal,  1./sumOfWeights());
      scale(_h_XpProtonTotal, 1./sumOfWeights());
      scale(_h_XpPiPlusL,     1. / *_weightLight);
      scale(_h_XpPiPlusC,     1. / *_weightCharm);
      scale(_h_XpPiPlusB,     1. / *_weightBottom);
      scale(_h_XpKPlusL,      1. / *_weightLight);
      scale(_h_XpKPlusC,      1. / *_weightCharm);
      scale(_h_XpKPlusB,      1. / *_weightBottom);
      scale(_h_XpProtonL,     1. / *_weightLight);
      scale(_h_XpProtonC,     1. / *_weightCharm);
      scale(_h_XpProtonB,     1. / *_weightBottom);

      scale(_h_XpChargedL, 1. / *_weightLight);
      scale(_h_XpChargedC, 1. / *_weightCharm);
      scale(_h_XpChargedB, 1. / *_weightBottom);

      scale(_h_NPiPlusL, 1. / *_weightLight);
      scale(_h_NPiPlusC, 1. / *_weightCharm);
      scale(_h_NPiPlusB, 1. / *_weightBottom);
      scale(_h_NKPlusL,  1. / *_weightLight);
      scale(_h_NKPlusC,  1. / *_weightCharm);
      scale(_h_NKPlusB,  1. / *_weightBottom);
      scale(_h_NProtonL, 1. / *_weightLight);
      scale(_h_NProtonC, 1. / *_weightCharm);
      scale(_h_NProtonB, 1. / *_weightBottom);

      // Paper suggests this should be 0.5/weight but it has to be 1.0 to get normalisations right...
      scale(_h_RPiPlus,  1. / *_weightLight);
      scale(_h_RPiMinus, 1. / *_weightLight);
      scale(_h_RKPlus,   1. / *_weightLight);
      scale(_h_RKMinus,  1. / *_weightLight);
      scale(_h_RProton,  1. / *_weightLight);
      scale(_h_RPBar,    1. / *_weightLight);

      // convert ratio to %
      _s_PiM_PiP->scaleY(100.);
      _s_KM_KP  ->scaleY(100.);
      _s_Pr_PBar->scaleY(100.);
    }

    //@}


  private:

    Scatter2DPtr tmp1, tmp2, tmp3, tmp4, tmp5;

    /// Multiplicities
    CounterPtr _weightedTotalChargedPartNumLight, _weightedTotalChargedPartNumCharm, _weightedTotalChargedPartNumBottom;

    /// Weights
    CounterPtr _weightLight, _weightCharm, _weightBottom;

    // Histograms
    /// @{
    Histo1DPtr _h_PCharged;
    Histo1DPtr _h_XpPiPlus, _h_XpKPlus, _h_XpProton;
    Histo1DPtr _h_XpPiPlusTotal, _h_XpKPlusTotal, _h_XpProtonTotal;
    Histo1DPtr _h_XpPiPlusL, _h_XpPiPlusC, _h_XpPiPlusB;
    Histo1DPtr _h_XpKPlusL, _h_XpKPlusC, _h_XpKPlusB;
    Histo1DPtr _h_XpProtonL, _h_XpProtonC, _h_XpProtonB;
    Histo1DPtr _h_XpChargedL, _h_XpChargedC, _h_XpChargedB;
    Histo1DPtr _h_NPiPlusL, _h_NPiPlusC, _h_NPiPlusB;
    Histo1DPtr _h_NKPlusL, _h_NKPlusC, _h_NKPlusB;
    Histo1DPtr _h_NProtonL, _h_NProtonC, _h_NProtonB;
    Histo1DPtr _h_RPiPlus, _h_RPiMinus, _h_RKPlus;
    Histo1DPtr _h_RKMinus, _h_RProton, _h_RPBar;
    Scatter2DPtr _s_PiM_PiP, _s_KM_KP, _s_Pr_PBar;
    /// @}

  };



  RIVET_DECLARE_ALIASED_PLUGIN(SLD_2004_S5693039, SLD_2004_I630327);

}