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(-2.37, 2.37, 16.0*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     int getEtaBin(double eta_w) const {
00065       double eta = fabs(eta_w);
00066 
00067       int v_iter=0;
00068       for (v_iter=0; v_iter < (int)_eta_bins_areaoffset.size()-1; v_iter++) {
00069         if (inRange(eta, _eta_bins_areaoffset[v_iter], _eta_bins_areaoffset[v_iter+1]))
00070           break;
00071       }
00072       return v_iter;
00073     }
00074 
00075 
00076     /// Perform the per-event analysis
00077     void analyze(const Event& event) {
00078 
00079       const double weight = event.weight();
00080 
00081       ///
00082       /// require at least 2 photons in final state
00083       ///
00084       Particles photons = applyProjection<IdentifiedFinalState>(event, "Photon").particlesByPt();
00085       if (photons.size() < 2) {
00086         vetoEvent;
00087       }
00088 
00089       ///
00090       /// compute the median energy density
00091       ///
00092       std::vector<double> _ptDensity;
00093       std::vector< std::vector<double> > ptDensities;
00094       std::vector<double> emptyVec;
00095       ptDensities.assign(_eta_bins_areaoffset.size()-1, emptyVec);
00096 
00097       const fastjet::ClusterSequenceArea* clust_seq_area = applyProjection<FastJets>(event, "KtJetsD05").clusterSeqArea();
00098       foreach (const fastjet::PseudoJet& jet, applyProjection<FastJets>(event, "KtJetsD05").pseudoJets(0.0*GeV)) {
00099         double eta = fabs(jet.eta());
00100         double pt = fabs(jet.perp());
00101 
00102         /// get the cluster sequence
00103         double area = clust_seq_area->area(jet);
00104 
00105         if(area > 10e-4 && fabs(eta) < _eta_bins_areaoffset[_eta_bins_areaoffset.size()-1]) {
00106           ptDensities.at(getEtaBin(fabs(eta))).push_back(pt/area);
00107         }
00108       }
00109 
00110       for(int b=0; b<(int)_eta_bins_areaoffset.size()-1; b++) {
00111         double median = 0.0;
00112         if(ptDensities[b].size() > 0) {
00113             std::sort(ptDensities[b].begin(), ptDensities[b].end());
00114             int nDens = ptDensities[b].size();
00115             if( nDens%2 == 0 )
00116               median = (ptDensities[b][nDens/2] + ptDensities[b][(nDens-2)/2]) / 2;
00117             else
00118               median = ptDensities[b][(nDens-1)/2];
00119           }
00120         _ptDensity.push_back(median);
00121       }
00122 
00123 
00124       ///
00125       /// Loop over photons and fill vector of isolated ones
00126       ///
00127       Particles isolated_photons;
00128       foreach (const Particle& photon, photons) {
00129 
00130         ///
00131         /// remove photons in crack
00132         ///
00133         double eta_P = photon.eta();
00134         if (fabs(eta_P)>=1.37 && fabs(eta_P)<1.52) continue;
00135 
00136         double phi_P = photon.momentum().phi();
00137 
00138         ///
00139         /// compute isolation
00140         ///
00141 
00142         /// std EtCone
00143         Particles fs = applyProjection<FinalState>(event, "FS").particles();
00144         FourMomentum mom_in_EtCone;
00145         foreach (const Particle& p, fs) {
00146           /// check if it's in the cone of .4
00147           if (deltaR(eta_P, phi_P, p.eta(), p.momentum().phi()) >= 0.4) continue;
00148 
00149           /// check if it's in the 5x7 central core
00150           if (fabs(eta_P-p.eta()) < .025*7.0*0.5 &&
00151               fabs(phi_P-p.momentum().phi()) < (M_PI/128.)*5.0*0.5) continue;
00152 
00153           mom_in_EtCone += p.momentum();
00154         }
00155 
00156         /// now figure out the correction (area*density)
00157         double EtCone_area = M_PI*.4*.4 - (7.0*.025)*(5.0*M_PI/128.);
00158         double correction = _ptDensity[getEtaBin(eta_P)]*EtCone_area;
00159 
00160         /// shouldn't need to subtract photon
00161         /// note: using expected cut at hadron/particle level, not cut at reco level
00162         if (mom_in_EtCone.Et()-correction > 4.0*GeV) {
00163           continue;
00164         }
00165 
00166         /// add photon to list of isolated ones
00167         isolated_photons.push_back(photon);
00168       }
00169 
00170       ///
00171       /// require at least two isolated photons
00172       ///
00173       if (isolated_photons.size() < 2) {
00174         vetoEvent;
00175       }
00176 
00177       ///
00178       /// select leading pT pair
00179       ///
00180       std::sort(isolated_photons.begin(), isolated_photons.end(), cmpMomByPt);
00181       FourMomentum y1=isolated_photons[0].momentum();
00182       FourMomentum y2=isolated_photons[1].momentum();
00183 
00184       ///
00185       /// require the two photons to be separated (dR>0.4)
00186       ///
00187       if (deltaR(y1, y2)<0.4) {
00188         vetoEvent;
00189       }
00190 
00191       FourMomentum yy = y1+y2;
00192       double Myy = yy.mass()/GeV;
00193       double pTyy = yy.pT()/GeV;
00194       double dPhiyy = deltaPhi(y1.phi(), y2.phi());
00195 
00196       _h_M->fill(Myy, weight);
00197       _h_pT->fill(pTyy, weight);
00198       _h_dPhi->fill(dPhiyy, weight);
00199     }
00200 
00201 
00202     /// Normalise histograms etc., after the run
00203     void finalize() {
00204       scale(_h_M, crossSection()/sumOfWeights());
00205       scale(_h_pT, crossSection()/sumOfWeights());
00206       scale(_h_dPhi, crossSection()/sumOfWeights());
00207     }
00208 
00209   private:
00210 
00211     Histo1DPtr _h_M;
00212     Histo1DPtr _h_pT;
00213     Histo1DPtr _h_dPhi;
00214 
00215     fastjet::AreaDefinition* _area_def;
00216 
00217     std::vector<double> _eta_bins_areaoffset;
00218   };
00219 
00220 
00221   // The hook for the plugin system
00222   DECLARE_RIVET_PLUGIN(ATLAS_2011_S9120807);
00223 
00224 }