CMS_2011_S9120041.cc
Go to the documentation of this file.
00001 // -*- C++ -*- 00002 00003 #include "Rivet/Analysis.hh" 00004 #include "Rivet/RivetYODA.hh" 00005 #include "Rivet/Tools/Logging.hh" 00006 #include "Rivet/Projections/FinalState.hh" 00007 #include "Rivet/Projections/ChargedFinalState.hh" 00008 #include "Rivet/Projections/FastJets.hh" 00009 00010 using namespace std; 00011 00012 namespace Rivet { 00013 00014 // UE charged particles vs. leading jet 00015 class CMS_2011_S9120041 : public Analysis { 00016 public: 00017 00018 /// Constructor 00019 CMS_2011_S9120041() : Analysis("CMS_2011_S9120041") {} 00020 00021 00022 void init() { 00023 const ChargedFinalState cfs(-2.0, 2.0, 500*MeV); 00024 addProjection(cfs, "CFS"); 00025 00026 const ChargedFinalState cfsforjet(-2.5, 2.5, 500*MeV); 00027 addProjection(cfsforjet, "CFSforjet"); 00028 00029 const FastJets jetpro(cfsforjet, FastJets::SISCONE, 0.5); 00030 addProjection(jetpro, "Jets"); 00031 00032 if (fuzzyEquals(sqrtS(), 7.0*TeV)) { 00033 _h_Nch_vs_pT = bookProfile1D(1, 1, 1); // Nch vs. pT_max 00034 _h_Sum_vs_pT = bookProfile1D(2, 1, 1); // sum(pT) vs. pT_max 00035 _h_pT3_Nch = bookHisto1D(5, 1, 1); // transverse Nch, pT_max > 3GeV 00036 _h_pT3_Sum = bookHisto1D(6, 1, 1); // transverse sum(pT), pT_max > 3GeV 00037 _h_pT3_pT = bookHisto1D(7, 1, 1); // transverse pT, pT_max > 3GeV 00038 _h_pT20_Nch = bookHisto1D(8, 1, 1); // transverse Nch, pT_max > 20GeV 00039 _h_pT20_Sum = bookHisto1D(9, 1, 1); // transverse sum(pT), pT_max > 20GeV 00040 _h_pT20_pT = bookHisto1D(10, 1, 1); // transverse pT, pT_max > 20GeV 00041 } 00042 00043 if (fuzzyEquals(sqrtS(), 0.9*TeV)) { 00044 _h_Nch_vs_pT = bookProfile1D(3, 1, 1); // Nch vs. pT_max 00045 _h_Sum_vs_pT = bookProfile1D(4, 1, 1); // sum(pT) vs. pT_max 00046 _h_pT3_Nch = bookHisto1D(11, 1, 1); // transverse Nch, pT_max > 3GeV 00047 _h_pT3_Sum = bookHisto1D(12, 1, 1); // transverse sum(pT), pT_max > 3GeV 00048 _h_pT3_pT = bookHisto1D(13, 1, 1); // transverse pT, pT_max > 3GeV 00049 } 00050 00051 sumOfWeights3 = 0.0; 00052 sumOfWeights20 = 0.0; 00053 00054 _nch_tot_pT3 = 0.0; 00055 _nch_tot_pT20 = 0.0; 00056 } 00057 00058 00059 /// Perform the per-event analysis 00060 void analyze(const Event& event) { 00061 const double weight = event.weight(); 00062 00063 Jets jets = applyProjection<FastJets>(event, "Jets").jetsByPt(1.0*GeV); 00064 if (jets.size() < 1) vetoEvent; 00065 00066 FourMomentum p_lead = jets[0].momentum(); 00067 const double philead = p_lead.phi(); 00068 const double pTlead = p_lead.pT(); 00069 00070 ParticleVector particles = applyProjection<ChargedFinalState>(event, "CFS").particlesByPt(); 00071 00072 int nTransverse = 0; 00073 double ptSumTransverse = 0.; 00074 foreach (const Particle& p, particles) { 00075 double dphi = fabs(deltaPhi(philead, p.momentum().phi())); 00076 if (dphi>PI/3. && dphi<PI*2./3.) { // Transverse region 00077 nTransverse++; 00078 00079 const double pT = p.momentum().pT()/GeV; 00080 ptSumTransverse += pT; 00081 00082 if (pTlead > 3.0*GeV) _h_pT3_pT->fill(pT, weight); 00083 if (fuzzyEquals(sqrtS(), 7.0*TeV) && pTlead > 20.0*GeV) _h_pT20_pT->fill(pT, weight); 00084 } 00085 } 00086 00087 const double area = 8./3. * PI; 00088 _h_Nch_vs_pT->fill(pTlead/GeV, 1./area*nTransverse, weight); 00089 _h_Sum_vs_pT->fill(pTlead/GeV, 1./area*ptSumTransverse, weight); 00090 if(pTlead > 3.0*GeV) { 00091 _h_pT3_Nch->fill(nTransverse, weight); 00092 _h_pT3_Sum->fill(ptSumTransverse, weight); 00093 sumOfWeights3 += weight; 00094 _nch_tot_pT3 += weight*nTransverse; 00095 } 00096 if (fuzzyEquals(sqrtS(), 7.0*TeV) && pTlead > 20.0*GeV) { 00097 _h_pT20_Nch->fill(nTransverse, weight); 00098 _h_pT20_Sum->fill(ptSumTransverse, weight); 00099 sumOfWeights20 += weight; 00100 _nch_tot_pT20 += weight*nTransverse; 00101 } 00102 } 00103 00104 00105 00106 /// Normalise histograms etc., after the run 00107 void finalize() { 00108 normalize(_h_pT3_Nch); 00109 normalize(_h_pT3_Sum); 00110 if (sumOfWeights3 != 0.0) normalize(_h_pT3_pT, _nch_tot_pT3 / sumOfWeights3); 00111 00112 if (fuzzyEquals(sqrtS(), 7.0*TeV)) { 00113 normalize(_h_pT20_Nch); 00114 normalize(_h_pT20_Sum); 00115 if (sumOfWeights20 != 0.0) normalize(_h_pT20_pT, _nch_tot_pT20 / sumOfWeights20); 00116 } 00117 } 00118 00119 00120 00121 private: 00122 00123 double sumOfWeights3; 00124 double sumOfWeights20; 00125 00126 double _nch_tot_pT3; 00127 double _nch_tot_pT20; 00128 00129 Profile1DPtr _h_Nch_vs_pT; 00130 Profile1DPtr _h_Sum_vs_pT; 00131 Histo1DPtr _h_pT3_Nch; 00132 Histo1DPtr _h_pT3_Sum; 00133 Histo1DPtr _h_pT3_pT; 00134 Histo1DPtr _h_pT20_Nch; 00135 Histo1DPtr _h_pT20_Sum; 00136 Histo1DPtr _h_pT20_pT; 00137 00138 }; 00139 00140 00141 // This global object acts as a hook for the plugin system 00142 DECLARE_RIVET_PLUGIN(CMS_2011_S9120041); 00143 } 00144 Generated on Fri Dec 21 2012 15:03:40 for The Rivet MC analysis system by ![]() |