ATLAS_2010_S8914702.cc
Go to the documentation of this file.
00001 // -*- C++ -*- 00002 #include <iostream> 00003 #include <sstream> 00004 #include <string> 00005 00006 #include "Rivet/Analysis.hh" 00007 #include "Rivet/Projections/FinalState.hh" 00008 00009 #include "Rivet/Projections/LeadingParticlesFinalState.hh" 00010 #include "Rivet/Jet.hh" 00011 #include "Rivet/Projections/FastJets.hh" 00012 00013 #include "fastjet/internal/base.hh" 00014 #include "fastjet/JetDefinition.hh" 00015 #include "fastjet/AreaDefinition.hh" 00016 #include "fastjet/ClusterSequence.hh" 00017 #include "fastjet/ClusterSequenceArea.hh" 00018 #include "fastjet/PseudoJet.hh" 00019 00020 namespace Rivet { 00021 00022 00023 class ATLAS_2010_S8914702 : public Analysis { 00024 public: 00025 00026 /// Constructor 00027 ATLAS_2010_S8914702() 00028 : Analysis("ATLAS_2010_S8914702") 00029 { 00030 _eta_bins.push_back( 0.00); 00031 _eta_bins.push_back( 0.60); 00032 _eta_bins.push_back( 1.37); 00033 _eta_bins.push_back( 1.52); 00034 _eta_bins.push_back( 1.81); 00035 00036 _eta_bins_areaoffset.push_back(0.0); 00037 _eta_bins_areaoffset.push_back(1.5); 00038 _eta_bins_areaoffset.push_back(3.0); 00039 } 00040 00041 00042 /// Book histograms and initialise projections before the run 00043 void init() { 00044 FinalState fs; 00045 addProjection(fs, "FS"); 00046 00047 FastJets fj(fs, FastJets::KT, 0.5); 00048 _area_def = new fastjet::AreaDefinition(fastjet::VoronoiAreaSpec()); 00049 fj.useJetArea(_area_def); 00050 addProjection(fj, "KtJetsD05"); 00051 00052 LeadingParticlesFinalState photonfs(FinalState(-1.81, 1.81, 15.0*GeV)); 00053 photonfs.addParticleId(PID::PHOTON); 00054 addProjection(photonfs, "LeadingPhoton"); 00055 00056 int hist_bin = 0; 00057 for (int i = 0; i < (int)_eta_bins.size()-1; ++i) { 00058 if (fabs(_eta_bins[i] - 1.37) < .0001) continue; 00059 _h_Et_photon[i] = bookHisto1D(1, 1, hist_bin+1); 00060 hist_bin += 1; 00061 } 00062 } 00063 00064 00065 int getEtaBin(double eta_w, bool area_eta) const { 00066 double eta = fabs(eta_w); 00067 int v_iter = 0; 00068 if (!area_eta) { 00069 for (v_iter=0; v_iter < (int)_eta_bins.size()-1; ++v_iter) { 00070 if (eta >= _eta_bins.at(v_iter) && eta < _eta_bins.at(v_iter+1)) break; 00071 } 00072 } else { 00073 for (v_iter=0; v_iter < (int)_eta_bins_areaoffset.size()-1; ++v_iter) { 00074 if (eta >= _eta_bins_areaoffset.at(v_iter) && eta < _eta_bins_areaoffset.at(v_iter+1)) break; 00075 } 00076 } 00077 return v_iter; 00078 } 00079 00080 00081 /// Perform the per-event analysis 00082 void analyze(const Event& event) { 00083 const double weight = event.weight(); 00084 00085 Particles photons = applyProjection<LeadingParticlesFinalState>(event, "LeadingPhoton").particles(); 00086 if (photons.size() != 1) { 00087 vetoEvent; 00088 } 00089 00090 FourMomentum leadingPhoton = photons[0].momentum(); 00091 double eta_P = leadingPhoton.eta(); 00092 double phi_P = leadingPhoton.phi(); 00093 00094 if(fabs(eta_P)>=1.37 && fabs(eta_P)<1.52){ 00095 vetoEvent; 00096 } 00097 00098 int eta_bin = getEtaBin(eta_P,false); 00099 00100 Particles fs = applyProjection<FinalState>(event, "FS").particles(); 00101 FourMomentum mom_in_EtCone; 00102 foreach (const Particle& p, fs) { 00103 // check if it's in the cone of .4 00104 if (deltaR(eta_P, phi_P, p.eta(), p.momentum().phi()) >= 0.4) continue; 00105 00106 // check if it's in the 5x7 central core 00107 if (fabs(eta_P-p.eta()) < .025*7.0*0.5 && 00108 fabs(phi_P-p.momentum().phi()) < (PI/128.)*5.0*0.5) continue; 00109 mom_in_EtCone += p.momentum(); 00110 } 00111 MSG_DEBUG("Done with initial EtCone."); 00112 00113 // Now compute the median energy density 00114 _ptDensity.clear(); 00115 _sigma.clear(); 00116 _Njets.clear(); 00117 00118 std::vector< std::vector<double> > ptDensities; 00119 std::vector<double> emptyVec; 00120 ptDensities.assign(_eta_bins_areaoffset.size()-1,emptyVec); 00121 00122 const fastjet::ClusterSequenceArea* clust_seq_area = applyProjection<FastJets>(event, "KtJetsD05").clusterSeqArea(); 00123 foreach (const fastjet::PseudoJet& jet, applyProjection<FastJets>(event, "KtJetsD05").pseudoJets(0.0*GeV)) { 00124 //const double y = fabs(jet.rapidity()); 00125 const double eta = fabs(jet.eta()); 00126 const double pt = fabs(jet.perp()); 00127 00128 // Get the cluster sequence 00129 double area = clust_seq_area->area(jet); 00130 00131 if (area > 10e-4 && fabs(eta)<_eta_bins_areaoffset[_eta_bins_areaoffset.size()-1]) { 00132 ptDensities.at(getEtaBin(fabs(eta),true)).push_back(pt/area); 00133 } 00134 } 00135 00136 for (int b = 0; b < (int)_eta_bins_areaoffset.size()-1; ++b) { 00137 double median = 0.0; 00138 double sigma = 0.0; 00139 int Njets = 0; 00140 if (ptDensities[b].size() > 0) { 00141 std::sort(ptDensities[b].begin(), ptDensities[b].end()); 00142 int nDens = ptDensities[b].size(); 00143 if (nDens % 2 == 0) { 00144 median = (ptDensities[b][nDens/2]+ptDensities[b][(nDens-2)/2])/2; 00145 } else { 00146 median = ptDensities[b][(nDens-1)/2]; 00147 } 00148 sigma = ptDensities[b][(int)(.15865*nDens)]; 00149 Njets = nDens; 00150 } 00151 _ptDensity.push_back(median); 00152 _sigma.push_back(sigma); 00153 _Njets.push_back(Njets); 00154 } 00155 00156 // Now figure out the correction 00157 float EtCone_area = PI*.4*.4 - (7.0*.025)*(5.0*PI/128.); 00158 float correction = _ptDensity[getEtaBin(eta_P,true)]*EtCone_area; 00159 MSG_DEBUG("Jet area correction done."); 00160 00161 // Shouldn't need to subtract photon 00162 // NB. Using expected cut at hadron/particle level, not cut at reco level 00163 if (mom_in_EtCone.Et() - correction/*-leadingPhoton.Et()*/ > 4.0*GeV) { 00164 vetoEvent; 00165 } 00166 MSG_DEBUG("Passed isolation cut."); 00167 00168 _h_Et_photon[eta_bin]->fill(leadingPhoton.Et(), weight); 00169 } 00170 00171 00172 /// Normalise histograms etc., after the run 00173 void finalize() { 00174 for (int i = 0; i < (int)_eta_bins.size()-1; ++i) { 00175 if (fabs(_eta_bins[i] - 1.37) < .0001) continue; 00176 scale(_h_Et_photon[i], crossSection()/sumOfWeights()); 00177 } 00178 } 00179 00180 00181 private: 00182 00183 Histo1DPtr _h_Et_photon[6]; 00184 00185 fastjet::AreaDefinition* _area_def; 00186 00187 std::vector<float> _eta_bins; 00188 std::vector<float> _eta_bins_areaoffset; 00189 00190 std::vector<float> _ptDensity; 00191 std::vector<float> _sigma; 00192 std::vector<float> _Njets; 00193 }; 00194 00195 00196 00197 // The hook for the plugin system 00198 DECLARE_RIVET_PLUGIN(ATLAS_2010_S8914702); 00199 00200 } Generated on Fri Oct 25 2013 12:41:43 for The Rivet MC analysis system by ![]() |