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