Rivet Analyses Reference

CMS_2011_S8978280

$K_s$, $\Lambda$, and Cascade$-$ transverse momentum and rapidity spectra at 900 and 7000 GeV.
Experiment: CMS (LHC)
Inspire ID: 890166
Status: VALIDATED
Authors:
  • Kevin Stenson
References:Beams: p+ p+
Beam energies: (450.0, 450.0); (3500.0, 3500.0) GeV
Run details:
  • Non-single-diffractive (NSD) events only. Should include double-diffractive (DD) events and non-diffractive (ND) events but NOT single-diffractive (SD) events. For example, in Pythia6 the SD processes to be turned off are 92 and 93, and in Pythia8 the SD processes are 103 and 104 (also called SoftQCD:singleDiffractive). Beam energy must be specified as analysis option "ENERGY" (900 or 7000) when rivet-merging samples.

The spectra of $K_S$, $\Lambda$, and Cascade- particles were measured versus transverse-momentum ($p_\perp$) and rapidity ($y$) in proton-proton collisions at center-of-mass energies 900 and 7000 GeV. The production is normalized to all non-single-diffractive (NSD) events using corrections for trigger and selection efficiency, acceptance, and branching ratios. The results cover a rapidity range of $|y| < 2$ and a $p_\perp$ range from 0 to 10 GeV ($K_S$ and $\Lambda$) and 0 to 6 GeV (Cascade-). Antiparticles are included in all measurements so only the sums of $\Lambda$ and $\bar{\Lambda}$, and Cascade$-$ and anti-Cascade$-$ are given. The rapidity distributions are shown versus $|y|$ but normalized to a unit of $y$. Ratios of $\Lambda$/$K_S$ and Cascade$-$/$\Lambda$ production versus $p_\perp$ and $|y|$ are also given, with somewhat smaller systematic uncertainties than obtained from taking the ratio of the individual distributions.' The data were corrected according to the SD/DD/ND content of the CMS trigger, as predicted by PYTHIA6. The uncertainties connected with correct or uncorrect modelling of diffraction were included in the systematic errors. Beam energy must be specified as analysis option "ENERGY" when rivet-merging samples.

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

namespace Rivet {


  /// @brief CMS strange particle spectra (Ks, Lambda, Cascade) in pp at 900 and 7000 GeV
  ///
  /// @author Kevin Stenson
  class CMS_2011_S8978280 : public Analysis {
  public:

    RIVET_DEFAULT_ANALYSIS_CTOR(CMS_2011_S8978280);


    void init() {
      UnstableParticles ufs(Cuts::absrap < 2);
      declare(ufs, "UFS");
      int beamEnergy = -1;
      if (isCompatibleWithSqrtS(900.))  beamEnergy = 1;
      else if (isCompatibleWithSqrtS(7000.))  beamEnergy = 2;
      else {
        MSG_WARNING("Could not decipher beam energy. For rivet-merge set -a CMS_2011_S8978280:energy=OPT, where OPT is 900 or 7000 (GeV is implied).");
      }
      
      // Particle distributions versus rapidity and transverse momentum
      if (beamEnergy == 1){
        book(_h_dNKshort_dy  ,1, 1, 1);
        book(_h_dNKshort_dpT ,2, 1, 1);
        book(_h_dNLambda_dy  ,3, 1, 1);
        book(_h_dNLambda_dpT ,4, 1, 1);
        book(_h_dNXi_dy      ,5, 1, 1);
        book(_h_dNXi_dpT     ,6, 1, 1);
        //
        book(_h_LampT_KpT , 7, 1, 1);
        book(_h_XipT_LampT, 8, 1, 1);
        book(_h_Lamy_Ky   , 9, 1, 1);
        book(_h_Xiy_Lamy  , 10, 1, 1);

      } else if (beamEnergy == 2){
        book(_h_dNKshort_dy  ,1, 1, 2);
        book(_h_dNKshort_dpT ,2, 1, 2);
        book(_h_dNLambda_dy  ,3, 1, 2);
        book(_h_dNLambda_dpT ,4, 1, 2);
        book(_h_dNXi_dy      ,5, 1, 2);
        book(_h_dNXi_dpT     ,6, 1, 2);
        //
        book(_h_LampT_KpT , 7, 1, 2);
        book(_h_XipT_LampT, 8, 1, 2);
        book(_h_Lamy_Ky   , 9, 1, 2);
        book(_h_Xiy_Lamy  , 10, 1, 2);
      } else {
        MSG_WARNING("Could not initialize properly.");
      }
    }


    void analyze(const Event& event) {

      const UnstableParticles& parts = apply<UnstableParticles>(event, "UFS");
      for (const Particle& p : parts.particles()) {
        switch (p.abspid()) {
        case PID::K0S:
          _h_dNKshort_dy->fill(p.absrap());
          _h_dNKshort_dpT->fill(p.pT()/GeV);
          break;

        case PID::LAMBDA:
          // Lambda should not have Cascade or Omega ancestors since they should not decay. But just in case...
          if ( !( p.hasAncestor(3322) || p.hasAncestor(-3322) || p.hasAncestor(3312) || p.hasAncestor(-3312) || p.hasAncestor(3334) || p.hasAncestor(-3334) ) ) {
            _h_dNLambda_dy->fill(p.absrap());
            _h_dNLambda_dpT->fill(p.pT()/GeV);
          }
          break;

        case PID::XIMINUS:
          // Cascade should not have Omega ancestors since it should not decay.  But just in case...
          if ( !( p.hasAncestor(3334) || p.hasAncestor(-3334) ) ) {
            _h_dNXi_dy->fill(p.absrap());
            _h_dNXi_dpT->fill(p.pT()/GeV);
          }
          break;
        }

      }
    }


    void finalize() {
      divide(_h_dNLambda_dpT,_h_dNKshort_dpT, _h_LampT_KpT);
      divide(_h_dNXi_dpT,_h_dNLambda_dpT, _h_XipT_LampT);
      divide(_h_dNLambda_dy,_h_dNKshort_dy, _h_Lamy_Ky);
      divide(_h_dNXi_dy,_h_dNLambda_dy, _h_Xiy_Lamy);
      const double normpT = 1.0/sumOfWeights();
      const double normy = 0.5*normpT; // Accounts for using |y| instead of y
      scale(_h_dNKshort_dy, normy);
      scale(_h_dNKshort_dpT, normpT);
      scale(_h_dNLambda_dy, normy);
      scale(_h_dNLambda_dpT, normpT);
      scale(_h_dNXi_dy, normy);
      scale(_h_dNXi_dpT, normpT);
    }


  private:

    /// @name Particle distributions versus rapidity and transverse momentum
    /// @{
    Histo1DPtr _h_dNKshort_dy, _h_dNKshort_dpT, _h_dNLambda_dy, _h_dNLambda_dpT, _h_dNXi_dy, _h_dNXi_dpT;
    Scatter2DPtr _h_LampT_KpT, _h_XipT_LampT, _h_Lamy_Ky, _h_Xiy_Lamy;
    /// @}

  };



  RIVET_DECLARE_ALIASED_PLUGIN(CMS_2011_S8978280, CMS_2011_I890166);

}