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/RivetYODA.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 = bookHisto1D(1, 1, 1); 00059 _h_central_opp_cross_section = bookHisto1D(2, 1, 1); 00060 _h_forward_same_cross_section = bookHisto1D(3, 1, 1); 00061 _h_forward_opp_cross_section = bookHisto1D(4, 1, 1); 00062 00063 _h_cen_opp_same = bookScatter2D(5, 1, 1); 00064 _h_fwd_opp_same = bookScatter2D(8, 1, 1); 00065 _h_cen_same_fwd_same = bookScatter2D(6, 1, 1); 00066 _h_cen_opp_fwd_same = bookScatter2D(7, 1, 1); 00067 _h_cen_same_fwd_opp = bookScatter2D(9, 1, 1); 00068 _h_cen_opp_fwd_opp = bookScatter2D(10, 1, 1); 00069 00070 } 00071 00072 00073 00074 /// Do the analysis 00075 void analyze(const Event& event) { 00076 const double weight = event.weight(); 00077 00078 // Get the photon 00079 const FinalState& photonfs = applyProjection<FinalState>(event, "LeadingPhoton"); 00080 if (photonfs.particles().size() != 1) { 00081 vetoEvent; 00082 } 00083 const FourMomentum photon = photonfs.particles().front().momentum(); 00084 00085 // Isolate photon by ensuring that a 0.4 cone around it contains less than 7% of the photon's energy 00086 double egamma = photon.E(); 00087 double eta_P = photon.pseudorapidity(); 00088 double phi_P = photon.azimuthalAngle(); 00089 double econe = 0.0; 00090 foreach (const Particle& p, applyProjection<FinalState>(event, "JetFS").particles()) { 00091 if (deltaR(eta_P, phi_P, 00092 p.momentum().pseudorapidity(), p.momentum().azimuthalAngle()) < 0.4) { 00093 econe += p.momentum().E(); 00094 // Veto as soon as E_cone gets larger 00095 if (econe/egamma > 0.07) { 00096 MSG_DEBUG("Vetoing event because photon is insufficiently isolated"); 00097 vetoEvent; 00098 } 00099 } 00100 } 00101 00102 Jets jets = applyProjection<FastJets>(event, "Jets").jetsByPt(15.0*GeV); 00103 if (jets.size()==0) { 00104 vetoEvent; 00105 } 00106 FourMomentum leadingJet = jets[0].momentum(); 00107 if (deltaR(eta_P, phi_P, leadingJet.eta(), leadingJet.phi())<0.7) { 00108 vetoEvent; 00109 } 00110 00111 int photon_jet_sign = sign( leadingJet.rapidity() * photon.rapidity() ); 00112 00113 // Veto if leading jet is outside plotted rapidity regions 00114 const double abs_y1 = fabs(leadingJet.rapidity()); 00115 if (inRange(abs_y1, 0.8, 1.5) || abs_y1 > 2.5) { 00116 MSG_DEBUG("Leading jet falls outside acceptance range; |y1| = " 00117 << abs_y1); 00118 vetoEvent; 00119 } 00120 00121 // Fill histos 00122 if (fabs(leadingJet.rapidity()) < 0.8) { 00123 if (photon_jet_sign >= 1) { 00124 _h_central_same_cross_section->fill(photon.pT(), weight); 00125 } else { 00126 _h_central_opp_cross_section->fill(photon.pT(), weight); 00127 } 00128 } else if (inRange( fabs(leadingJet.rapidity()), 1.5, 2.5)) { 00129 if (photon_jet_sign >= 1) { 00130 _h_forward_same_cross_section->fill(photon.pT(), weight); 00131 } else { 00132 _h_forward_opp_cross_section->fill(photon.pT(), weight); 00133 } 00134 } 00135 00136 } 00137 00138 00139 00140 /// Finalize 00141 void finalize() { 00142 const double lumi_gen = sumOfWeights()/crossSection(); 00143 const double dy_photon = 2.0; 00144 const double dy_jet_central = 1.6; 00145 const double dy_jet_forward = 2.0; 00146 00147 00148 // Cross-section ratios (6 plots) 00149 // Central/central and forward/forward ratios 00150 00151 divide(_h_central_opp_cross_section,_h_central_same_cross_section, 00152 _h_cen_opp_same); 00153 00154 divide(_h_forward_opp_cross_section,_h_forward_same_cross_section, 00155 _h_fwd_opp_same); 00156 00157 00158 // Central/forward ratio combinations 00159 00160 divide(_h_central_same_cross_section,_h_forward_same_cross_section, 00161 _h_cen_same_fwd_same); 00162 00163 divide(_h_central_opp_cross_section,_h_forward_same_cross_section, 00164 _h_cen_opp_fwd_same); 00165 00166 divide(_h_central_same_cross_section,_h_forward_opp_cross_section, 00167 _h_cen_same_fwd_opp); 00168 00169 divide(_h_central_opp_cross_section,_h_forward_opp_cross_section, 00170 _h_cen_opp_fwd_opp); 00171 00172 // Use generator cross section for remaining histograms 00173 // Each of these needs the additional factor 2 because the 00174 // y_photon * y_jet requirement reduces the corresponding 2D "bin width" 00175 // by a factor 1/2. 00176 scale(_h_central_same_cross_section, 2.0/lumi_gen * 1.0/dy_photon * 1.0/dy_jet_central); 00177 scale(_h_central_opp_cross_section, 2.0/lumi_gen * 1.0/dy_photon * 1.0/dy_jet_central); 00178 scale(_h_forward_same_cross_section, 2.0/lumi_gen * 1.0/dy_photon * 1.0/dy_jet_forward); 00179 scale(_h_forward_opp_cross_section, 2.0/lumi_gen * 1.0/dy_photon * 1.0/dy_jet_forward); 00180 } 00181 00182 //@} 00183 00184 private: 00185 00186 /// @name Histograms 00187 //@{ 00188 Histo1DPtr _h_central_same_cross_section; 00189 Histo1DPtr _h_central_opp_cross_section; 00190 Histo1DPtr _h_forward_same_cross_section; 00191 Histo1DPtr _h_forward_opp_cross_section; 00192 00193 Scatter2DPtr _h_cen_opp_same; 00194 Scatter2DPtr _h_fwd_opp_same; 00195 Scatter2DPtr _h_cen_same_fwd_same; 00196 Scatter2DPtr _h_cen_opp_fwd_same; 00197 Scatter2DPtr _h_cen_same_fwd_opp; 00198 Scatter2DPtr _h_cen_opp_fwd_opp; 00199 //@} 00200 00201 }; 00202 00203 00204 00205 // The hook for the plugin system 00206 DECLARE_RIVET_PLUGIN(D0_2008_S7719523); 00207 00208 } Generated on Fri Dec 21 2012 15:03:40 for The Rivet MC analysis system by ![]() |