rivet is hosted by Hepforge, IPPP Durham
ATLAS_2016_I1457605.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/PromptFinalState.hh"
00005 #include "Rivet/Projections/LeadingParticlesFinalState.hh"
00006 #include "Rivet/Projections/FastJets.hh"
00007 
00008 namespace Rivet {
00009 
00010 
00011   /// Inclusive isolated prompt photon analysis with 2012 LHC data
00012   class ATLAS_2016_I1457605 : public Analysis {
00013   public:
00014 
00015     /// Constructor
00016     DEFAULT_RIVET_ANALYSIS_CTOR(ATLAS_2016_I1457605);
00017 
00018     /// Book histograms and initialise projections before the run
00019     void init() {
00020 
00021       FinalState fs;
00022       addProjection(fs, "FS");
00023 
00024       // Consider the final state jets for the energy density calculation
00025       FastJets fj(fs, FastJets::KT, 0.5);
00026       fj.useJetArea(new fastjet::AreaDefinition(fastjet::VoronoiAreaSpec()));
00027       addProjection(fj, "KtJetsD05");
00028 
00029       // Consider the leading pt photon with |eta| < 2.37 and pT > 25 GeV
00030       LeadingParticlesFinalState photonfs(PromptFinalState(FinalState(Cuts::abseta < 2.37 && Cuts::pT > 25*GeV)));
00031       photonfs.addParticleId(PID::PHOTON);
00032       addProjection(photonfs, "LeadingPhoton");
00033 
00034       // Book the dsigma/dEt (in eta bins) histograms
00035       for (size_t i = 0; i < _eta_bins.size() - 1; ++i) {
00036         if (fuzzyEquals(_eta_bins[i], 1.37)) continue; // skip this bin
00037         int offset = i > 2? 0 : 1;
00038         _h_Et_photon[i] = bookHisto1D(i + offset, 1, 1);
00039       }
00040 
00041     }
00042 
00043 
00044     /// Return eta bin for either dsigma/dET histogram (area_eta=false) or energy density correction (area_eta=true)
00045     size_t _getEtaBin(double eta_w, bool area_eta) const {
00046       const double eta = fabs(eta_w);
00047       if (!area_eta) {
00048         return binIndex(eta, _eta_bins);
00049       } else {
00050         return binIndex(eta, _eta_bins_areaoffset);
00051       }
00052     }
00053 
00054 
00055     /// Perform the per-event analysis
00056     void analyze(const Event& event) {
00057       // Retrieve leading photon
00058       Particles photons = applyProjection<LeadingParticlesFinalState>(event, "LeadingPhoton").particles();
00059       if (photons.size() < 1)  vetoEvent;
00060       const Particle& leadingPhoton = photons[0];
00061 
00062       // Veto events with photon in ECAL crack
00063       if (inRange(leadingPhoton.abseta(), 1.37, 1.56)) vetoEvent;
00064 
00065       // Compute isolation energy in cone of radius .4 around photon (all particles)
00066       FourMomentum mom_in_EtCone;
00067       Particles fs = applyProjection<FinalState>(event, "FS").particles();
00068       foreach (const Particle& p, fs) {
00069         // Check if it's outside the cone of 0.4
00070         if (deltaR(leadingPhoton, p) >= 0.4) continue;
00071         // Except muons or neutrinos
00072         if (PID::isNeutrino(p.abspid()) || p.abspid() == PID::MUON) continue;
00073         // Increment isolation energy
00074         mom_in_EtCone += p.momentum();
00075       }
00076       // Remove the photon energy from the isolation
00077       mom_in_EtCone -= leadingPhoton.momentum();
00078 
00079       // Get the area-filtered jet inputs for computing median energy density, etc.
00080       vector<double> ptDensity;
00081       vector< vector<double> > ptDensities(_eta_bins_areaoffset.size()-1);
00082       const FastJets& fast_jets = applyProjection<FastJets>(event, "KtJetsD05");
00083       const auto clust_seq_area = fast_jets.clusterSeqArea();
00084       foreach (const Jet& jet, fast_jets.jets()) {
00085         const double area = clust_seq_area->area(jet);
00086         if (area > 1e-3 && jet.abseta() < _eta_bins_areaoffset.back())
00087           ptDensities.at( _getEtaBin(jet.abseta(), true) ) += jet.pT()/area;
00088       }
00089       // Compute the median energy density, etc.
00090       for (size_t b = 0; b < _eta_bins_areaoffset.size()-1; ++b) {
00091         const int njets = ptDensities[b].size();
00092         ptDensity += (njets > 0) ? median(ptDensities[b]) : 0.0;
00093       }
00094       // Compute the isolation energy correction (cone area*energy density)
00095       const double etCone_area = PI * sqr(0.4);
00096       const double correction = ptDensity[_getEtaBin(leadingPhoton.abseta(), true)] * etCone_area;
00097 
00098       // Apply isolation cut on area-corrected value
00099       // cut is Etiso < 4.8GeV + 4.2E-03 * Et_gamma.
00100       if (mom_in_EtCone.Et() - correction > 4.8*GeV + 0.0042*leadingPhoton.Et()) vetoEvent;
00101 
00102       // Fill histograms
00103       const size_t eta_bin = _getEtaBin(leadingPhoton.abseta(), false);
00104       _h_Et_photon[eta_bin]->fill(leadingPhoton.Et(), event.weight());
00105     }
00106 
00107 
00108     /// Normalise histograms etc., after the run
00109     void finalize() {
00110       double sf = crossSection() / (picobarn * sumOfWeights());
00111       for (size_t i = 0; i < _eta_bins.size()-1; ++i) {
00112         if (fuzzyEquals(_eta_bins[i], 1.37)) continue;
00113         scale(_h_Et_photon[i], sf);
00114       }
00115     }
00116 
00117 
00118   private:
00119 
00120     Histo1DPtr _h_Et_photon[5];
00121 
00122     const vector<double> _eta_bins = {0.00, 0.60, 1.37, 1.56, 1.81, 2.37 };
00123     const vector<double> _eta_bins_areaoffset = {0.0, 1.5, 3.0};
00124 
00125   };
00126 
00127 
00128   DECLARE_RIVET_PLUGIN(ATLAS_2016_I1457605);
00129 
00130 }