rivet is hosted by Hepforge, IPPP Durham
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       _h_dNKshort_dy  = bookHisto1D(1, 1, 1);
00024       _h_dNKshort_dpT = bookHisto1D(2, 1, 1);
00025       _h_dNLambda_dy  = bookHisto1D(3, 1, 1);
00026       _h_dNLambda_dpT = bookHisto1D(4, 1, 1);
00027       _h_dNXi_dy      = bookHisto1D(5, 1, 1);
00028       _h_dNXi_dpT     = bookHisto1D(6, 1, 1);
00029 
00030       _h_LampT_KpT  = bookScatter2D(7, 1, 1);
00031       _h_XipT_LampT     = bookScatter2D(8, 1, 1);
00032       _h_Lamy_Ky    = bookScatter2D(9, 1, 1);
00033       _h_Xiy_Lamy   = bookScatter2D(10, 1, 1);
00034 
00035     }
00036 
00037 
00038     void analyze(const Event& event) {
00039       const double weight = event.weight();
00040 
00041       const UnstableFinalState& parts = applyProjection<UnstableFinalState>(event, "UFS");
00042 
00043       foreach (const Particle& p, parts.particles()) {
00044         const double pT = p.pT();
00045         const double y = fabs(p.rapidity());
00046         const PdgId pid = abs(p.pdgId());
00047 
00048         if (y < 2.0) {
00049           switch (pid) {
00050           case PID::K0S:
00051               _h_dNKshort_dy->fill(y, weight);
00052               _h_dNKshort_dpT->fill(pT, weight);
00053             break;
00054           case PID::LAMBDA:
00055             // Lambda should not have Cascade or Omega ancestors since they should not decay. But just in case...
00056             if ( !( p.hasAncestor(3322) || p.hasAncestor(-3322) || p.hasAncestor(3312) || p.hasAncestor(-3312) || p.hasAncestor(3334) || p.hasAncestor(-3334) ) ) {
00057               _h_dNLambda_dy->fill(y, weight);
00058               _h_dNLambda_dpT->fill(pT, weight);
00059             }
00060             break;
00061           case PID::XIMINUS:
00062             // Cascade should not have Omega ancestors since it should not decay.  But just in case...
00063             if ( !( p.hasAncestor(3334) || p.hasAncestor(-3334) ) ) {
00064               _h_dNXi_dy->fill(y, weight);
00065               _h_dNXi_dpT->fill(pT, weight);
00066             }
00067             break;
00068           }
00069         }
00070       }
00071     }
00072 
00073 
00074     void finalize() {
00075 
00076       divide(_h_dNLambda_dpT,_h_dNKshort_dpT,
00077          _h_LampT_KpT);
00078 
00079       divide(_h_dNXi_dpT,_h_dNLambda_dpT,
00080          _h_XipT_LampT);
00081 
00082       divide(_h_dNLambda_dy,_h_dNKshort_dy,
00083          _h_Lamy_Ky);
00084 
00085       divide(_h_dNXi_dy,_h_dNLambda_dy,
00086          _h_Xiy_Lamy);
00087 
00088       double normpT = 1.0/sumOfWeights();
00089       double normy = 0.5*normpT; // Accounts for using |y| instead of y
00090       scale(_h_dNKshort_dy, normy);
00091       scale(_h_dNKshort_dpT, normpT);
00092       scale(_h_dNLambda_dy, normy);
00093       scale(_h_dNLambda_dpT, normpT);
00094       scale(_h_dNXi_dy, normy);
00095       scale(_h_dNXi_dpT, normpT);
00096     }
00097 
00098 
00099   private:
00100 
00101     // Particle distributions versus rapidity and transverse momentum
00102     Histo1DPtr _h_dNKshort_dy;
00103     Histo1DPtr _h_dNKshort_dpT;
00104     Histo1DPtr _h_dNLambda_dy;
00105     Histo1DPtr _h_dNLambda_dpT;
00106     Histo1DPtr _h_dNXi_dy;
00107     Histo1DPtr _h_dNXi_dpT;
00108 
00109     Scatter2DPtr _h_LampT_KpT;
00110     Scatter2DPtr _h_XipT_LampT;
00111     Scatter2DPtr _h_Lamy_Ky;
00112     Scatter2DPtr _h_Xiy_Lamy;
00113 
00114   };
00115 
00116 
00117 
00118   // The hook for the plugin system
00119   DECLARE_RIVET_PLUGIN(CMS_2011_S8978280);
00120 
00121 }