00001
00002 #include "Rivet/Analysis.hh"
00003 #include "Rivet/RivetAIDA.hh"
00004 #include "Rivet/Tools/Logging.hh"
00005 #include "Rivet/Projections/FastJets.hh"
00006 #include "Rivet/Projections/FinalState.hh"
00007
00008 namespace Rivet {
00009
00010 namespace {
00011
00012
00013
00014
00015
00016
00017
00018 double calc_BZ(const vector<fastjet::PseudoJet>& jets) {
00019 assert(jets.size() == 4);
00020 Vector3 p12 = cross( momentum3(jets[0]), momentum3(jets[1]));
00021 Vector3 p34 = cross( momentum3(jets[2]), momentum3(jets[3]));
00022 return dot(p12,p34) / (p12.mod()*p34.mod());
00023 }
00024
00025
00026
00027
00028 double calc_KSW(const vector<fastjet::PseudoJet>& jets) {
00029 assert(jets.size() == 4);
00030 Vector3 p13 = cross( momentum3(jets[0]), momentum3(jets[2]));
00031 Vector3 p24 = cross( momentum3(jets[1]), momentum3(jets[3]));
00032 Vector3 p14 = cross( momentum3(jets[0]), momentum3(jets[3]));
00033 Vector3 p23 = cross( momentum3(jets[1]), momentum3(jets[2]));
00034 return cos (0.5*( acos (dot(p14,p23) / (p14.mod()*p23.mod())) +
00035 acos (dot(p13,p24) / (p13.mod()*p24.mod())) ));
00036 }
00037
00038
00039
00040
00041 double calc_NR(const vector<fastjet::PseudoJet>& jets) {
00042 assert(jets.size() == 4);
00043 Vector3 p12 = momentum3(jets[0]) - momentum3(jets[1]);
00044 Vector3 p34 = momentum3(jets[2]) - momentum3(jets[3]);
00045 return dot(p12,p34) / (p12.mod()*p34.mod());
00046 }
00047
00048
00049
00050 double calc_ALPHA34(const vector<fastjet::PseudoJet>& jets) {
00051 assert(jets.size() == 4);
00052 Vector3 p3 = momentum3(jets[2]);
00053 Vector3 p4 = momentum3(jets[3]);
00054 return dot(p3,p4) / (p3.mod()*p4.mod());
00055 }
00056
00057
00058
00059 }
00060
00061
00062 class OPAL_2001_S4553896 : public Analysis {
00063 public:
00064
00065
00066
00067
00068
00069 OPAL_2001_S4553896()
00070 : Analysis("OPAL_2001_S4553896")
00071 {
00072 setBeams(ELECTRON, POSITRON);
00073 setNeedsCrossSection(false);
00074 }
00075
00076
00077
00078
00079 public:
00080
00081
00082
00083
00084
00085 void init() {
00086
00087
00088 const FinalState fs;
00089 addProjection(fs, "FS");
00090 addProjection(FastJets(fs, FastJets::DURHAM, 0.7), "Jets");
00091
00092
00093
00094 _h_BZ = bookHistogram1D(3, 1, 1);
00095 _h_KSW = bookHistogram1D(4, 1, 1);
00096 _h_NR = bookHistogram1D(5, 1, 1);
00097 _h_ALPHA34 = bookHistogram1D(6, 1, 1);
00098 }
00099
00100
00101
00102 void analyze(const Event& event) {
00103 const double weight = event.weight();
00104
00105
00106 if (applyProjection<FinalState>(event, "FS").particles().size() < 2) {
00107 vetoEvent;
00108 }
00109
00110 const FastJets& fastjets = applyProjection<FastJets>(event, "Jets");
00111 if (fastjets.clusterSeq()) {
00112 vector<fastjet::PseudoJet> jets;
00113 foreach (const fastjet::PseudoJet& jet,
00114 fastjet::sorted_by_E(fastjets.clusterSeq()->exclusive_jets_ycut(0.008))) {
00115 if (jet.E()>3.0*GeV) jets.push_back(jet);
00116 }
00117 if (jets.size() == 4) {
00118 _h_BZ->fill(fabs(calc_BZ(jets)), weight);
00119 _h_KSW->fill(calc_KSW(jets), weight);
00120 _h_NR->fill(fabs(calc_NR(jets)), weight);
00121 _h_ALPHA34->fill(calc_ALPHA34(jets), weight);
00122 }
00123 }
00124
00125
00126 }
00127
00128
00129
00130 void finalize() {
00131
00132
00133 normalize(_h_BZ);
00134 normalize(_h_KSW);
00135 normalize(_h_NR);
00136 normalize(_h_ALPHA34);
00137
00138 }
00139
00140
00141
00142
00143 private:
00144
00145
00146
00147
00148 AIDA::IHistogram1D *_h_BZ;
00149 AIDA::IHistogram1D *_h_KSW;
00150 AIDA::IHistogram1D *_h_NR;
00151 AIDA::IHistogram1D *_h_ALPHA34;
00152
00153
00154 };
00155
00156
00157
00158
00159 AnalysisBuilder<OPAL_2001_S4553896> plugin_OPAL_2001_S4553896;
00160
00161
00162 }