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