rivet is hosted by Hepforge, IPPP Durham
ATLAS_2012_I1093738.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/RivetYODA.hh"
00008 #include "Rivet/Tools/Logging.hh"
00009 #include "Rivet/Projections/FinalState.hh"
00010 #include "Rivet/Projections/LeadingParticlesFinalState.hh"
00011 #include "Rivet/Projections/VetoedFinalState.hh"
00012 #include "Rivet/Jet.hh"
00013 #include "Rivet/Projections/FastJets.hh"
00014 
00015 #include "fastjet/internal/base.hh"
00016 #include "fastjet/JetDefinition.hh"
00017 #include "fastjet/AreaDefinition.hh"
00018 #include "fastjet/ClusterSequence.hh"
00019 #include "fastjet/ClusterSequenceArea.hh"
00020 #include "fastjet/PseudoJet.hh"
00021 
00022 #define MYDEBUG if(false) getLog() << __LINE__
00023 
00024 
00025 namespace Rivet {
00026 
00027   /// @brief Measurement of isolated gamma + jet + X differential cross-sections
00028   ///
00029   /// Inclusive isolated gamma + jet cross-sections, differential in pT(gamma), for
00030   /// various photon and jet rapidity configurations.
00031   ///
00032   /// @author Giovanni Marchiori
00033 
00034   class ATLAS_2012_I1093738 : public Analysis {
00035   public:
00036 
00037     // Constructor
00038     ATLAS_2012_I1093738()
00039       : Analysis("ATLAS_2012_I1093738")
00040     {
00041       _eta_bins_ph.push_back(0.0);
00042       _eta_bins_ph.push_back(1.37);
00043       _eta_bins_ph.push_back(1.52);
00044       _eta_bins_ph.push_back(2.37);
00045 
00046       _eta_bins_jet.push_back(0.0);
00047       _eta_bins_jet.push_back(1.2);
00048       _eta_bins_jet.push_back(2.8);
00049       _eta_bins_jet.push_back(4.4);
00050 
00051       _eta_bins_areaoffset.push_back(0.0);
00052       _eta_bins_areaoffset.push_back(1.5);
00053       _eta_bins_areaoffset.push_back(3.0);
00054     }
00055 
00056   public:
00057 
00058     // Book histograms and initialise projections before the run
00059     void init() {
00060       MYDEBUG << "Entering init." << std::endl;
00061 
00062       // Final state
00063       FinalState fs;
00064       addProjection(fs, "FS");
00065 
00066       // Voronoi eta-phi tassellation with KT jets, for ambient energy density calculation
00067       FastJets fj(fs, FastJets::KT, 0.5);
00068       _area_def = new fastjet::AreaDefinition(fastjet::VoronoiAreaSpec());
00069       fj.useJetArea(_area_def);
00070       addProjection(fj, "KtJetsD05");
00071 
00072       // Leading photon
00073       LeadingParticlesFinalState photonfs(FinalState(-1.37, 1.37, 25.0*GeV));
00074       photonfs.addParticleId(PHOTON);
00075       addProjection(photonfs, "LeadingPhoton");
00076 
00077       // FS excluding the leading photon
00078       VetoedFinalState vfs(fs);
00079       vfs.addVetoOnThisFinalState(photonfs);
00080       addProjection(vfs, "JetFS");
00081 
00082       // Jets
00083       FastJets jetpro(vfs, FastJets::ANTIKT, 0.4);
00084       //FastJets jetpro(vfs, FastJets::ANTIKT, 0.6);
00085       jetpro.useInvisibles();
00086       addProjection(jetpro, "Jets");
00087 
00088 
00089       MYDEBUG << "\t... Booking Histograms " << std::endl;
00090       // use autobook option, histo d??-x??-y??
00091       _h_phbarrel_jetcentral_SS = bookHisto1D(1, 1, 1);
00092       _h_phbarrel_jetmedium_SS  = bookHisto1D(2, 1, 1);
00093       _h_phbarrel_jetforward_SS = bookHisto1D(3, 1, 1);
00094 
00095       _h_phbarrel_jetcentral_OS = bookHisto1D(4, 1, 1);
00096       _h_phbarrel_jetmedium_OS  = bookHisto1D(5, 1, 1);
00097       _h_phbarrel_jetforward_OS = bookHisto1D(6, 1, 1);
00098 
00099       MYDEBUG << "Exiting init." << std::endl;
00100     }
00101 
00102 
00103     int getEtaBin(double eta_w, int what) const {
00104       double eta = fabs(eta_w);
00105 
00106       int v_iter=0;
00107       if (what==0) {
00108         for(v_iter=0; v_iter < (int)_eta_bins_ph.size()-1; v_iter++){
00109           if(eta >= _eta_bins_ph.at(v_iter) && eta < _eta_bins_ph.at(v_iter+1))
00110             break;
00111         }
00112       }
00113       else if (what==1) {
00114         for (v_iter=0; v_iter < (int)_eta_bins_jet.size()-1; v_iter++){
00115           if(eta >= _eta_bins_jet.at(v_iter) && eta < _eta_bins_jet.at(v_iter+1))
00116             break;
00117         }
00118       }
00119       else {
00120         for(v_iter=0; v_iter < (int)_eta_bins_areaoffset.size()-1; v_iter++){
00121           if(eta >= _eta_bins_areaoffset.at(v_iter) && eta < _eta_bins_areaoffset.at(v_iter+1))
00122             break;
00123         }
00124       }
00125 
00126       return v_iter;
00127     }
00128 
00129 
00130     // Perform the per-event analysis
00131     void analyze(const Event& event) {
00132 
00133       MYDEBUG << "Entering Analyze." << std::endl;
00134 
00135       const double weight = event.weight();
00136 
00137       MYDEBUG << "...First projections." << std::endl;
00138 
00139 
00140       // Get the photon
00141       const FinalState& photonfs = applyProjection<FinalState>(event, "LeadingPhoton");
00142       if (photonfs.particles().size() < 1) {
00143         MYDEBUG << "...Going to veto event(1)." << std::endl;
00144         vetoEvent;
00145       }
00146       MYDEBUG << "...Didn't veto event(1)." << std::endl;
00147 
00148       const FourMomentum photon = photonfs.particles().front().momentum();
00149       double eta_P = photon.eta();
00150       double phi_P = photon.phi();
00151 
00152       // Get the jet
00153       Jets jets = applyProjection<FastJets>(event, "Jets").jetsByPt(20.0*GeV);
00154       if (jets.size()==0) {
00155         MYDEBUG << "...Going to veto event(3)." << std::endl;
00156         vetoEvent;
00157       }
00158       MYDEBUG << "...Didn't veto event(3)." << std::endl;
00159       FourMomentum leadingJet = jets[0].momentum();
00160 
00161       // Require jet separated from photon
00162       if (deltaR(eta_P, phi_P, leadingJet.eta(), leadingJet.phi())<1.0) {
00163         MYDEBUG << "...Going to veto event(4)." << std::endl;
00164         vetoEvent;
00165       }
00166       MYDEBUG << "...Didn't veto event(4)." << std::endl;
00167 
00168       // Veto if leading jet is outside plotted rapidity regions
00169       const double abs_y1 = fabs(leadingJet.rapidity());
00170       if (abs_y1 > 4.4) {
00171         MYDEBUG << "...Going to veto event(5)." << std::endl;
00172         vetoEvent;
00173       }
00174       MYDEBUG << "...Didn't veto event(5)." << std::endl;
00175 
00176 
00177       // compute the median event energy density
00178       const unsigned int skipnhardjets = 0;
00179       _ptDensity.clear();
00180       _sigma.clear();
00181       _Njets.clear();
00182       std::vector< std::vector<double> > ptDensities;
00183       std::vector<double> emptyVec;
00184       ptDensities.assign(_eta_bins_areaoffset.size()-1,emptyVec);
00185 
00186       const fastjet::ClusterSequenceArea* clust_seq_area = applyProjection<FastJets>(event, "KtJetsD05").clusterSeqArea();
00187       foreach (const fastjet::PseudoJet& jet, applyProjection<FastJets>(event, "KtJetsD05").pseudoJets(0.0*GeV)) {
00188         double eta = fabs(jet.eta());
00189         double pt = fabs(jet.perp());
00190 
00191         /// get the cluster sequence
00192         double area = clust_seq_area->area(jet);
00193 
00194         if(area > 10e-4 && fabs(eta)<_eta_bins_areaoffset[_eta_bins_areaoffset.size()-1]){
00195           ptDensities.at(getEtaBin(fabs(eta),2)).push_back(pt/area);
00196         }
00197       }
00198 
00199       for(int b=0; b<(int)_eta_bins_areaoffset.size()-1; b++){
00200         double median = 0.0;
00201         double sigma = 0.0;
00202         int Njets = 0;
00203         if(ptDensities[b].size() > skipnhardjets)
00204           {
00205             std::sort(ptDensities[b].begin(), ptDensities[b].end());
00206             int nDens = ptDensities[b].size() - skipnhardjets;
00207             if( nDens%2 == 0 )
00208               median = (ptDensities[b][nDens/2]+ptDensities[b][(nDens-2)/2])/2;
00209             else
00210               median = ptDensities[b][(nDens-1)/2];
00211             sigma = ptDensities[b][(int)(.15865*nDens)];
00212             Njets = nDens;
00213           }
00214         _ptDensity.push_back(median);
00215         _sigma.push_back(sigma);
00216         _Njets.push_back(Njets);
00217       }
00218 
00219 
00220       // compute photon isolation
00221 
00222       // std EtCone
00223       ParticleVector fs = applyProjection<FinalState>(event, "JetFS").particles();
00224       FourMomentum mom_in_EtCone;
00225       float iso_dR = 0.4;
00226       float cluster_eta_width = 0.25*7.0;
00227       float cluster_phi_width = (PI/128.)*5.0;
00228       foreach (const Particle& p, fs) {
00229         // check if it's in the cone of .4
00230         if (deltaR(eta_P, phi_P, p.momentum().eta(), p.momentum().phi()) >= iso_dR) continue;
00231 
00232         // check if it's in the 5x7 central core
00233         if (fabs(eta_P-p.momentum().eta()) < cluster_eta_width*0.5 &&
00234             fabs(phi_P-p.momentum().phi()) < cluster_phi_width*0.5) continue;
00235 
00236         mom_in_EtCone += p.momentum();
00237       }
00238       MYDEBUG << "...Done with initial EtCone." << std::endl;
00239 
00240       // now figure out the correction (area*density)
00241       float EtCone_area = PI*iso_dR*iso_dR - cluster_eta_width*cluster_phi_width;
00242       float correction = _ptDensity[getEtaBin(eta_P,2)]*EtCone_area;
00243       MYDEBUG << "...Done with jet-area correction." << std::endl;
00244 
00245       // require photon to be isolated
00246       if(mom_in_EtCone.Et()-correction > 4.0*GeV){
00247         MYDEBUG << "...Going to veto event(6)." << std::endl;
00248         vetoEvent;
00249       }
00250       MYDEBUG << "...Didn't veto event(6)." << std::endl;
00251 
00252 
00253       int photon_jet_sign = sign( leadingJet.rapidity() * photon.rapidity() );
00254 
00255       // Fill histos
00256       float abs_jet_rapidity = fabs(leadingJet.rapidity());
00257       float photon_pt = photon.pT()/GeV;
00258       float abs_photon_eta = fabs(photon.eta());
00259 
00260       if (abs_photon_eta<1.37) {
00261 
00262         if (abs_jet_rapidity < 1.2) {
00263 
00264           if (photon_jet_sign >= 1) {
00265             _h_phbarrel_jetcentral_SS->fill(photon_pt, weight);
00266           } else {
00267             _h_phbarrel_jetcentral_OS->fill(photon_pt, weight);
00268           }
00269 
00270         } else if (abs_jet_rapidity < 2.8) {
00271 
00272           if (photon_jet_sign >= 1) {
00273             _h_phbarrel_jetmedium_SS->fill(photon_pt, weight);
00274           } else {
00275             _h_phbarrel_jetmedium_OS->fill(photon_pt, weight);
00276           }
00277 
00278         } else if (abs_jet_rapidity < 4.4) {
00279 
00280           if (photon_jet_sign >= 1) {
00281             _h_phbarrel_jetforward_SS->fill(photon_pt, weight);
00282           } else {
00283             _h_phbarrel_jetforward_OS->fill(photon_pt, weight);
00284           }
00285         }
00286 
00287       }
00288 
00289       MYDEBUG << "...Done with analyze." << std::endl;
00290     }
00291 
00292 
00293     /// Normalise histograms etc., after the run
00294     void finalize() {
00295 
00296       MYDEBUG << "In finalize." << std::endl;
00297 
00298       scale(_h_phbarrel_jetcentral_SS, crossSection()/sumOfWeights());
00299       scale(_h_phbarrel_jetcentral_OS, crossSection()/sumOfWeights());
00300       scale(_h_phbarrel_jetmedium_SS, crossSection()/sumOfWeights());
00301       scale(_h_phbarrel_jetmedium_OS, crossSection()/sumOfWeights());
00302       scale(_h_phbarrel_jetforward_SS, crossSection()/sumOfWeights());
00303       scale(_h_phbarrel_jetforward_OS, crossSection()/sumOfWeights());
00304 
00305       MYDEBUG << "Done with finalize." << std::endl;
00306     }
00307 
00308   private:
00309 
00310     Histo1DPtr _h_phbarrel_jetcentral_SS;
00311     Histo1DPtr _h_phbarrel_jetmedium_SS;
00312     Histo1DPtr _h_phbarrel_jetforward_SS;
00313 
00314     Histo1DPtr _h_phbarrel_jetcentral_OS;
00315     Histo1DPtr _h_phbarrel_jetmedium_OS;
00316     Histo1DPtr _h_phbarrel_jetforward_OS;
00317 
00318     fastjet::AreaDefinition* _area_def;
00319 
00320     std::vector<float> _eta_bins_ph;
00321     std::vector<float> _eta_bins_jet;
00322     std::vector<float> _eta_bins_areaoffset;
00323 
00324     std::vector<float> _ptDensity;
00325     std::vector<float> _sigma;
00326     std::vector<float> _Njets;
00327   };
00328 
00329   // The hook for the plugin system
00330   DECLARE_RIVET_PLUGIN(ATLAS_2012_I1093738);
00331 
00332 }