CMS_2013_I1258128.cc
Go to the documentation of this file.
00001 #include "Rivet/Analysis.hh" 00002 #include "Rivet/Tools/BinnedHistogram.hh" 00003 #include "Rivet/Projections/FinalState.hh" 00004 #include "Rivet/Projections/FastJets.hh" 00005 #include "Rivet/Projections/ZFinder.hh" 00006 #include "Rivet/Projections/Thrust.hh" 00007 #include "Rivet/Projections/LeadingParticlesFinalState.hh" 00008 00009 namespace Rivet { 00010 00011 00012 /// CMS Z rapidity measurement 00013 class CMS_2013_I1258128 : public Analysis { 00014 public: 00015 00016 // Constructor 00017 CMS_2013_I1258128() 00018 : Analysis("CMS_2013_I1258128") 00019 { } 00020 00021 00022 void init() { 00023 // Full final state 00024 const FinalState fs(Cuts::abseta < 5); 00025 addProjection(fs, "FS"); 00026 00027 // Z finders for electrons and muons 00028 Cut cuts = Cuts::abseta < 2.1 && Cuts::pT > 20*GeV; 00029 const ZFinder zfe(fs, cuts, PID::ELECTRON, 76*GeV, 106*GeV); 00030 const ZFinder zfm(fs, cuts, PID::MUON, 76*GeV, 106*GeV); 00031 addProjection(zfe, "ZFE"); 00032 addProjection(zfm, "ZFM"); 00033 00034 // Try to get the leading photon 00035 LeadingParticlesFinalState photonfs(FinalState(-2.5, 2.5, 40.0*GeV)); 00036 photonfs.addParticleId(PID::PHOTON); 00037 addProjection(photonfs, "LeadingPhoton"); 00038 00039 // Jets 00040 const FastJets jets(fs, FastJets::ANTIKT, 0.5); 00041 addProjection(jets, "JETS"); 00042 00043 // Histograms 00044 _hist1YZ = bookHisto1D(1, 1, 1); 00045 _hist1YJet = bookHisto1D(2, 1, 1); 00046 _hist1YSum = bookHisto1D(3, 1, 1); 00047 _hist1YDif = bookHisto1D(4, 1, 1); 00048 _hist2YPhoton = bookHisto1D(5, 1, 1); 00049 _hist2YJet = bookHisto1D(6, 1, 1); 00050 _hist2YSum = bookHisto1D(7, 1, 1); 00051 _hist2YDif = bookHisto1D(8, 1, 1); 00052 } 00053 00054 00055 void makeZCut(const Event& event) { 00056 // Apply the Z finders and veto if no Z found 00057 const ZFinder& zfe = applyProjection<ZFinder>(event, "ZFE"); 00058 const ZFinder& zfm = applyProjection<ZFinder>(event, "ZFM"); 00059 if (zfe.empty() && zfm.empty()) vetoEvent; 00060 00061 // Choose the Z candidate 00062 const ParticleVector& z = (!zfm.empty()) ? zfm.bosons() : zfe.bosons(); 00063 const ParticleVector& clusteredConstituents = (!zfm.empty()) ? zfm.constituents() : zfe.constituents(); 00064 00065 // Insist that the Z is in a high-pT (boosted) regime 00066 if (z[0].pT() < 40*GeV) return; 00067 00068 // Build the jets 00069 const FastJets& jetfs = applyProjection<FastJets>(event, "JETS"); 00070 Jets jets = jetfs.jetsByPt(Cuts::pT > 30*GeV && Cuts::abseta < 2.4); 00071 if (jets.empty()) return; 00072 00073 // Clean the jets against the lepton candidates with a DeltaR cut of 0.5 00074 vector<const Jet*> cleanedJets; 00075 foreach (const Jet& j, jets) { 00076 bool isolated = true; 00077 foreach (const Particle& p, clusteredConstituents) { 00078 if (deltaR(p, j) < 0.5) { 00079 isolated = false; 00080 break; 00081 } 00082 } 00083 if (isolated) cleanedJets.push_back(&j); 00084 } 00085 // Require exactly 1 isolated jet 00086 if (cleanedJets.size() != 1) return; 00087 00088 // Fill histos 00089 const double weight = event.weight(); 00090 const double yz = z[0].rapidity(); 00091 const double yjet = cleanedJets[0]->momentum().rapidity(); 00092 _hist1YZ->fill(fabs(yz), weight); 00093 _hist1YJet->fill(fabs(yjet), weight); 00094 _hist1YSum->fill(0.5*fabs(yz + yjet), weight); 00095 _hist1YDif->fill(0.5*fabs(yz - yjet), weight); 00096 } 00097 00098 00099 void makePhotonCut(const Event& event) { 00100 // Get the photon 00101 const FinalState& photonfs = applyProjection<FinalState>(event, "LeadingPhoton"); 00102 if (photonfs.particles().size() < 1) return; 00103 const Particle& photon = photonfs.particles().front(); 00104 if (photon.pT() < 40*GeV) return; 00105 if (fabs(photon.eta()) > 1.4442 ) return; 00106 00107 // Build the jets 00108 const FastJets& jetfs = applyProjection<FastJets>(event, "JETS"); 00109 Jets jets = jetfs.jetsByPt(Cuts::pT > 30*GeV && Cuts::abseta < 2.4); 00110 if (jets.empty()) return; 00111 00112 // Clean the jets against the photon candidate with a DeltaR cut of 0.5 00113 vector<const Jet*> cleanedJets; 00114 foreach (const Jet& j, jets) 00115 if (deltaR(photon, j) > 0.5) 00116 cleanedJets.push_back(&j); 00117 // Require exactly 1 jet 00118 if (cleanedJets.size() != 1) return; 00119 00120 // Fill histos 00121 const double weight = event.weight(); 00122 const double ypho = photon.rapidity(); 00123 const double yjet = cleanedJets[0]->momentum().rapidity(); 00124 _hist2YPhoton->fill(fabs(ypho), weight); 00125 _hist2YJet->fill(fabs(yjet), weight); 00126 _hist2YSum->fill(0.5*fabs(ypho + yjet), weight); 00127 _hist2YDif->fill(0.5*fabs(ypho - yjet), weight); 00128 } 00129 00130 00131 void analyze(const Event& event) { 00132 makeZCut(event); 00133 makePhotonCut(event); 00134 } 00135 00136 00137 void finalize() { 00138 normalizeByContents(_hist1YZ); 00139 normalizeByContents(_hist1YJet); 00140 normalizeByContents(_hist1YSum); 00141 normalizeByContents(_hist1YDif); 00142 normalizeByContents(_hist2YPhoton); 00143 normalizeByContents(_hist2YJet); 00144 normalizeByContents(_hist2YSum); 00145 normalizeByContents(_hist2YDif); 00146 } 00147 00148 00149 // The CMS normalization in this analysis is that the sum over bin contents 00150 // is equal to 1. This function normalizes to area = area*bin_width. / 00151 // @note This is a strange definition... why? 00152 void normalizeByContents(Histo1DPtr h) { 00153 normalize(h, h->bin(0).xWidth()); 00154 } 00155 00156 00157 private: 00158 00159 Histo1DPtr _hist1YZ, _hist1YJet, _hist1YSum, _hist1YDif; 00160 Histo1DPtr _hist2YPhoton, _hist2YJet, _hist2YSum, _hist2YDif; 00161 00162 }; 00163 00164 00165 // Plugin system hook 00166 DECLARE_RIVET_PLUGIN(CMS_2013_I1258128); 00167 00168 } Generated on Tue Mar 24 2015 17:35:26 for The Rivet MC analysis system by ![]() |