rivet is hosted by Hepforge, IPPP Durham
ATLAS_2011_S9120807.cc
Go to the documentation of this file.
00001 // -*- C++ -*-
00002 #include "Rivet/Analysis.hh"
00003 #include "Rivet/Projections/FinalState.hh"
00004 #include "Rivet/Projections/IdentifiedFinalState.hh"
00005 #include "Rivet/Projections/FastJets.hh"
00006 
00007 namespace Rivet {
00008 
00009 
00010   /// @brief Measurement of isolated diphoton + X differential cross-sections
00011   ///
00012   /// Inclusive isolated gamma gamma cross-sections, differential in M(gg), pT(gg), dphi(gg)
00013   ///
00014   /// @author Giovanni Marchiori
00015   class ATLAS_2011_S9120807 : public Analysis {
00016   public:
00017 
00018     /// Constructor
00019     ATLAS_2011_S9120807()
00020       : Analysis("ATLAS_2011_S9120807")
00021     {    }
00022 
00023 
00024     /// Book histograms and initialise projections before the run
00025     void init() {
00026       FinalState fs;
00027       declare(fs, "FS");
00028 
00029       FastJets fj(fs, FastJets::KT, 0.5);
00030       fj.useJetArea(new fastjet::AreaDefinition(fastjet::VoronoiAreaSpec()));
00031       declare(fj, "KtJetsD05");
00032 
00033       IdentifiedFinalState photonfs(Cuts::abseta < 2.37 && Cuts::pT > 16*GeV);
00034       photonfs.acceptId(PID::PHOTON);
00035       declare(photonfs, "Photon");
00036 
00037       _h_M    = bookHisto1D(1, 1, 1);
00038       _h_pT   = bookHisto1D(2, 1, 1);
00039       _h_dPhi = bookHisto1D(3, 1, 1);
00040     }
00041 
00042 
00043     size_t getEtaBin(double eta) const {
00044       const double aeta = fabs(eta);
00045       return binIndex(aeta, _eta_bins_areaoffset);
00046     }
00047 
00048 
00049     /// Perform the per-event analysis
00050     void analyze(const Event& event) {
00051       // Require at least 2 photons in final state
00052       Particles photons = apply<IdentifiedFinalState>(event, "Photon").particlesByPt();
00053       if (photons.size() < 2) vetoEvent;
00054 
00055       // Compute jet pT densities
00056       vector<double> ptDensity;
00057       vector< vector<double> > ptDensities(_eta_bins_areaoffset.size()-1);
00058       const shared_ptr<fastjet::ClusterSequenceArea> clust_seq_area = apply<FastJets>(event, "KtJetsD05").clusterSeqArea();
00059       for (const Jet& jet : apply<FastJets>(event, "KtJetsD05").jets()) {
00060         const double area = clust_seq_area->area(jet); // .pseudojet() called implicitly
00061         /// @todo Should be 1e-4?
00062         if (area > 10e-4 && jet.abseta() < _eta_bins_areaoffset.back()) {
00063           ptDensities.at(getEtaBin(jet.abseta())).push_back(jet.pT()/area);
00064         }
00065       }
00066 
00067       // Compute the median energy density
00068       for (size_t b = 0; b < _eta_bins_areaoffset.size()-1; ++b) {
00069         ptDensity += ptDensities[b].empty() ? 0 : median(ptDensities[b]);
00070       }
00071 
00072       // Loop over photons and fill vector of isolated ones
00073       Particles isolated_photons;
00074       for (const Particle& photon : photons) {
00075 
00076         // Remove photons in crack
00077         if (inRange(photon.abseta(), 1.37, 1.52)) continue;
00078 
00079         // Standard ET cone isolation
00080         const Particles& fs = apply<FinalState>(event, "FS").particles();
00081         FourMomentum mom_in_EtCone;
00082         for (const Particle& p : fs) {
00083           // Check if it's in the cone of .4
00084           if (deltaR(photon, p) >= 0.4) continue;
00085           // Veto if it's in the 5x7 central core
00086           if (fabs(deltaEta(photon, p)) < 0.025*5.0*0.5 &&
00087               fabs(deltaPhi(photon, p)) < (M_PI/128.)*7.0*0.5) continue;
00088           // Increment isolation cone ET sum
00089           mom_in_EtCone += p.momentum();
00090         }
00091 
00092         // Now figure out the correction (area*density)
00093         const double ETCONE_AREA = M_PI*.4*.4 - (7.0*.025)*(5.0*M_PI/128.);
00094         const double correction = ptDensity[getEtaBin(photon.abseta())] * ETCONE_AREA;
00095 
00096         // Shouldn't need to subtract photon
00097         // NB. Using expected cut at hadron/particle level, not cut at reco level
00098         if (mom_in_EtCone.Et() - correction > 4.0*GeV) continue;
00099 
00100         // Add to list of isolated photons
00101         isolated_photons.push_back(photon);
00102       }
00103 
00104       // Require at least two isolated photons
00105       if (isolated_photons.size() < 2) vetoEvent;
00106 
00107       // Select leading pT pair
00108       std::sort(isolated_photons.begin(), isolated_photons.end(), cmpMomByPt);
00109       const FourMomentum y1 = isolated_photons[0].momentum();
00110       const FourMomentum y2 = isolated_photons[1].momentum();
00111 
00112       // Require the two photons to be separated (dR>0.4)
00113       if (deltaR(y1, y2) < 0.4) vetoEvent;
00114 
00115       const FourMomentum yy = y1 + y2;
00116       const double Myy = yy.mass()/GeV;
00117       const double pTyy = yy.pT()/GeV;
00118       const double dPhiyy = deltaPhi(y1.phi(), y2.phi());
00119 
00120       const double weight = event.weight();
00121       _h_M->fill(Myy, weight);
00122       _h_pT->fill(pTyy, weight);
00123       _h_dPhi->fill(dPhiyy, weight);
00124     }
00125 
00126 
00127     /// Normalise histograms etc., after the run
00128     void finalize() {
00129       scale(_h_M, crossSection()/sumOfWeights());
00130       scale(_h_pT, crossSection()/sumOfWeights());
00131       scale(_h_dPhi, crossSection()/sumOfWeights());
00132     }
00133 
00134 
00135   private:
00136 
00137     Histo1DPtr _h_M, _h_pT, _h_dPhi;
00138     const vector<double> _eta_bins_areaoffset = {0.0, 1.5, 3.0};
00139 
00140   };
00141 
00142 
00143   // The hook for the plugin system
00144   DECLARE_RIVET_PLUGIN(ATLAS_2011_S9120807);
00145 
00146 
00147 }