ATLAS_2010_S8914702.cc
Go to the documentation of this file.
00001 // -*- C++ -*- 00002 #include "Rivet/Analysis.hh" 00003 #include "Rivet/Projections/FinalState.hh" 00004 #include "Rivet/Projections/LeadingParticlesFinalState.hh" 00005 #include "Rivet/Projections/FastJets.hh" 00006 00007 namespace Rivet { 00008 00009 00010 class ATLAS_2010_S8914702 : public Analysis { 00011 public: 00012 00013 /// Constructor 00014 ATLAS_2010_S8914702() 00015 : Analysis("ATLAS_2010_S8914702") 00016 { } 00017 00018 00019 /// Book histograms and initialise projections before the run 00020 void init() { 00021 FinalState fs; 00022 declare(fs, "FS"); 00023 00024 FastJets fj(fs, FastJets::KT, 0.5); 00025 fj.useJetArea(new fastjet::AreaDefinition(fastjet::VoronoiAreaSpec())); 00026 declare(fj, "KtJetsD05"); 00027 00028 LeadingParticlesFinalState photonfs(FinalState(Cuts::abseta < 1.81 && Cuts::pT > 15*GeV)); 00029 photonfs.addParticleId(PID::PHOTON); 00030 declare(photonfs, "LeadingPhoton"); 00031 00032 size_t hist_bin = 0; 00033 for (size_t i = 0; i < _eta_bins.size()-1; ++i) { 00034 if (fabs(_eta_bins[i] - 1.37) < .0001) continue; 00035 _h_Et_photon[i] = bookHisto1D(1, 1, hist_bin+1); 00036 hist_bin += 1; 00037 } 00038 } 00039 00040 00041 size_t getEtaBin(double eta, bool area_eta) const { 00042 return (!area_eta) ? binIndex(fabs(eta), _eta_bins) : binIndex(fabs(eta), _eta_bins_areaoffset); 00043 } 00044 00045 00046 /// Perform the per-event analysis 00047 void analyze(const Event& event) { 00048 Particles photons = apply<LeadingParticlesFinalState>(event, "LeadingPhoton").particles(); 00049 if (photons.size() != 1) vetoEvent; 00050 00051 const Particle& leadingPhoton = photons[0]; 00052 if (inRange(leadingPhoton.abseta(), 1.37, 1.52)) vetoEvent; 00053 const int eta_bin = getEtaBin(leadingPhoton.abseta(), false); 00054 00055 const Particles& fs = apply<FinalState>(event, "FS").particles(); 00056 FourMomentum mom_in_EtCone; 00057 for (const Particle& p : fs) { 00058 // Check if it's in the cone of .4 00059 if (deltaR(leadingPhoton, p) >= 0.4) continue; 00060 // Check if it's in the 5x7 central core 00061 if (fabs(deltaEta(leadingPhoton, p)) < .025*5.0*0.5 && 00062 fabs(deltaPhi(leadingPhoton, p)) < (PI/128.)*7.0*0.5) continue; 00063 // Increment 00064 mom_in_EtCone += p.momentum(); 00065 } 00066 MSG_DEBUG("Done with initial Et cone"); 00067 00068 // Get the jet pT densities 00069 vector< vector<double> > ptDensities(_eta_bins_areaoffset.size()-1); 00070 FastJets fastjets = apply<FastJets>(event, "KtJetsD05"); 00071 const shared_ptr<fastjet::ClusterSequenceArea> clust_seq_area = fastjets.clusterSeqArea(); 00072 for (const Jet& jet : fastjets.jets()) { 00073 const double area = clust_seq_area->area(jet); //< Implicit call to pseudojet() 00074 if (area > 1e-4 && jet.abseta() < _eta_bins_areaoffset.back()) { 00075 ptDensities.at(getEtaBin(jet.abseta(), true)) += jet.pT()/area; 00076 } 00077 } 00078 00079 // Now compute the median energy densities 00080 vector<double> ptDensity; 00081 for (size_t b = 0; b < _eta_bins_areaoffset.size()-1; ++b) { 00082 ptDensity += ptDensities[b].empty() ? 0 : median(ptDensities[b]); 00083 } 00084 00085 // Now figure out the correction 00086 const double ETCONE_AREA = PI*.4*.4 - (7.0*.025)*(5.0*PI/128.); 00087 const double correction = ptDensity[getEtaBin(leadingPhoton.abseta(), true)]*ETCONE_AREA; 00088 MSG_DEBUG("Jet area correction done"); 00089 00090 // Shouldn't need to subtract photon 00091 // NB. Using expected cut at hadron/particle level, not cut at reco level 00092 if (mom_in_EtCone.Et() - correction/*-leadingPhoton.Et()*/ > 4.0*GeV) vetoEvent; 00093 MSG_DEBUG("Passed isolation cut"); 00094 00095 // Fill histogram 00096 _h_Et_photon[eta_bin]->fill(leadingPhoton.Et()/GeV, event.weight()); 00097 } 00098 00099 00100 /// Normalise histograms etc., after the run 00101 void finalize() { 00102 for (size_t i = 0; i < _eta_bins.size()-1; ++i) { 00103 if (fabs(_eta_bins[i] - 1.37) < .0001) continue; 00104 scale(_h_Et_photon[i], crossSection()/sumOfWeights()); 00105 } 00106 } 00107 00108 00109 private: 00110 00111 Histo1DPtr _h_Et_photon[6]; 00112 00113 const vector<double> _eta_bins = {0.00, 0.60, 1.37, 1.52, 1.81}; 00114 const vector<double> _eta_bins_areaoffset = {0.0, 1.5, 3.0}; 00115 00116 }; 00117 00118 00119 00120 // The hook for the plugin system 00121 DECLARE_RIVET_PLUGIN(ATLAS_2010_S8914702); 00122 00123 } Generated on Tue Dec 13 2016 16:32:34 for The Rivet MC analysis system by ![]() |