ATLAS_2011_S9120807.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/RivetYODA.hh" 00008 #include "Rivet/Tools/Logging.hh" 00009 #include "Rivet/Projections/FinalState.hh" 00010 00011 #include "Rivet/Projections/IdentifiedFinalState.hh" 00012 #include "Rivet/Jet.hh" 00013 #include "Rivet/Projections/FastJets.hh" 00014 00015 #include "fastjet/internal/base.hh" 00016 #include "fastjet/JetDefinition.hh" 00017 #include "fastjet/AreaDefinition.hh" 00018 #include "fastjet/ClusterSequence.hh" 00019 #include "fastjet/ClusterSequenceArea.hh" 00020 #include "fastjet/PseudoJet.hh" 00021 00022 namespace Rivet { 00023 00024 00025 /// @brief Measurement of isolated diphoton + X differential cross-sections 00026 /// 00027 /// Inclusive isolated gamma gamma cross-sections, differential in M(gg), pT(gg), 00028 /// dphi(gg) 00029 /// 00030 /// @author Giovanni Marchiori 00031 class ATLAS_2011_S9120807 : public Analysis { 00032 public: 00033 00034 /// Constructor 00035 ATLAS_2011_S9120807() 00036 : Analysis("ATLAS_2011_S9120807") 00037 { 00038 _eta_bins_areaoffset.push_back(0.0); 00039 _eta_bins_areaoffset.push_back(1.5); 00040 _eta_bins_areaoffset.push_back(3.0); 00041 } 00042 00043 00044 public: 00045 00046 /// Book histograms and initialise projections before the run 00047 void init() { 00048 FinalState fs; 00049 addProjection(fs, "FS"); 00050 00051 FastJets fj(fs, FastJets::KT, 0.5); 00052 _area_def = new fastjet::AreaDefinition(fastjet::VoronoiAreaSpec()); 00053 fj.useJetArea(_area_def); 00054 addProjection(fj, "KtJetsD05"); 00055 00056 IdentifiedFinalState photonfs(-2.37, 2.37, 16.0*GeV); 00057 photonfs.acceptId(PHOTON); 00058 addProjection(photonfs, "Photon"); 00059 00060 _h_M = bookHisto1D(1, 1, 1); 00061 _h_pT = bookHisto1D(2, 1, 1); 00062 _h_dPhi = bookHisto1D(3, 1, 1); 00063 } 00064 00065 00066 int getEtaBin(double eta_w) const { 00067 double eta = fabs(eta_w); 00068 00069 int v_iter=0; 00070 for (v_iter=0; v_iter < (int)_eta_bins_areaoffset.size()-1; v_iter++) { 00071 if (inRange(eta, _eta_bins_areaoffset[v_iter], _eta_bins_areaoffset[v_iter+1])) 00072 break; 00073 } 00074 return v_iter; 00075 } 00076 00077 00078 /// Perform the per-event analysis 00079 void analyze(const Event& event) { 00080 00081 const double weight = event.weight(); 00082 00083 /// 00084 /// require at least 2 photons in final state 00085 /// 00086 ParticleVector photons = applyProjection<IdentifiedFinalState>(event, "Photon").particlesByPt(); 00087 if (photons.size() < 2) { 00088 vetoEvent; 00089 } 00090 00091 /// 00092 /// compute the median energy density 00093 /// 00094 std::vector<double> _ptDensity; 00095 std::vector< std::vector<double> > ptDensities; 00096 std::vector<double> emptyVec; 00097 ptDensities.assign(_eta_bins_areaoffset.size()-1, emptyVec); 00098 00099 const fastjet::ClusterSequenceArea* clust_seq_area = applyProjection<FastJets>(event, "KtJetsD05").clusterSeqArea(); 00100 foreach (const fastjet::PseudoJet& jet, applyProjection<FastJets>(event, "KtJetsD05").pseudoJets(0.0*GeV)) { 00101 double eta = fabs(jet.eta()); 00102 double pt = fabs(jet.perp()); 00103 00104 /// get the cluster sequence 00105 double area = clust_seq_area->area(jet); 00106 00107 if(area > 10e-4 && fabs(eta) < _eta_bins_areaoffset[_eta_bins_areaoffset.size()-1]) { 00108 ptDensities.at(getEtaBin(fabs(eta))).push_back(pt/area); 00109 } 00110 } 00111 00112 for(int b=0; b<(int)_eta_bins_areaoffset.size()-1; b++) { 00113 double median = 0.0; 00114 if(ptDensities[b].size() > 0) { 00115 std::sort(ptDensities[b].begin(), ptDensities[b].end()); 00116 int nDens = ptDensities[b].size(); 00117 if( nDens%2 == 0 ) 00118 median = (ptDensities[b][nDens/2] + ptDensities[b][(nDens-2)/2]) / 2; 00119 else 00120 median = ptDensities[b][(nDens-1)/2]; 00121 } 00122 _ptDensity.push_back(median); 00123 } 00124 00125 00126 /// 00127 /// Loop over photons and fill vector of isolated ones 00128 /// 00129 ParticleVector isolated_photons; 00130 foreach (const Particle& photon, photons) { 00131 00132 /// 00133 /// remove photons in crack 00134 /// 00135 double eta_P = photon.momentum().eta(); 00136 if (fabs(eta_P)>=1.37 && fabs(eta_P)<1.52) continue; 00137 00138 double phi_P = photon.momentum().phi(); 00139 00140 /// 00141 /// compute isolation 00142 /// 00143 00144 /// std EtCone 00145 ParticleVector fs = applyProjection<FinalState>(event, "FS").particles(); 00146 FourMomentum mom_in_EtCone; 00147 foreach (const Particle& p, fs) { 00148 /// check if it's in the cone of .4 00149 if (deltaR(eta_P, phi_P, p.momentum().eta(), p.momentum().phi()) >= 0.4) continue; 00150 00151 /// check if it's in the 5x7 central core 00152 if (fabs(eta_P-p.momentum().eta()) < .025*7.0*0.5 && 00153 fabs(phi_P-p.momentum().phi()) < (M_PI/128.)*5.0*0.5) continue; 00154 00155 mom_in_EtCone += p.momentum(); 00156 } 00157 00158 /// now figure out the correction (area*density) 00159 double EtCone_area = M_PI*.4*.4 - (7.0*.025)*(5.0*M_PI/128.); 00160 double correction = _ptDensity[getEtaBin(eta_P)]*EtCone_area; 00161 00162 /// shouldn't need to subtract photon 00163 /// note: using expected cut at hadron/particle level, not cut at reco level 00164 if (mom_in_EtCone.Et()-correction > 4.0*GeV) { 00165 continue; 00166 } 00167 00168 /// add photon to list of isolated ones 00169 isolated_photons.push_back(photon); 00170 } 00171 00172 /// 00173 /// require at least two isolated photons 00174 /// 00175 if (isolated_photons.size() < 2) { 00176 vetoEvent; 00177 } 00178 00179 /// 00180 /// select leading pT pair 00181 /// 00182 std::sort(isolated_photons.begin(), isolated_photons.end(), cmpParticleByPt); 00183 FourMomentum y1=isolated_photons[0].momentum(); 00184 FourMomentum y2=isolated_photons[1].momentum(); 00185 00186 /// 00187 /// require the two photons to be separated (dR>0.4) 00188 /// 00189 if (deltaR(y1, y2)<0.4) { 00190 vetoEvent; 00191 } 00192 00193 FourMomentum yy = y1+y2; 00194 double Myy = yy.mass()/GeV; 00195 double pTyy = yy.pT()/GeV; 00196 double dPhiyy = deltaPhi(y1.phi(), y2.phi()); 00197 00198 _h_M->fill(Myy, weight); 00199 _h_pT->fill(pTyy, weight); 00200 _h_dPhi->fill(dPhiyy, weight); 00201 } 00202 00203 00204 /// Normalise histograms etc., after the run 00205 void finalize() { 00206 scale(_h_M, crossSection()/sumOfWeights()); 00207 scale(_h_pT, crossSection()/sumOfWeights()); 00208 scale(_h_dPhi, crossSection()/sumOfWeights()); 00209 } 00210 00211 private: 00212 00213 Histo1DPtr _h_M; 00214 Histo1DPtr _h_pT; 00215 Histo1DPtr _h_dPhi; 00216 00217 fastjet::AreaDefinition* _area_def; 00218 00219 std::vector<double> _eta_bins_areaoffset; 00220 }; 00221 00222 00223 // The hook for the plugin system 00224 DECLARE_RIVET_PLUGIN(ATLAS_2011_S9120807); 00225 00226 } Generated on Fri Dec 21 2012 15:03:38 for The Rivet MC analysis system by ![]() |