00001
00002 #include "Rivet/Analysis.hh"
00003 #include "Rivet/RivetAIDA.hh"
00004 #include "Rivet/Projections/ChargedFinalState.hh"
00005 #include "Rivet/Projections/TriggerCDFRun0Run1.hh"
00006
00007 namespace Rivet {
00008
00009
00010
00011 class CDF_1988_S1865951 : public Analysis {
00012 public:
00013
00014
00015 CDF_1988_S1865951()
00016 : Analysis("CDF_1988_S1865951")
00017 {
00018 setBeams(PROTON, ANTIPROTON);
00019 _sumWTrig = 0;
00020 }
00021
00022
00023
00024
00025
00026
00027 void init() {
00028
00029 addProjection(TriggerCDFRun0Run1(), "Trigger");
00030 const ChargedFinalState cfs(-1.0, 1.0, 0.4*GeV);
00031 addProjection(cfs, "CFS");
00032
00033
00034 if (fuzzyEquals(sqrtS()/GeV, 1800, 1E-3)) {
00035 _hist_pt = bookHistogram1D(1, 1, 1);
00036 } else if (fuzzyEquals(sqrtS()/GeV, 630, 1E-3)) {
00037 _hist_pt = bookHistogram1D(2, 1, 1);
00038 }
00039 }
00040
00041
00042
00043 void analyze(const Event& event) {
00044
00045 const bool trigger = applyProjection<TriggerCDFRun0Run1>(event, "Trigger").minBiasDecision();
00046 if (!trigger) vetoEvent;
00047 const double weight = event.weight();
00048 _sumWTrig += weight;
00049
00050 const FinalState& trackfs = applyProjection<ChargedFinalState>(event, "CFS");
00051 foreach (Particle p, trackfs.particles()) {
00052 const double pt = p.momentum().pT()/GeV;
00053
00054 const double eff_weight = weight/(2*2*TWOPI*pt);
00055 _hist_pt->fill(pt, eff_weight);
00056 }
00057 }
00058
00059
00060
00061 void finalize() {
00062 scale(_hist_pt, crossSectionPerEvent()/millibarn);
00063 }
00064
00065
00066
00067
00068 private:
00069
00070
00071
00072 double _sumWTrig;
00073
00074
00075
00076
00077 AIDA::IHistogram1D* _hist_pt;
00078
00079
00080 };
00081
00082
00083
00084
00085 AnalysisBuilder<CDF_1988_S1865951> plugin_CDF_1988_S1865951;
00086
00087 }