ATLAS_2010_S8914702.cc
Go to the documentation of this file.00001
00002 #include <iostream>
00003 #include <sstream>
00004 #include <string>
00005
00006 #include "Rivet/Analysis.hh"
00007 #include "Rivet/RivetAIDA.hh"
00008 #include "Rivet/Tools/Logging.hh"
00009 #include "Rivet/Projections/FinalState.hh"
00010
00011 #include "Rivet/Projections/LeadingParticlesFinalState.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 class ATLAS_2010_S8914702 : public Analysis {
00026 public:
00027
00028
00029 ATLAS_2010_S8914702()
00030 : Analysis("ATLAS_2010_S8914702")
00031 {
00032 setNeedsCrossSection(true);
00033
00034 _eta_bins.push_back( 0.00);
00035 _eta_bins.push_back( 0.60);
00036 _eta_bins.push_back( 1.37);
00037 _eta_bins.push_back( 1.52);
00038 _eta_bins.push_back( 1.81);
00039
00040 _eta_bins_areaoffset.push_back(0.0);
00041 _eta_bins_areaoffset.push_back(1.5);
00042 _eta_bins_areaoffset.push_back(3.0);
00043 }
00044
00045
00046
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 LeadingParticlesFinalState photonfs(FinalState(-1.81, 1.81, 15.0*GeV));
00057 photonfs.addParticleId(PHOTON);
00058 addProjection(photonfs, "LeadingPhoton");
00059
00060 int hist_bin = 0;
00061 for (int i = 0; i < (int)_eta_bins.size()-1; ++i) {
00062 if (fabs(_eta_bins[i] - 1.37) < .0001) continue;
00063 _h_Et_photon[i] = bookHistogram1D(1, 1, hist_bin+1);
00064 hist_bin += 1;
00065 }
00066 }
00067
00068
00069 int getEtaBin(double eta_w, bool area_eta) const {
00070 double eta = fabs(eta_w);
00071 int v_iter = 0;
00072 if (!area_eta) {
00073 for (v_iter=0; v_iter < (int)_eta_bins.size()-1; ++v_iter) {
00074 if (eta >= _eta_bins.at(v_iter) && eta < _eta_bins.at(v_iter+1)) break;
00075 }
00076 } else {
00077 for (v_iter=0; v_iter < (int)_eta_bins_areaoffset.size()-1; ++v_iter) {
00078 if (eta >= _eta_bins_areaoffset.at(v_iter) && eta < _eta_bins_areaoffset.at(v_iter+1)) break;
00079 }
00080 }
00081 return v_iter;
00082 }
00083
00084
00085
00086 void analyze(const Event& event) {
00087 const double weight = event.weight();
00088
00089 ParticleVector photons = applyProjection<LeadingParticlesFinalState>(event, "LeadingPhoton").particles();
00090 if (photons.size() != 1) {
00091 vetoEvent;
00092 }
00093
00094 FourMomentum leadingPhoton = photons[0].momentum();
00095 double eta_P = leadingPhoton.eta();
00096 double phi_P = leadingPhoton.phi();
00097
00098 if(fabs(eta_P)>=1.37 && fabs(eta_P)<1.52){
00099 vetoEvent;
00100 }
00101
00102 int eta_bin = getEtaBin(eta_P,false);
00103
00104 ParticleVector fs = applyProjection<FinalState>(event, "FS").particles();
00105 FourMomentum mom_in_EtCone;
00106 foreach (const Particle& p, fs) {
00107
00108 if (deltaR(eta_P, phi_P, p.momentum().eta(), p.momentum().phi()) >= 0.4) continue;
00109
00110
00111 if (fabs(eta_P-p.momentum().eta()) < .025*7.0*0.5 &&
00112 fabs(phi_P-p.momentum().phi()) < (PI/128.)*5.0*0.5) continue;
00113 mom_in_EtCone += p.momentum();
00114 }
00115 MSG_DEBUG("Done with initial EtCone.");
00116
00117
00118 _ptDensity.clear();
00119 _sigma.clear();
00120 _Njets.clear();
00121
00122 std::vector< std::vector<double> > ptDensities;
00123 std::vector<double> emptyVec;
00124 ptDensities.assign(_eta_bins_areaoffset.size()-1,emptyVec);
00125
00126 const fastjet::ClusterSequenceArea* clust_seq_area = applyProjection<FastJets>(event, "KtJetsD05").clusterSeqArea();
00127 foreach (const fastjet::PseudoJet& jet, applyProjection<FastJets>(event, "KtJetsD05").pseudoJets(0.0*GeV)) {
00128
00129 const double eta = fabs(jet.eta());
00130 const double pt = fabs(jet.perp());
00131
00132
00133 double area = clust_seq_area->area(jet);
00134
00135 if (area > 10e-4 && fabs(eta)<_eta_bins_areaoffset[_eta_bins_areaoffset.size()-1]) {
00136 ptDensities.at(getEtaBin(fabs(eta),true)).push_back(pt/area);
00137 }
00138 }
00139
00140 for (int b = 0; b < (int)_eta_bins_areaoffset.size()-1; ++b) {
00141 double median = 0.0;
00142 double sigma = 0.0;
00143 int Njets = 0;
00144 if (ptDensities[b].size() > 0) {
00145 std::sort(ptDensities[b].begin(), ptDensities[b].end());
00146 int nDens = ptDensities[b].size();
00147 if (nDens % 2 == 0) {
00148 median = (ptDensities[b][nDens/2]+ptDensities[b][(nDens-2)/2])/2;
00149 } else {
00150 median = ptDensities[b][(nDens-1)/2];
00151 }
00152 sigma = ptDensities[b][(int)(.15865*nDens)];
00153 Njets = nDens;
00154 }
00155 _ptDensity.push_back(median);
00156 _sigma.push_back(sigma);
00157 _Njets.push_back(Njets);
00158 }
00159
00160
00161 float EtCone_area = PI*.4*.4 - (7.0*.025)*(5.0*PI/128.);
00162 float correction = _ptDensity[getEtaBin(eta_P,true)]*EtCone_area;
00163 MSG_DEBUG("Jet area correction done.");
00164
00165
00166
00167 if (mom_in_EtCone.Et() - correction > 4.0*GeV) {
00168 vetoEvent;
00169 }
00170 MSG_DEBUG("Passed isolation cut.");
00171
00172 _h_Et_photon[eta_bin]->fill(leadingPhoton.Et(), weight);
00173 }
00174
00175
00176
00177 void finalize() {
00178 for (int i = 0; i < (int)_eta_bins.size()-1; ++i) {
00179 if (fabs(_eta_bins[i] - 1.37) < .0001) continue;
00180 scale(_h_Et_photon[i], crossSection()/sumOfWeights());
00181 }
00182 }
00183
00184
00185 private:
00186
00187 AIDA::IHistogram1D *_h_Et_photon[6];
00188
00189 fastjet::AreaDefinition* _area_def;
00190
00191 std::vector<float> _eta_bins;
00192 std::vector<float> _eta_bins_areaoffset;
00193
00194 std::vector<float> _ptDensity;
00195 std::vector<float> _sigma;
00196 std::vector<float> _Njets;
00197 };
00198
00199
00200
00201
00202 AnalysisBuilder<ATLAS_2010_S8914702> plugin_ATLAS_2010_S8914702;
00203
00204 }