00001
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
00014
00015
00016
00017
00018
00019
00020 class D0_2008_S7719523 : public Analysis {
00021
00022 public:
00023
00024
00025
00026
00027
00028 D0_2008_S7719523()
00029 : Analysis("D0_2008_S7719523")
00030 {
00031 setBeams(PROTON, ANTIPROTON);
00032 setNeedsCrossSection(true);
00033 }
00034
00035
00036
00037
00038
00039
00040
00041
00042 void init() {
00043
00044 FinalState fs;
00045 addProjection(fs, "FS");
00046
00047
00048 LeadingParticlesFinalState photonfs(FinalState(-1.0, 1.0, 30.0*GeV));
00049 photonfs.addParticleId(PHOTON);
00050 addProjection(photonfs, "LeadingPhoton");
00051
00052
00053 VetoedFinalState vfs(fs);
00054 vfs.addVetoOnThisFinalState(photonfs);
00055 addProjection(vfs, "JetFS");
00056
00057
00058 FastJets jetpro(vfs, FastJets::D0ILCONE, 0.7);
00059 addProjection(jetpro, "Jets");
00060
00061
00062 _h_central_same_cross_section = bookHistogram1D(1, 1, 1);
00063 _h_central_opp_cross_section = bookHistogram1D(2, 1, 1);
00064 _h_forward_same_cross_section = bookHistogram1D(3, 1, 1);
00065 _h_forward_opp_cross_section = bookHistogram1D(4, 1, 1);
00066 }
00067
00068
00069
00070
00071 void analyze(const Event& event) {
00072 const double weight = event.weight();
00073
00074
00075 const FinalState& photonfs = applyProjection<FinalState>(event, "LeadingPhoton");
00076 if (photonfs.particles().size() != 1) {
00077 vetoEvent;
00078 }
00079 const FourMomentum photon = photonfs.particles().front().momentum();
00080
00081
00082 double egamma = photon.E();
00083 double eta_P = photon.pseudorapidity();
00084 double phi_P = photon.azimuthalAngle();
00085 double econe = 0.0;
00086 foreach (const Particle& p, applyProjection<FinalState>(event, "JetFS").particles()) {
00087 if (deltaR(eta_P, phi_P,
00088 p.momentum().pseudorapidity(), p.momentum().azimuthalAngle()) < 0.4) {
00089 econe += p.momentum().E();
00090
00091 if (econe/egamma > 0.07) {
00092 getLog() << Log::DEBUG << "Vetoing event because photon is insufficiently isolated" << endl;
00093 vetoEvent;
00094 }
00095 }
00096 }
00097
00098 Jets jets = applyProjection<FastJets>(event, "Jets").jetsByPt(15.0*GeV);
00099 if (jets.size()==0) {
00100 vetoEvent;
00101 }
00102 FourMomentum leadingJet = jets[0].momentum();
00103 if (deltaR(eta_P, phi_P, leadingJet.eta(), leadingJet.phi())<0.7) {
00104 vetoEvent;
00105 }
00106
00107 int photon_jet_sign = sign( leadingJet.rapidity() * photon.rapidity() );
00108
00109
00110 const double abs_y1 = fabs(leadingJet.rapidity());
00111 if (inRange(abs_y1, 0.8, 1.5) || abs_y1 > 2.5) {
00112 getLog() << Log::DEBUG << "Leading jet falls outside acceptance range; |y1| = "
00113 << abs_y1 << endl;
00114 vetoEvent;
00115 }
00116
00117
00118 if (fabs(leadingJet.rapidity()) < 0.8) {
00119 if (photon_jet_sign >= 1) {
00120 _h_central_same_cross_section->fill(photon.pT(), weight);
00121 } else {
00122 _h_central_opp_cross_section->fill(photon.pT(), weight);
00123 }
00124 } else if (inRange( fabs(leadingJet.rapidity()), 1.5, 2.5)) {
00125 if (photon_jet_sign >= 1) {
00126 _h_forward_same_cross_section->fill(photon.pT(), weight);
00127 } else {
00128 _h_forward_opp_cross_section->fill(photon.pT(), weight);
00129 }
00130 }
00131
00132 }
00133
00134
00135
00136
00137 void finalize() {
00138 const double lumi_gen = sumOfWeights()/crossSection();
00139 const double dy_photon = 2.0;
00140 const double dy_jet_central = 1.6;
00141 const double dy_jet_forward = 2.0;
00142
00143
00144
00145 AIDA::IHistogramFactory& hf = histogramFactory();
00146 const string dir = histoDir();
00147
00148 hf.divide(dir + "/d05-x01-y01", *_h_central_opp_cross_section, *_h_central_same_cross_section);
00149 hf.divide(dir + "/d08-x01-y01", *_h_forward_opp_cross_section, *_h_forward_same_cross_section);
00150
00151
00152 hf.divide(dir + "/d06-x01-y01", *_h_central_same_cross_section,
00153 *_h_forward_same_cross_section)->scale(dy_jet_forward/dy_jet_central, 1);
00154 hf.divide(dir + "/d07-x01-y01", *_h_central_opp_cross_section,
00155 *_h_forward_same_cross_section)->scale(dy_jet_forward/dy_jet_central, 1);
00156 hf.divide(dir + "/d09-x01-y01", *_h_central_same_cross_section,
00157 *_h_forward_opp_cross_section)->scale(dy_jet_forward/dy_jet_central, 1);
00158 hf.divide(dir + "/d10-x01-y01", *_h_central_opp_cross_section,
00159 *_h_forward_opp_cross_section)->scale(dy_jet_forward/dy_jet_central, 1);
00160
00161
00162
00163
00164
00165 scale(_h_central_same_cross_section, 2.0/lumi_gen * 1.0/dy_photon * 1.0/dy_jet_central);
00166 scale(_h_central_opp_cross_section, 2.0/lumi_gen * 1.0/dy_photon * 1.0/dy_jet_central);
00167 scale(_h_forward_same_cross_section, 2.0/lumi_gen * 1.0/dy_photon * 1.0/dy_jet_forward);
00168 scale(_h_forward_opp_cross_section, 2.0/lumi_gen * 1.0/dy_photon * 1.0/dy_jet_forward);
00169 }
00170
00171
00172
00173 private:
00174
00175
00176
00177 AIDA::IHistogram1D* _h_central_same_cross_section;
00178 AIDA::IHistogram1D* _h_central_opp_cross_section;
00179 AIDA::IHistogram1D* _h_forward_same_cross_section;
00180 AIDA::IHistogram1D* _h_forward_opp_cross_section;
00181
00182
00183 };
00184
00185
00186
00187
00188 AnalysisBuilder<D0_2008_S7719523> plugin_D0_2008_S7719523;
00189
00190 }