CMS_2011_S8978280.cc
Go to the documentation of this file.
00001 // -*- C++ -*- 00002 #include "Rivet/Analysis.hh" 00003 #include "Rivet/Projections/UnstableFinalState.hh" 00004 00005 namespace Rivet { 00006 00007 00008 /// @brief CMS strange particle spectra (Ks, Lambda, Cascade) in pp at 900 and 7000 GeV 00009 /// @author Kevin Stenson 00010 class CMS_2011_S8978280 : public Analysis { 00011 public: 00012 00013 /// Constructor 00014 CMS_2011_S8978280() : Analysis("CMS_2011_S8978280") {} 00015 00016 00017 void init() { 00018 // Need wide range of eta because cut on rapidity not pseudorapidity 00019 UnstableFinalState ufs(-8.0, 8.0, 0.0*GeV); 00020 addProjection(ufs, "UFS"); 00021 00022 // Particle distributions versus rapidity and transverse momentum 00023 if (fuzzyEquals(sqrtS()/GeV, 900.0*GeV)){ 00024 _h_dNKshort_dy = bookHisto1D(1, 1, 1); 00025 _h_dNKshort_dpT = bookHisto1D(2, 1, 1); 00026 _h_dNLambda_dy = bookHisto1D(3, 1, 1); 00027 _h_dNLambda_dpT = bookHisto1D(4, 1, 1); 00028 _h_dNXi_dy = bookHisto1D(5, 1, 1); 00029 _h_dNXi_dpT = bookHisto1D(6, 1, 1); 00030 00031 _h_LampT_KpT = bookScatter2D(7, 1, 1); 00032 _h_XipT_LampT = bookScatter2D(8, 1, 1); 00033 _h_Lamy_Ky = bookScatter2D(9, 1, 1); 00034 _h_Xiy_Lamy = bookScatter2D(10, 1, 1); 00035 } else if (fuzzyEquals(sqrtS()/GeV, 7000.0*GeV)){ 00036 _h_dNKshort_dy = bookHisto1D(1, 1, 2); 00037 _h_dNKshort_dpT = bookHisto1D(2, 1, 2); 00038 _h_dNLambda_dy = bookHisto1D(3, 1, 2); 00039 _h_dNLambda_dpT = bookHisto1D(4, 1, 2); 00040 _h_dNXi_dy = bookHisto1D(5, 1, 2); 00041 _h_dNXi_dpT = bookHisto1D(6, 1, 2); 00042 00043 _h_LampT_KpT = bookScatter2D(7, 1, 2); 00044 _h_XipT_LampT = bookScatter2D(8, 1, 2); 00045 _h_Lamy_Ky = bookScatter2D(9, 1, 2); 00046 _h_Xiy_Lamy = bookScatter2D(10, 1, 2); 00047 } 00048 } 00049 00050 00051 void analyze(const Event& event) { 00052 const double weight = event.weight(); 00053 00054 const UnstableFinalState& parts = applyProjection<UnstableFinalState>(event, "UFS"); 00055 00056 foreach (const Particle& p, parts.particles()) { 00057 const double pT = p.pT(); 00058 const double y = p.absrap(); 00059 const PdgId pid = p.abspid(); 00060 00061 if (y < 2.0) { 00062 switch (pid) { 00063 case PID::K0S: 00064 _h_dNKshort_dy->fill(y, weight); 00065 _h_dNKshort_dpT->fill(pT, weight); 00066 break; 00067 case PID::LAMBDA: 00068 // Lambda should not have Cascade or Omega ancestors since they should not decay. But just in case... 00069 if ( !( p.hasAncestor(3322) || p.hasAncestor(-3322) || p.hasAncestor(3312) || p.hasAncestor(-3312) || p.hasAncestor(3334) || p.hasAncestor(-3334) ) ) { 00070 _h_dNLambda_dy->fill(y, weight); 00071 _h_dNLambda_dpT->fill(pT, weight); 00072 } 00073 break; 00074 case PID::XIMINUS: 00075 // Cascade should not have Omega ancestors since it should not decay. But just in case... 00076 if ( !( p.hasAncestor(3334) || p.hasAncestor(-3334) ) ) { 00077 _h_dNXi_dy->fill(y, weight); 00078 _h_dNXi_dpT->fill(pT, weight); 00079 } 00080 break; 00081 } 00082 } 00083 } 00084 } 00085 00086 00087 void finalize() { 00088 00089 divide(_h_dNLambda_dpT,_h_dNKshort_dpT, 00090 _h_LampT_KpT); 00091 00092 divide(_h_dNXi_dpT,_h_dNLambda_dpT, 00093 _h_XipT_LampT); 00094 00095 divide(_h_dNLambda_dy,_h_dNKshort_dy, 00096 _h_Lamy_Ky); 00097 00098 divide(_h_dNXi_dy,_h_dNLambda_dy, 00099 _h_Xiy_Lamy); 00100 00101 double normpT = 1.0/sumOfWeights(); 00102 double normy = 0.5*normpT; // Accounts for using |y| instead of y 00103 scale(_h_dNKshort_dy, normy); 00104 scale(_h_dNKshort_dpT, normpT); 00105 scale(_h_dNLambda_dy, normy); 00106 scale(_h_dNLambda_dpT, normpT); 00107 scale(_h_dNXi_dy, normy); 00108 scale(_h_dNXi_dpT, normpT); 00109 } 00110 00111 00112 private: 00113 00114 // Particle distributions versus rapidity and transverse momentum 00115 Histo1DPtr _h_dNKshort_dy; 00116 Histo1DPtr _h_dNKshort_dpT; 00117 Histo1DPtr _h_dNLambda_dy; 00118 Histo1DPtr _h_dNLambda_dpT; 00119 Histo1DPtr _h_dNXi_dy; 00120 Histo1DPtr _h_dNXi_dpT; 00121 00122 Scatter2DPtr _h_LampT_KpT; 00123 Scatter2DPtr _h_XipT_LampT; 00124 Scatter2DPtr _h_Lamy_Ky; 00125 Scatter2DPtr _h_Xiy_Lamy; 00126 00127 }; 00128 00129 00130 00131 // The hook for the plugin system 00132 DECLARE_RIVET_PLUGIN(CMS_2011_S8978280); 00133 00134 } Generated on Tue Sep 30 2014 19:45:43 for The Rivet MC analysis system by ![]() |