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() 00015 : Analysis("CMS_2011_S8978280") 00016 { } 00017 00018 00019 void init() { 00020 UnstableFinalState ufs(Cuts::absrap < 2); 00021 addProjection(ufs, "UFS"); 00022 00023 // Particle distributions versus rapidity and transverse momentum 00024 if (fuzzyEquals(sqrtS()/GeV, 900*GeV)){ 00025 _h_dNKshort_dy = bookHisto1D(1, 1, 1); 00026 _h_dNKshort_dpT = bookHisto1D(2, 1, 1); 00027 _h_dNLambda_dy = bookHisto1D(3, 1, 1); 00028 _h_dNLambda_dpT = bookHisto1D(4, 1, 1); 00029 _h_dNXi_dy = bookHisto1D(5, 1, 1); 00030 _h_dNXi_dpT = bookHisto1D(6, 1, 1); 00031 // 00032 _h_LampT_KpT = bookScatter2D(7, 1, 1); 00033 _h_XipT_LampT = bookScatter2D(8, 1, 1); 00034 _h_Lamy_Ky = bookScatter2D(9, 1, 1); 00035 _h_Xiy_Lamy = bookScatter2D(10, 1, 1); 00036 00037 } else if (fuzzyEquals(sqrtS()/GeV, 7000*GeV)){ 00038 _h_dNKshort_dy = bookHisto1D(1, 1, 2); 00039 _h_dNKshort_dpT = bookHisto1D(2, 1, 2); 00040 _h_dNLambda_dy = bookHisto1D(3, 1, 2); 00041 _h_dNLambda_dpT = bookHisto1D(4, 1, 2); 00042 _h_dNXi_dy = bookHisto1D(5, 1, 2); 00043 _h_dNXi_dpT = bookHisto1D(6, 1, 2); 00044 // 00045 _h_LampT_KpT = bookScatter2D(7, 1, 2); 00046 _h_XipT_LampT = bookScatter2D(8, 1, 2); 00047 _h_Lamy_Ky = bookScatter2D(9, 1, 2); 00048 _h_Xiy_Lamy = bookScatter2D(10, 1, 2); 00049 } 00050 } 00051 00052 00053 void analyze(const Event& event) { 00054 const double weight = event.weight(); 00055 00056 const UnstableFinalState& parts = applyProjection<UnstableFinalState>(event, "UFS"); 00057 foreach (const Particle& p, parts.particles()) { 00058 switch (p.abspid()) { 00059 case PID::K0S: 00060 _h_dNKshort_dy->fill(p.absrap(), weight); 00061 _h_dNKshort_dpT->fill(p.pT(), weight); 00062 break; 00063 00064 case PID::LAMBDA: 00065 // Lambda should not have Cascade or Omega ancestors since they should not decay. But just in case... 00066 if ( !( p.hasAncestor(3322) || p.hasAncestor(-3322) || p.hasAncestor(3312) || p.hasAncestor(-3312) || p.hasAncestor(3334) || p.hasAncestor(-3334) ) ) { 00067 _h_dNLambda_dy->fill(p.absrap(), weight); 00068 _h_dNLambda_dpT->fill(p.pT(), weight); 00069 } 00070 break; 00071 00072 case PID::XIMINUS: 00073 // Cascade should not have Omega ancestors since it should not decay. But just in case... 00074 if ( !( p.hasAncestor(3334) || p.hasAncestor(-3334) ) ) { 00075 _h_dNXi_dy->fill(p.absrap(), weight); 00076 _h_dNXi_dpT->fill(p.pT(), weight); 00077 } 00078 break; 00079 } 00080 00081 } 00082 } 00083 00084 00085 void finalize() { 00086 divide(_h_dNLambda_dpT,_h_dNKshort_dpT, _h_LampT_KpT); 00087 divide(_h_dNXi_dpT,_h_dNLambda_dpT, _h_XipT_LampT); 00088 divide(_h_dNLambda_dy,_h_dNKshort_dy, _h_Lamy_Ky); 00089 divide(_h_dNXi_dy,_h_dNLambda_dy, _h_Xiy_Lamy); 00090 const double normpT = 1.0/sumOfWeights(); 00091 const double normy = 0.5*normpT; // Accounts for using |y| instead of y 00092 scale(_h_dNKshort_dy, normy); 00093 scale(_h_dNKshort_dpT, normpT); 00094 scale(_h_dNLambda_dy, normy); 00095 scale(_h_dNLambda_dpT, normpT); 00096 scale(_h_dNXi_dy, normy); 00097 scale(_h_dNXi_dpT, normpT); 00098 } 00099 00100 00101 private: 00102 00103 // Particle distributions versus rapidity and transverse momentum 00104 Histo1DPtr _h_dNKshort_dy, _h_dNKshort_dpT, _h_dNLambda_dy, _h_dNLambda_dpT, _h_dNXi_dy, _h_dNXi_dpT; 00105 Scatter2DPtr _h_LampT_KpT, _h_XipT_LampT, _h_Lamy_Ky, _h_Xiy_Lamy; 00106 00107 }; 00108 00109 00110 00111 // The hook for the plugin system 00112 DECLARE_RIVET_PLUGIN(CMS_2011_S8978280); 00113 00114 } Generated on Thu Mar 10 2016 08:29:49 for The Rivet MC analysis system by ![]() |