00001
00002 #include "Rivet/Analysis.hh"
00003 #include "Rivet/Tools/Logging.hh"
00004 #include "Rivet/Projections/FinalState.hh"
00005 #include "Rivet/Projections/ZFinder.hh"
00006 #include "Rivet/Projections/FastJets.hh"
00007 #include "Rivet/RivetAIDA.hh"
00008
00009 namespace Rivet {
00010
00011
00012 class D0_2009_S8202443 : public Analysis {
00013
00014 public:
00015
00016
00017
00018
00019 D0_2009_S8202443() : Analysis("D0_2009_S8202443"),
00020 _sum_of_weights(0.0), _sum_of_weights_constrained(0.0)
00021 {
00022 setBeams(PROTON, ANTIPROTON);
00023 }
00024
00025
00026
00027
00028
00029
00030
00031
00032 void init() {
00033
00034 vector<pair<double, double> > etaRanges;
00035 etaRanges.push_back(make_pair(-2.5, -1.5));
00036 etaRanges.push_back(make_pair(-1.1, 1.1));
00037 etaRanges.push_back(make_pair(1.5, 2.5));
00038 ZFinder zfinder_constrained(etaRanges, 25.0*GeV, ELECTRON,
00039 65.0*GeV, 115.0*GeV, 0.2);
00040 addProjection(zfinder_constrained, "ZFinderConstrained");
00041 FastJets conefinder_constrained(zfinder_constrained.remainingFinalState(),
00042 FastJets::D0ILCONE, 0.5);
00043 addProjection(conefinder_constrained, "ConeFinderConstrained");
00044
00045
00046 ZFinder zfinder(FinalState(), ELECTRON, 65.0*GeV, 115.0*GeV, 0.2);
00047 addProjection(zfinder, "ZFinder");
00048 FastJets conefinder(zfinder.remainingFinalState(), FastJets::D0ILCONE, 0.5);
00049 addProjection(conefinder, "ConeFinder");
00050
00051 _h_jet1_pT_constrained = bookHistogram1D(1, 1, 1);
00052 _h_jet2_pT_constrained = bookHistogram1D(3, 1, 1);
00053 _h_jet3_pT_constrained = bookHistogram1D(5, 1, 1);
00054 _h_jet1_pT = bookHistogram1D(2, 1, 1);
00055 _h_jet2_pT = bookHistogram1D(4, 1, 1);
00056 _h_jet3_pT = bookHistogram1D(6, 1, 1);
00057 }
00058
00059
00060
00061
00062 void analyze(const Event& e) {
00063 double weight = e.weight();
00064
00065
00066 const ZFinder& zfinder = applyProjection<ZFinder>(e, "ZFinder");
00067 if (zfinder.particles().size()==1) {
00068 _sum_of_weights += weight;
00069 const JetAlg& jetpro = applyProjection<JetAlg>(e, "ConeFinder");
00070 const Jets& jets = jetpro.jetsByPt(20.0*GeV);
00071 Jets jets_cut;
00072 foreach (const Jet& j, jets) {
00073 if (fabs(j.momentum().pseudorapidity()) < 2.5) {
00074 jets_cut.push_back(j);
00075 }
00076 }
00077
00078 if (jets_cut.size()>0) {
00079 _h_jet1_pT->fill(jets_cut[0].momentum().pT()/GeV, weight);
00080 }
00081 if (jets_cut.size()>1) {
00082 _h_jet2_pT->fill(jets_cut[1].momentum().pT()/GeV, weight);
00083 }
00084 if (jets_cut.size()>2) {
00085 _h_jet3_pT->fill(jets_cut[2].momentum().pT()/GeV, weight);
00086 }
00087 }
00088 else {
00089 getLog() << Log::DEBUG << "no unique lepton pair found." << endl;
00090 }
00091
00092
00093
00094 const ZFinder& zfinder_constrained = applyProjection<ZFinder>(e, "ZFinderConstrained");
00095 if (zfinder_constrained.particles().size()==1) {
00096 _sum_of_weights_constrained += weight;
00097 const JetAlg& jetpro = applyProjection<JetAlg>(e, "ConeFinderConstrained");
00098 const Jets& jets = jetpro.jetsByPt(20.0*GeV);
00099 Jets jets_cut;
00100 foreach (const Jet& j, jets) {
00101 if (fabs(j.momentum().pseudorapidity()) < 2.5) {
00102 jets_cut.push_back(j);
00103 }
00104 }
00105
00106 if (jets_cut.size()>0) {
00107 _h_jet1_pT_constrained->fill(jets_cut[0].momentum().pT()/GeV, weight);
00108 }
00109 if (jets_cut.size()>1) {
00110 _h_jet2_pT_constrained->fill(jets_cut[1].momentum().pT()/GeV, weight);
00111 }
00112 if (jets_cut.size()>2) {
00113 _h_jet3_pT_constrained->fill(jets_cut[2].momentum().pT()/GeV, weight);
00114 }
00115 }
00116 else {
00117 getLog() << Log::DEBUG << "no unique lepton pair found." << endl;
00118 vetoEvent;
00119 }
00120 }
00121
00122
00123
00124
00125 void finalize() {
00126 scale(_h_jet1_pT, 1.0/_sum_of_weights);
00127 scale(_h_jet2_pT, 1.0/_sum_of_weights);
00128 scale(_h_jet3_pT, 1.0/_sum_of_weights);
00129 scale(_h_jet1_pT_constrained, 1.0/_sum_of_weights_constrained);
00130 scale(_h_jet2_pT_constrained, 1.0/_sum_of_weights_constrained);
00131 scale(_h_jet3_pT_constrained, 1.0/_sum_of_weights_constrained);
00132 }
00133
00134
00135
00136
00137 private:
00138
00139
00140
00141 AIDA::IHistogram1D * _h_jet1_pT;
00142 AIDA::IHistogram1D * _h_jet2_pT;
00143 AIDA::IHistogram1D * _h_jet3_pT;
00144 AIDA::IHistogram1D * _h_jet1_pT_constrained;
00145 AIDA::IHistogram1D * _h_jet2_pT_constrained;
00146 AIDA::IHistogram1D * _h_jet3_pT_constrained;
00147
00148
00149 double _sum_of_weights, _sum_of_weights_constrained;
00150
00151 };
00152
00153
00154
00155
00156 AnalysisBuilder<D0_2009_S8202443> plugin_D0_2009_S8202443;
00157
00158 }