D0_2008_S7719523.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 #include "Rivet/Analysis.hh"
00003 #include "Rivet/Tools/Logging.hh"
00004 #include "Rivet/Projections/FinalState.hh"
00005 #include "Rivet/Projections/LeadingParticlesFinalState.hh"
00006 #include "Rivet/Projections/VetoedFinalState.hh"
00007 #include "Rivet/Projections/FastJets.hh"
00008 #include "Rivet/RivetAIDA.hh"
00009 
00010 namespace Rivet {
00011 
00012 
00013   /// @brief Measurement of isolated gamma + jet + X differential cross-sections
00014   ///
00015   /// Inclusive isolated gamma + jet cross-sections, differential in pT(gamma), for
00016   /// various photon and jet rapidity bins.
00017   ///
00018   /// @author Andy Buckley
00019   /// @author Gavin Hesketh
00020   class D0_2008_S7719523 : public Analysis {
00021   public:
00022 
00023     /// @name Constructors etc.
00024     //@{
00025 
00026     /// Constructor
00027     D0_2008_S7719523()
00028       : Analysis("D0_2008_S7719523")
00029     {    }
00030 
00031     //@}
00032 
00033 
00034     /// @name Analysis methods
00035     //@{
00036 
00037     /// Set up projections and book histograms
00038     void init() {
00039       // General FS
00040       FinalState fs;
00041       addProjection(fs, "FS");
00042 
00043       // Get leading photon
00044       LeadingParticlesFinalState photonfs(FinalState(-1.0, 1.0, 30.0*GeV));
00045       photonfs.addParticleId(PHOTON);
00046       addProjection(photonfs, "LeadingPhoton");
00047 
00048       // FS excluding the leading photon
00049       VetoedFinalState vfs(fs);
00050       vfs.addVetoOnThisFinalState(photonfs);
00051       addProjection(vfs, "JetFS");
00052 
00053       // Jets
00054       FastJets jetpro(vfs, FastJets::D0ILCONE, 0.7);
00055       addProjection(jetpro, "Jets");
00056 
00057       // Histograms
00058       _h_central_same_cross_section = bookHistogram1D(1, 1, 1);
00059       _h_central_opp_cross_section  = bookHistogram1D(2, 1, 1);
00060       _h_forward_same_cross_section = bookHistogram1D(3, 1, 1);
00061       _h_forward_opp_cross_section  = bookHistogram1D(4, 1, 1);
00062     }
00063 
00064 
00065 
00066     /// Do the analysis
00067     void analyze(const Event& event) {
00068       const double weight = event.weight();
00069 
00070       // Get the photon
00071       const FinalState& photonfs = applyProjection<FinalState>(event, "LeadingPhoton");
00072       if (photonfs.particles().size() != 1) {
00073         vetoEvent;
00074       }
00075       const FourMomentum photon = photonfs.particles().front().momentum();
00076 
00077       // Isolate photon by ensuring that a 0.4 cone around it contains less than 7% of the photon's energy
00078       double egamma = photon.E();
00079       double eta_P = photon.pseudorapidity();
00080       double phi_P = photon.azimuthalAngle();
00081       double econe = 0.0;
00082       foreach (const Particle& p, applyProjection<FinalState>(event, "JetFS").particles()) {
00083         if (deltaR(eta_P, phi_P,
00084                    p.momentum().pseudorapidity(), p.momentum().azimuthalAngle()) < 0.4) {
00085           econe += p.momentum().E();
00086           // Veto as soon as E_cone gets larger
00087           if (econe/egamma > 0.07) {
00088             MSG_DEBUG("Vetoing event because photon is insufficiently isolated");
00089             vetoEvent;
00090           }
00091         }
00092       }
00093 
00094       Jets jets = applyProjection<FastJets>(event, "Jets").jetsByPt(15.0*GeV);
00095       if (jets.size()==0) {
00096         vetoEvent;
00097       }
00098       FourMomentum leadingJet = jets[0].momentum();
00099       if (deltaR(eta_P, phi_P, leadingJet.eta(), leadingJet.phi())<0.7) {
00100         vetoEvent;
00101       }
00102 
00103       int photon_jet_sign = sign( leadingJet.rapidity() * photon.rapidity() );
00104 
00105       // Veto if leading jet is outside plotted rapidity regions
00106       const double abs_y1 = fabs(leadingJet.rapidity());
00107       if (inRange(abs_y1, 0.8, 1.5) || abs_y1 > 2.5) {
00108         MSG_DEBUG("Leading jet falls outside acceptance range; |y1| = "
00109                   << abs_y1);
00110         vetoEvent;
00111       }
00112 
00113       // Fill histos
00114       if (fabs(leadingJet.rapidity()) < 0.8) {
00115         if (photon_jet_sign >= 1) {
00116           _h_central_same_cross_section->fill(photon.pT(), weight);
00117         } else {
00118           _h_central_opp_cross_section->fill(photon.pT(), weight);
00119         }
00120       } else if (inRange( fabs(leadingJet.rapidity()), 1.5, 2.5)) {
00121         if (photon_jet_sign >= 1) {
00122           _h_forward_same_cross_section->fill(photon.pT(), weight);
00123         } else {
00124           _h_forward_opp_cross_section->fill(photon.pT(), weight);
00125         }
00126       }
00127 
00128     }
00129 
00130 
00131 
00132     /// Finalize
00133     void finalize() {
00134       const double lumi_gen = sumOfWeights()/crossSection();
00135       const double dy_photon = 2.0;
00136       const double dy_jet_central = 1.6;
00137       const double dy_jet_forward = 2.0;
00138 
00139       // Cross-section ratios (6 plots)
00140       // Central/central and forward/forward ratios
00141       AIDA::IHistogramFactory& hf = histogramFactory();
00142       const string dir = histoDir();
00143 
00144       hf.divide(dir + "/d05-x01-y01", *_h_central_opp_cross_section, *_h_central_same_cross_section);
00145       hf.divide(dir + "/d08-x01-y01", *_h_forward_opp_cross_section, *_h_forward_same_cross_section);
00146 
00147       // Central/forward ratio combinations
00148       hf.divide(dir + "/d06-x01-y01", *_h_central_same_cross_section,
00149                 *_h_forward_same_cross_section)->scale(dy_jet_forward/dy_jet_central, 1);
00150       hf.divide(dir + "/d07-x01-y01", *_h_central_opp_cross_section,
00151                 *_h_forward_same_cross_section)->scale(dy_jet_forward/dy_jet_central, 1);
00152       hf.divide(dir + "/d09-x01-y01", *_h_central_same_cross_section,
00153                 *_h_forward_opp_cross_section)->scale(dy_jet_forward/dy_jet_central, 1);
00154       hf.divide(dir + "/d10-x01-y01", *_h_central_opp_cross_section,
00155                 *_h_forward_opp_cross_section)->scale(dy_jet_forward/dy_jet_central, 1);
00156 
00157       // Use generator cross section for remaining histograms
00158       // Each of these needs the additional factor 2 because the
00159       // y_photon * y_jet requirement reduces the corresponding 2D "bin width"
00160       // by a factor 1/2.
00161       scale(_h_central_same_cross_section, 2.0/lumi_gen * 1.0/dy_photon * 1.0/dy_jet_central);
00162       scale(_h_central_opp_cross_section, 2.0/lumi_gen * 1.0/dy_photon * 1.0/dy_jet_central);
00163       scale(_h_forward_same_cross_section, 2.0/lumi_gen * 1.0/dy_photon * 1.0/dy_jet_forward);
00164       scale(_h_forward_opp_cross_section, 2.0/lumi_gen * 1.0/dy_photon * 1.0/dy_jet_forward);
00165     }
00166 
00167     //@}
00168 
00169   private:
00170 
00171     /// @name Histograms
00172     //@{
00173     AIDA::IHistogram1D* _h_central_same_cross_section;
00174     AIDA::IHistogram1D* _h_central_opp_cross_section;
00175     AIDA::IHistogram1D* _h_forward_same_cross_section;
00176     AIDA::IHistogram1D* _h_forward_opp_cross_section;
00177     //@}
00178 
00179   };
00180 
00181 
00182 
00183   // The hook for the plugin system
00184   DECLARE_RIVET_PLUGIN(D0_2008_S7719523);
00185 
00186 }