D0_2008_S7719523.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/VetoedFinalState.hh" 00006 #include "Rivet/Projections/FastJets.hh" 00007 00008 namespace Rivet { 00009 00010 00011 // A local scope function for division, handling the div-by-zero case 00012 /// @todo Why isn't the math divide() function being found? 00013 namespace { 00014 inline double _safediv(double a, double b, double result_if_err) { 00015 return (b != 0) ? a/b : result_if_err; 00016 } 00017 } 00018 00019 00020 /// @brief Measurement of isolated gamma + jet + X differential cross-sections 00021 /// 00022 /// Inclusive isolated gamma + jet cross-sections, differential in pT(gamma), for 00023 /// various photon and jet rapidity bins. 00024 /// 00025 /// @author Andy Buckley 00026 /// @author Gavin Hesketh 00027 class D0_2008_S7719523 : public Analysis { 00028 public: 00029 00030 /// @name Constructors etc. 00031 //@{ 00032 00033 /// Constructor 00034 D0_2008_S7719523() 00035 : Analysis("D0_2008_S7719523") 00036 { } 00037 00038 //@} 00039 00040 00041 /// @name Analysis methods 00042 //@{ 00043 00044 /// Set up projections and book histograms 00045 void init() { 00046 // General FS 00047 FinalState fs; 00048 addProjection(fs, "FS"); 00049 00050 // Get leading photon 00051 LeadingParticlesFinalState photonfs(FinalState(-1.0, 1.0, 30.0*GeV)); 00052 photonfs.addParticleId(PID::PHOTON); 00053 addProjection(photonfs, "LeadingPhoton"); 00054 00055 // FS excluding the leading photon 00056 VetoedFinalState vfs(fs); 00057 vfs.addVetoOnThisFinalState(photonfs); 00058 addProjection(vfs, "JetFS"); 00059 00060 // Jets 00061 FastJets jetpro(vfs, FastJets::D0ILCONE, 0.7); 00062 addProjection(jetpro, "Jets"); 00063 00064 // Histograms 00065 _h_central_same_cross_section = bookHisto1D(1, 1, 1); 00066 _h_central_opp_cross_section = bookHisto1D(2, 1, 1); 00067 _h_forward_same_cross_section = bookHisto1D(3, 1, 1); 00068 _h_forward_opp_cross_section = bookHisto1D(4, 1, 1); 00069 00070 // Ratio histos to be filled by divide() 00071 _h_cen_opp_same = bookScatter2D(5, 1, 1); 00072 _h_fwd_opp_same = bookScatter2D(8, 1, 1); 00073 // Ratio histos to be filled manually, since the num/denom inputs don't match 00074 _h_cen_same_fwd_same = bookScatter2D(6, 1, 1, true); 00075 _h_cen_opp_fwd_same = bookScatter2D(7, 1, 1, true); 00076 _h_cen_same_fwd_opp = bookScatter2D(9, 1, 1, true); 00077 _h_cen_opp_fwd_opp = bookScatter2D(10, 1, 1, true); 00078 } 00079 00080 00081 00082 /// Do the analysis 00083 void analyze(const Event& event) { 00084 const double weight = event.weight(); 00085 00086 // Get the photon 00087 const FinalState& photonfs = applyProjection<FinalState>(event, "LeadingPhoton"); 00088 if (photonfs.particles().size() != 1) { 00089 vetoEvent; 00090 } 00091 const FourMomentum photon = photonfs.particles().front().momentum(); 00092 00093 // Isolate photon by ensuring that a 0.4 cone around it contains less than 7% of the photon's energy 00094 double egamma = photon.E(); 00095 double eta_P = photon.eta(); 00096 double phi_P = photon.phi(); 00097 double econe = 0.0; 00098 foreach (const Particle& p, applyProjection<FinalState>(event, "JetFS").particles()) { 00099 if (deltaR(eta_P, phi_P, p.eta(), p.phi()) < 0.4) { 00100 econe += p.E(); 00101 // Veto as soon as E_cone gets larger 00102 if (econe/egamma > 0.07) { 00103 MSG_DEBUG("Vetoing event because photon is insufficiently isolated"); 00104 vetoEvent; 00105 } 00106 } 00107 } 00108 00109 Jets jets = applyProjection<FastJets>(event, "Jets").jetsByPt(15.0*GeV); 00110 if (jets.empty()) vetoEvent; 00111 00112 FourMomentum leadingJet = jets[0].momentum(); 00113 if (deltaR(eta_P, phi_P, leadingJet.eta(), leadingJet.phi()) < 0.7) { 00114 vetoEvent; 00115 } 00116 00117 int photon_jet_sign = sign( leadingJet.rapidity() * photon.rapidity() ); 00118 00119 // Veto if leading jet is outside plotted rapidity regions 00120 const double abs_y1 = fabs(leadingJet.rapidity()); 00121 if (inRange(abs_y1, 0.8, 1.5) || abs_y1 > 2.5) { 00122 MSG_DEBUG("Leading jet falls outside acceptance range; |y1| = " << abs_y1); 00123 vetoEvent; 00124 } 00125 00126 // Fill histos 00127 if (fabs(leadingJet.rapidity()) < 0.8) { 00128 Histo1DPtr h = (photon_jet_sign >= 1) ? _h_central_same_cross_section : _h_central_opp_cross_section; 00129 h->fill(photon.pT(), weight); 00130 } else if (inRange( fabs(leadingJet.rapidity()), 1.5, 2.5)) { 00131 Histo1DPtr h = (photon_jet_sign >= 1) ? _h_forward_same_cross_section : _h_forward_opp_cross_section; 00132 h->fill(photon.pT(), weight); 00133 } 00134 00135 } 00136 00137 00138 /// Finalize 00139 void finalize() { 00140 const double lumi_gen = sumOfWeights()/crossSection(); 00141 const double dy_photon = 2.0; 00142 const double dy_jet_central = 1.6; 00143 const double dy_jet_forward = 2.0; 00144 00145 // Cross-section ratios (6 plots) 00146 // Central/central and forward/forward ratios 00147 divide(_h_central_opp_cross_section, _h_central_same_cross_section, _h_cen_opp_same); 00148 divide(_h_forward_opp_cross_section, _h_forward_same_cross_section, _h_fwd_opp_same); 00149 // Central/forward ratio combinations 00150 /// @note The central/forward histo binnings are not the same! Hence the need to do these by hand :-( 00151 for (size_t i = 0; i < _h_cen_same_fwd_same->numPoints(); ++i) { 00152 const YODA::HistoBin1D& cen_same_bini = _h_central_same_cross_section->bin(i); 00153 const YODA::HistoBin1D& cen_opp_bini = _h_central_opp_cross_section->bin(i); 00154 const YODA::HistoBin1D& fwd_same_bini = _h_central_same_cross_section->bin(i); 00155 const YODA::HistoBin1D& fwd_opp_bini = _h_central_opp_cross_section->bin(i); 00156 _h_cen_same_fwd_same->point(i).setY(_safediv(cen_same_bini.sumW(), fwd_same_bini.sumW(), 0), 00157 add_quad(cen_same_bini.relErr(), fwd_same_bini.relErr())); 00158 _h_cen_opp_fwd_same->point(i).setY(_safediv(cen_opp_bini.sumW(), fwd_same_bini.sumW(), 0), 00159 add_quad(cen_opp_bini.relErr(), fwd_same_bini.relErr())); 00160 _h_cen_same_fwd_opp->point(i).setY(_safediv(cen_same_bini.sumW(), fwd_opp_bini.sumW(), 0), 00161 add_quad(cen_same_bini.relErr(), fwd_opp_bini.relErr())); 00162 _h_cen_opp_fwd_opp->point(i).setY(_safediv(cen_opp_bini.sumW(), fwd_opp_bini.sumW(), 0), 00163 add_quad(cen_opp_bini.relErr(), fwd_opp_bini.relErr())); 00164 } 00165 00166 // Use generator cross section for remaining histograms 00167 // Each of these needs the additional factor 2 because the 00168 // y_photon * y_jet requirement reduces the corresponding 2D "bin width" 00169 // by a factor 1/2. 00170 scale(_h_central_same_cross_section, 2.0/lumi_gen * 1.0/dy_photon * 1.0/dy_jet_central); 00171 scale(_h_central_opp_cross_section, 2.0/lumi_gen * 1.0/dy_photon * 1.0/dy_jet_central); 00172 scale(_h_forward_same_cross_section, 2.0/lumi_gen * 1.0/dy_photon * 1.0/dy_jet_forward); 00173 scale(_h_forward_opp_cross_section, 2.0/lumi_gen * 1.0/dy_photon * 1.0/dy_jet_forward); 00174 } 00175 00176 //@} 00177 00178 private: 00179 00180 /// @name Histograms 00181 //@{ 00182 Histo1DPtr _h_central_same_cross_section; 00183 Histo1DPtr _h_central_opp_cross_section; 00184 Histo1DPtr _h_forward_same_cross_section; 00185 Histo1DPtr _h_forward_opp_cross_section; 00186 00187 Scatter2DPtr _h_cen_opp_same; 00188 Scatter2DPtr _h_fwd_opp_same; 00189 Scatter2DPtr _h_cen_same_fwd_same; 00190 Scatter2DPtr _h_cen_opp_fwd_same; 00191 Scatter2DPtr _h_cen_same_fwd_opp; 00192 Scatter2DPtr _h_cen_opp_fwd_opp; 00193 //@} 00194 00195 }; 00196 00197 00198 00199 // The hook for the plugin system 00200 DECLARE_RIVET_PLUGIN(D0_2008_S7719523); 00201 00202 } Generated on Wed Oct 7 2015 12:09:12 for The Rivet MC analysis system by ![]() |