ATLAS_2012_I1093738.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 #include "Rivet/Projections/LeadingParticlesFinalState.hh" 00009 #include "Rivet/Projections/VetoedFinalState.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 /// @brief Measurement of isolated gamma + jet + X differential cross-sections 00024 /// 00025 /// Inclusive isolated gamma + jet cross-sections, differential in pT(gamma), for 00026 /// various photon and jet rapidity configurations. 00027 /// 00028 /// @author Giovanni Marchiori 00029 class ATLAS_2012_I1093738 : public Analysis { 00030 public: 00031 00032 // Constructor 00033 ATLAS_2012_I1093738() 00034 : Analysis("ATLAS_2012_I1093738") 00035 { 00036 _eta_bins_ph.push_back(0.0); 00037 _eta_bins_ph.push_back(1.37); 00038 _eta_bins_ph.push_back(1.52); 00039 _eta_bins_ph.push_back(2.37); 00040 00041 _eta_bins_jet.push_back(0.0); 00042 _eta_bins_jet.push_back(1.2); 00043 _eta_bins_jet.push_back(2.8); 00044 _eta_bins_jet.push_back(4.4); 00045 00046 _eta_bins_areaoffset.push_back(0.0); 00047 _eta_bins_areaoffset.push_back(1.5); 00048 _eta_bins_areaoffset.push_back(3.0); 00049 } 00050 00051 public: 00052 00053 // Book histograms and initialise projections before the run 00054 void init() { 00055 // Final state 00056 FinalState fs; 00057 addProjection(fs, "FS"); 00058 00059 // Voronoi eta-phi tassellation with KT jets, for ambient energy density calculation 00060 FastJets fj(fs, FastJets::KT, 0.5); 00061 _area_def = new fastjet::AreaDefinition(fastjet::VoronoiAreaSpec()); 00062 fj.useJetArea(_area_def); 00063 addProjection(fj, "KtJetsD05"); 00064 00065 // Leading photon 00066 LeadingParticlesFinalState photonfs(FinalState(-1.37, 1.37, 25.0*GeV)); 00067 photonfs.addParticleId(PID::PHOTON); 00068 addProjection(photonfs, "LeadingPhoton"); 00069 00070 // FS excluding the leading photon 00071 VetoedFinalState vfs(fs); 00072 vfs.addVetoOnThisFinalState(photonfs); 00073 addProjection(vfs, "JetFS"); 00074 00075 // Jets 00076 FastJets jetpro(vfs, FastJets::ANTIKT, 0.4); 00077 jetpro.useInvisibles(); 00078 addProjection(jetpro, "Jets"); 00079 00080 _h_phbarrel_jetcentral_SS = bookHisto1D(1, 1, 1); 00081 _h_phbarrel_jetmedium_SS = bookHisto1D(2, 1, 1); 00082 _h_phbarrel_jetforward_SS = bookHisto1D(3, 1, 1); 00083 00084 _h_phbarrel_jetcentral_OS = bookHisto1D(4, 1, 1); 00085 _h_phbarrel_jetmedium_OS = bookHisto1D(5, 1, 1); 00086 _h_phbarrel_jetforward_OS = bookHisto1D(6, 1, 1); 00087 } 00088 00089 00090 int getEtaBin(double eta_w, int what) const { 00091 double eta = fabs(eta_w); 00092 int v_iter = 0; 00093 if (what==0) { 00094 for (v_iter=0; v_iter < (int)_eta_bins_ph.size()-1; v_iter++){ 00095 if (eta >= _eta_bins_ph.at(v_iter) && eta < _eta_bins_ph.at(v_iter+1)) 00096 break; 00097 } 00098 } 00099 else if (what==1) { 00100 for (v_iter=0; v_iter < (int)_eta_bins_jet.size()-1; v_iter++){ 00101 if (eta >= _eta_bins_jet.at(v_iter) && eta < _eta_bins_jet.at(v_iter+1)) 00102 break; 00103 } 00104 } 00105 else { 00106 for (v_iter=0; v_iter < (int)_eta_bins_areaoffset.size()-1; v_iter++){ 00107 if (eta >= _eta_bins_areaoffset.at(v_iter) && eta < _eta_bins_areaoffset.at(v_iter+1)) 00108 break; 00109 } 00110 } 00111 return v_iter; 00112 } 00113 00114 00115 // Perform the per-event analysis 00116 void analyze(const Event& event) { 00117 const double weight = event.weight(); 00118 00119 // Get the photon 00120 const FinalState& photonfs = applyProjection<FinalState>(event, "LeadingPhoton"); 00121 if (photonfs.particles().size() < 1) vetoEvent; 00122 00123 const FourMomentum photon = photonfs.particles().front().momentum(); 00124 double eta_P = photon.eta(); 00125 double phi_P = photon.phi(); 00126 00127 // Get the jet 00128 Jets jets = applyProjection<FastJets>(event, "Jets").jetsByPt(20.0*GeV); 00129 if (jets.size()==0) { 00130 vetoEvent; 00131 } 00132 FourMomentum leadingJet = jets[0].momentum(); 00133 00134 // Require jet separated from photon 00135 if (deltaR(eta_P, phi_P, leadingJet.eta(), leadingJet.phi())<1.0) { 00136 vetoEvent; 00137 } 00138 00139 // Veto if leading jet is outside plotted rapidity regions 00140 const double abs_y1 = fabs(leadingJet.rapidity()); 00141 if (abs_y1 > 4.4) { 00142 vetoEvent; 00143 } 00144 00145 // compute the median event energy density 00146 const unsigned int skipnhardjets = 0; 00147 _ptDensity.clear(); 00148 _sigma.clear(); 00149 _Njets.clear(); 00150 std::vector< std::vector<double> > ptDensities; 00151 std::vector<double> emptyVec; 00152 ptDensities.assign(_eta_bins_areaoffset.size()-1,emptyVec); 00153 00154 const fastjet::ClusterSequenceArea* clust_seq_area = applyProjection<FastJets>(event, "KtJetsD05").clusterSeqArea(); 00155 foreach (const fastjet::PseudoJet& jet, applyProjection<FastJets>(event, "KtJetsD05").pseudoJets(0.0*GeV)) { 00156 double eta = fabs(jet.eta()); 00157 double pt = fabs(jet.perp()); 00158 double area = clust_seq_area->area(jet); 00159 if (area > 10e-4 && fabs(eta)<_eta_bins_areaoffset[_eta_bins_areaoffset.size()-1]) { 00160 ptDensities.at(getEtaBin(fabs(eta),2)).push_back(pt/area); 00161 } 00162 } 00163 00164 for (int b = 0; b<(int)_eta_bins_areaoffset.size()-1; b++) { 00165 double median = 0.0; 00166 double sigma = 0.0; 00167 int Njets = 0; 00168 if (ptDensities[b].size() > skipnhardjets) { 00169 std::sort(ptDensities[b].begin(), ptDensities[b].end()); 00170 int nDens = ptDensities[b].size() - skipnhardjets; 00171 if ( nDens%2 == 0 ) { 00172 median = (ptDensities[b][nDens/2]+ptDensities[b][(nDens-2)/2])/2; 00173 } else { 00174 median = ptDensities[b][(nDens-1)/2]; 00175 } 00176 sigma = ptDensities[b][(int)(.15865*nDens)]; 00177 Njets = nDens; 00178 } 00179 _ptDensity.push_back(median); 00180 _sigma.push_back(sigma); 00181 _Njets.push_back(Njets); 00182 } 00183 00184 00185 // compute photon isolation 00186 00187 // std EtCone 00188 Particles fs = applyProjection<FinalState>(event, "JetFS").particles(); 00189 FourMomentum mom_in_EtCone; 00190 float iso_dR = 0.4; 00191 float cluster_eta_width = 0.25*7.0; 00192 float cluster_phi_width = (PI/128.)*5.0; 00193 foreach (const Particle& p, fs) { 00194 // check if it's in the cone of .4 00195 if (deltaR(eta_P, phi_P, p.eta(), p.momentum().phi()) >= iso_dR) continue; 00196 00197 // check if it's in the 5x7 central core 00198 if (fabs(eta_P-p.eta()) < cluster_eta_width*0.5 && 00199 fabs(phi_P-p.momentum().phi()) < cluster_phi_width*0.5) continue; 00200 00201 mom_in_EtCone += p.momentum(); 00202 } 00203 00204 // now figure out the correction (area*density) 00205 float EtCone_area = PI*iso_dR*iso_dR - cluster_eta_width*cluster_phi_width; 00206 float correction = _ptDensity[getEtaBin(eta_P,2)]*EtCone_area; 00207 00208 // require photon to be isolated 00209 if (mom_in_EtCone.Et()-correction > 4.0*GeV) vetoEvent; 00210 00211 int photon_jet_sign = sign( leadingJet.rapidity() * photon.rapidity() ); 00212 00213 // Fill histos 00214 float abs_jet_rapidity = fabs(leadingJet.rapidity()); 00215 float photon_pt = photon.pT()/GeV; 00216 float abs_photon_eta = fabs(photon.eta()); 00217 00218 if (abs_photon_eta<1.37) { 00219 00220 if (abs_jet_rapidity < 1.2) { 00221 00222 if (photon_jet_sign >= 1) { 00223 _h_phbarrel_jetcentral_SS->fill(photon_pt, weight); 00224 } else { 00225 _h_phbarrel_jetcentral_OS->fill(photon_pt, weight); 00226 } 00227 00228 } else if (abs_jet_rapidity < 2.8) { 00229 00230 if (photon_jet_sign >= 1) { 00231 _h_phbarrel_jetmedium_SS->fill(photon_pt, weight); 00232 } else { 00233 _h_phbarrel_jetmedium_OS->fill(photon_pt, weight); 00234 } 00235 00236 } else if (abs_jet_rapidity < 4.4) { 00237 00238 if (photon_jet_sign >= 1) { 00239 _h_phbarrel_jetforward_SS->fill(photon_pt, weight); 00240 } else { 00241 _h_phbarrel_jetforward_OS->fill(photon_pt, weight); 00242 } 00243 } 00244 00245 } 00246 } 00247 00248 00249 /// Normalise histograms etc., after the run 00250 void finalize() { 00251 scale(_h_phbarrel_jetcentral_SS, crossSection()/sumOfWeights()); 00252 scale(_h_phbarrel_jetcentral_OS, crossSection()/sumOfWeights()); 00253 scale(_h_phbarrel_jetmedium_SS, crossSection()/sumOfWeights()); 00254 scale(_h_phbarrel_jetmedium_OS, crossSection()/sumOfWeights()); 00255 scale(_h_phbarrel_jetforward_SS, crossSection()/sumOfWeights()); 00256 scale(_h_phbarrel_jetforward_OS, crossSection()/sumOfWeights()); 00257 } 00258 00259 00260 private: 00261 00262 Histo1DPtr _h_phbarrel_jetcentral_SS; 00263 Histo1DPtr _h_phbarrel_jetmedium_SS; 00264 Histo1DPtr _h_phbarrel_jetforward_SS; 00265 00266 Histo1DPtr _h_phbarrel_jetcentral_OS; 00267 Histo1DPtr _h_phbarrel_jetmedium_OS; 00268 Histo1DPtr _h_phbarrel_jetforward_OS; 00269 00270 fastjet::AreaDefinition* _area_def; 00271 00272 std::vector<float> _eta_bins_ph; 00273 std::vector<float> _eta_bins_jet; 00274 std::vector<float> _eta_bins_areaoffset; 00275 00276 std::vector<float> _ptDensity; 00277 std::vector<float> _sigma; 00278 std::vector<float> _Njets; 00279 }; 00280 00281 00282 00283 // The hook for the plugin system 00284 DECLARE_RIVET_PLUGIN(ATLAS_2012_I1093738); 00285 00286 } Generated on Thu Feb 6 2014 17:38:42 for The Rivet MC analysis system by ![]() |