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
00013 class D0_2009_S8202443 : public Analysis {
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, true, true);
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(-MAXRAPIDITY, MAXRAPIDITY, 0.0*GeV, ELECTRON,
00047 65.0*GeV, 115.0*GeV, 0.2, true, true);
00048 addProjection(zfinder, "ZFinder");
00049 FastJets conefinder(zfinder.remainingFinalState(), FastJets::D0ILCONE, 0.5);
00050 addProjection(conefinder, "ConeFinder");
00051
00052 _h_jet1_pT_constrained = bookHistogram1D(1, 1, 1);
00053 _h_jet2_pT_constrained = bookHistogram1D(3, 1, 1);
00054 _h_jet3_pT_constrained = bookHistogram1D(5, 1, 1);
00055 _h_jet1_pT = bookHistogram1D(2, 1, 1);
00056 _h_jet2_pT = bookHistogram1D(4, 1, 1);
00057 _h_jet3_pT = bookHistogram1D(6, 1, 1);
00058 }
00059
00060
00061
00062
00063 void analyze(const Event& e) {
00064 double weight = e.weight();
00065
00066
00067 const ZFinder& zfinder = applyProjection<ZFinder>(e, "ZFinder");
00068 if (zfinder.bosons().size()==1) {
00069 _sum_of_weights += weight;
00070 const JetAlg& jetpro = applyProjection<JetAlg>(e, "ConeFinder");
00071 const Jets& jets = jetpro.jetsByPt(20.0*GeV);
00072 Jets jets_cut;
00073 foreach (const Jet& j, jets) {
00074 if (fabs(j.momentum().pseudorapidity()) < 2.5) {
00075 jets_cut.push_back(j);
00076 }
00077 }
00078
00079 if (jets_cut.size()>0) {
00080 _h_jet1_pT->fill(jets_cut[0].momentum().pT()/GeV, weight);
00081 }
00082 if (jets_cut.size()>1) {
00083 _h_jet2_pT->fill(jets_cut[1].momentum().pT()/GeV, weight);
00084 }
00085 if (jets_cut.size()>2) {
00086 _h_jet3_pT->fill(jets_cut[2].momentum().pT()/GeV, weight);
00087 }
00088 }
00089 else {
00090 getLog() << Log::DEBUG << "no unique lepton pair found." << endl;
00091 }
00092
00093
00094
00095 const ZFinder& zfinder_constrained = applyProjection<ZFinder>(e, "ZFinderConstrained");
00096 if (zfinder_constrained.bosons().size()==1) {
00097 _sum_of_weights_constrained += weight;
00098 const JetAlg& jetpro = applyProjection<JetAlg>(e, "ConeFinderConstrained");
00099 const Jets& jets = jetpro.jetsByPt(20.0*GeV);
00100 Jets jets_cut;
00101 foreach (const Jet& j, jets) {
00102 if (fabs(j.momentum().pseudorapidity()) < 2.5) {
00103 jets_cut.push_back(j);
00104 }
00105 }
00106
00107 if (jets_cut.size()>0) {
00108 _h_jet1_pT_constrained->fill(jets_cut[0].momentum().pT()/GeV, weight);
00109 }
00110 if (jets_cut.size()>1) {
00111 _h_jet2_pT_constrained->fill(jets_cut[1].momentum().pT()/GeV, weight);
00112 }
00113 if (jets_cut.size()>2) {
00114 _h_jet3_pT_constrained->fill(jets_cut[2].momentum().pT()/GeV, weight);
00115 }
00116 }
00117 else {
00118 getLog() << Log::DEBUG << "no unique lepton pair found." << endl;
00119 vetoEvent;
00120 }
00121 }
00122
00123
00124
00125
00126 void finalize() {
00127 scale(_h_jet1_pT, 1.0/_sum_of_weights);
00128 scale(_h_jet2_pT, 1.0/_sum_of_weights);
00129 scale(_h_jet3_pT, 1.0/_sum_of_weights);
00130 scale(_h_jet1_pT_constrained, 1.0/_sum_of_weights_constrained);
00131 scale(_h_jet2_pT_constrained, 1.0/_sum_of_weights_constrained);
00132 scale(_h_jet3_pT_constrained, 1.0/_sum_of_weights_constrained);
00133 }
00134
00135
00136
00137
00138 private:
00139
00140
00141
00142 AIDA::IHistogram1D * _h_jet1_pT;
00143 AIDA::IHistogram1D * _h_jet2_pT;
00144 AIDA::IHistogram1D * _h_jet3_pT;
00145 AIDA::IHistogram1D * _h_jet1_pT_constrained;
00146 AIDA::IHistogram1D * _h_jet2_pT_constrained;
00147 AIDA::IHistogram1D * _h_jet3_pT_constrained;
00148
00149
00150 double _sum_of_weights, _sum_of_weights_constrained;
00151
00152 };
00153
00154
00155
00156
00157 AnalysisBuilder<D0_2009_S8202443> plugin_D0_2009_S8202443;
00158
00159 }