Rivet analyses referenceCMS_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:
Beam energies: (450.0, 450.0); (3500.0, 3500.0) GeV Run details:
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// -*- C++ -*-
2#include "Rivet/Analysis.hh"
3#include "Rivet/Projections/UnstableParticles.hh"
4
5namespace Rivet {
6
7
8 /// @brief CMS strange particle spectra (Ks, Lambda, Cascade) in pp at 900 and 7000 GeV
9 ///
10 /// @author Kevin Stenson
11 class CMS_2011_S8978280 : public Analysis {
12 public:
13
14 RIVET_DEFAULT_ANALYSIS_CTOR(CMS_2011_S8978280);
15
16
17 void init() {
18 UnstableParticles ufs(Cuts::absrap < 2);
19 declare(ufs, "UFS");
20 int beamEnergy = -1;
21 if (isCompatibleWithSqrtS(900.)) beamEnergy = 1;
22 else if (isCompatibleWithSqrtS(7000.)) beamEnergy = 2;
23 else {
24 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).");
25 }
26
27 // Particle distributions versus rapidity and transverse momentum
28 if (beamEnergy == 1){
29 book(_h_dNKshort_dy ,1, 1, 1);
30 book(_h_dNKshort_dpT ,2, 1, 1);
31 book(_h_dNLambda_dy ,3, 1, 1);
32 book(_h_dNLambda_dpT ,4, 1, 1);
33 book(_h_dNXi_dy ,5, 1, 1);
34 book(_h_dNXi_dpT ,6, 1, 1);
35 //
36 book(_h_LampT_KpT , 7, 1, 1);
37 book(_h_XipT_LampT, 8, 1, 1);
38 book(_h_Lamy_Ky , 9, 1, 1);
39 book(_h_Xiy_Lamy , 10, 1, 1);
40
41 } else if (beamEnergy == 2){
42 book(_h_dNKshort_dy ,1, 1, 2);
43 book(_h_dNKshort_dpT ,2, 1, 2);
44 book(_h_dNLambda_dy ,3, 1, 2);
45 book(_h_dNLambda_dpT ,4, 1, 2);
46 book(_h_dNXi_dy ,5, 1, 2);
47 book(_h_dNXi_dpT ,6, 1, 2);
48 //
49 book(_h_LampT_KpT , 7, 1, 2);
50 book(_h_XipT_LampT, 8, 1, 2);
51 book(_h_Lamy_Ky , 9, 1, 2);
52 book(_h_Xiy_Lamy , 10, 1, 2);
53 } else {
54 MSG_WARNING("Could not initialize properly.");
55 }
56 }
57
58
59 void analyze(const Event& event) {
60
61 const UnstableParticles& parts = apply<UnstableParticles>(event, "UFS");
62 for (const Particle& p : parts.particles()) {
63 switch (p.abspid()) {
64 case PID::K0S:
65 _h_dNKshort_dy->fill(p.absrap());
66 _h_dNKshort_dpT->fill(p.pT()/GeV);
67 break;
68
69 case PID::LAMBDA:
70 // Lambda should not have Cascade or Omega ancestors since they should not decay. But just in case...
71 if ( !( p.hasAncestor(3322) || p.hasAncestor(-3322) || p.hasAncestor(3312) || p.hasAncestor(-3312) || p.hasAncestor(3334) || p.hasAncestor(-3334) ) ) {
72 _h_dNLambda_dy->fill(p.absrap());
73 _h_dNLambda_dpT->fill(p.pT()/GeV);
74 }
75 break;
76
77 case PID::XIMINUS:
78 // Cascade should not have Omega ancestors since it should not decay. But just in case...
79 if ( !( p.hasAncestor(3334) || p.hasAncestor(-3334) ) ) {
80 _h_dNXi_dy->fill(p.absrap());
81 _h_dNXi_dpT->fill(p.pT()/GeV);
82 }
83 break;
84 }
85
86 }
87 }
88
89
90 void finalize() {
91 divide(_h_dNLambda_dpT,_h_dNKshort_dpT, _h_LampT_KpT);
92 divide(_h_dNXi_dpT,_h_dNLambda_dpT, _h_XipT_LampT);
93 divide(_h_dNLambda_dy,_h_dNKshort_dy, _h_Lamy_Ky);
94 divide(_h_dNXi_dy,_h_dNLambda_dy, _h_Xiy_Lamy);
95 const double normpT = 1.0/sumOfWeights();
96 const double normy = 0.5*normpT; // Accounts for using |y| instead of y
97 scale(_h_dNKshort_dy, normy);
98 scale(_h_dNKshort_dpT, normpT);
99 scale(_h_dNLambda_dy, normy);
100 scale(_h_dNLambda_dpT, normpT);
101 scale(_h_dNXi_dy, normy);
102 scale(_h_dNXi_dpT, normpT);
103 }
104
105
106 private:
107
108 /// @name Particle distributions versus rapidity and transverse momentum
109 /// @{
110 Histo1DPtr _h_dNKshort_dy, _h_dNKshort_dpT, _h_dNLambda_dy, _h_dNLambda_dpT, _h_dNXi_dy, _h_dNXi_dpT;
111 Scatter2DPtr _h_LampT_KpT, _h_XipT_LampT, _h_Lamy_Ky, _h_Xiy_Lamy;
112 /// @}
113
114 };
115
116
117
118 RIVET_DECLARE_ALIASED_PLUGIN(CMS_2011_S8978280, CMS_2011_I890166);
119
120}
|