CMS_2011_S8978280.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 #include "Rivet/Analysis.hh"
00003 #include "Rivet/RivetAIDA.hh"
00004 #include "Rivet/Tools/Logging.hh"
00005 #include "Rivet/Projections/UnstableFinalState.hh"
00006 
00007 namespace Rivet {
00008 
00009 
00010   /// @brief CMS strange particle spectra (Ks, Lambda, Cascade) in pp at 900 and 7000 GeV
00011   /// @author Kevin Stenson
00012   class CMS_2011_S8978280 : public Analysis {
00013   public:
00014 
00015     /// Constructor
00016     CMS_2011_S8978280() : Analysis("CMS_2011_S8978280") {}
00017 
00018 
00019     void init() {
00020       // Need wide range of eta because cut on rapidity not pseudorapidity
00021       UnstableFinalState ufs(-8.0, 8.0, 0.0*GeV);
00022       addProjection(ufs, "UFS");
00023 
00024       // Particle distributions versus rapidity and transverse momentum
00025       // Only make histograms if the correct energy is used.
00026       if (fuzzyEquals(sqrtS()/GeV, 900)){
00027         _h_dNKshort_dy  = bookHistogram1D(1, 1, 1);
00028         _h_dNKshort_dpT = bookHistogram1D(2, 1, 1);
00029         _h_dNLambda_dy  = bookHistogram1D(3, 1, 1);
00030         _h_dNLambda_dpT = bookHistogram1D(4, 1, 1);
00031         _h_dNXi_dy      = bookHistogram1D(5, 1, 1);
00032         _h_dNXi_dpT     = bookHistogram1D(6, 1, 1);
00033       } else if (fuzzyEquals(sqrtS()/GeV, 7000)){
00034         _h_dNKshort_dy  = bookHistogram1D(1, 1, 2);
00035         _h_dNKshort_dpT = bookHistogram1D(2, 1, 2);
00036         _h_dNLambda_dy  = bookHistogram1D(3, 1, 2);
00037         _h_dNLambda_dpT = bookHistogram1D(4, 1, 2);
00038         _h_dNXi_dy      = bookHistogram1D(5, 1, 2);
00039         _h_dNXi_dpT     = bookHistogram1D(6, 1, 2);
00040       }
00041     }
00042 
00043 
00044     void analyze(const Event& event) {
00045       const double weight = event.weight();
00046 
00047       const UnstableFinalState& parts = applyProjection<UnstableFinalState>(event, "UFS");
00048 
00049       foreach (const Particle& p, parts.particles()) {
00050         const double pT = p.momentum().pT();
00051         const double y = fabs(p.momentum().rapidity());
00052         const PdgId pid = abs(p.pdgId());
00053 
00054         if (y < 2.0) {
00055           switch (pid) {
00056           case K0S:
00057               _h_dNKshort_dy->fill(y, weight);
00058               _h_dNKshort_dpT->fill(pT, weight);
00059             break;
00060           case LAMBDA:
00061             // Lambda should not have Cascade or Omega ancestors since they should not decay. But just in case...
00062             if ( !( p.hasAncestor(3322) || p.hasAncestor(-3322) || p.hasAncestor(3312) || p.hasAncestor(-3312) || p.hasAncestor(3334) || p.hasAncestor(-3334) ) ) {
00063               _h_dNLambda_dy->fill(y, weight);
00064               _h_dNLambda_dpT->fill(pT, weight);
00065             }
00066             break;
00067           case XIMINUS:
00068             // Cascade should not have Omega ancestors since it should not decay.  But just in case...
00069             if ( !( p.hasAncestor(3334) || p.hasAncestor(-3334) ) ) {
00070               _h_dNXi_dy->fill(y, weight);
00071               _h_dNXi_dpT->fill(pT, weight);
00072             }
00073             break;
00074           }
00075         }
00076       }
00077     }
00078 
00079 
00080     void finalize() {
00081       AIDA::IHistogramFactory& hf = histogramFactory();
00082       const string dir = histoDir();
00083 
00084       // Making the Lambda/Kshort and Xi/Lambda ratios vs pT and y
00085       if (fuzzyEquals(sqrtS()/GeV, 900)) {
00086         hf.divide(dir + "/d07-x01-y01",*_h_dNLambda_dpT, *_h_dNKshort_dpT);
00087         hf.divide(dir + "/d08-x01-y01",*_h_dNXi_dpT, *_h_dNLambda_dpT);
00088         hf.divide(dir + "/d09-x01-y01",*_h_dNLambda_dy, *_h_dNKshort_dy);
00089         hf.divide(dir + "/d10-x01-y01",*_h_dNXi_dy, *_h_dNLambda_dy);
00090       } else if (fuzzyEquals(sqrtS()/GeV, 7000)) {
00091         hf.divide(dir + "/d07-x01-y02",*_h_dNLambda_dpT, *_h_dNKshort_dpT);
00092         hf.divide(dir + "/d08-x01-y02",*_h_dNXi_dpT, *_h_dNLambda_dpT);
00093         hf.divide(dir + "/d09-x01-y02",*_h_dNLambda_dy, *_h_dNKshort_dy);
00094         hf.divide(dir + "/d10-x01-y02",*_h_dNXi_dy, *_h_dNLambda_dy);
00095       }
00096 
00097       double normpT = 1.0/sumOfWeights();
00098       double normy = 0.5*normpT; // Accounts for using |y| instead of y
00099       scale(_h_dNKshort_dy, normy);
00100       scale(_h_dNKshort_dpT, normpT);
00101       scale(_h_dNLambda_dy, normy);
00102       scale(_h_dNLambda_dpT, normpT);
00103       scale(_h_dNXi_dy, normy);
00104       scale(_h_dNXi_dpT, normpT);
00105     }
00106 
00107 
00108   private:
00109 
00110     // Particle distributions versus rapidity and transverse momentum
00111     AIDA::IHistogram1D *_h_dNKshort_dy;
00112     AIDA::IHistogram1D *_h_dNKshort_dpT;
00113     AIDA::IHistogram1D *_h_dNLambda_dy;
00114     AIDA::IHistogram1D *_h_dNLambda_dpT;
00115     AIDA::IHistogram1D *_h_dNXi_dy;
00116     AIDA::IHistogram1D *_h_dNXi_dpT;
00117 
00118   };
00119 
00120 
00121 
00122   // The hook for the plugin system
00123   DECLARE_RIVET_PLUGIN(CMS_2011_S8978280);
00124 
00125 }