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