CDF_1997_S3541940.cc
Go to the documentation of this file.
00001 // -*- C++ -*- 00002 #include "Rivet/Analysis.hh" 00003 #include "Rivet/Projections/FinalState.hh" 00004 #include "Rivet/Projections/FastJets.hh" 00005 00006 namespace Rivet { 00007 00008 00009 /// @brief CDF properties of 6-jet events with large 6-jet mass 00010 class CDF_1997_S3541940 : public Analysis { 00011 public: 00012 00013 CDF_1997_S3541940() 00014 : Analysis("CDF_1997_S3541940") 00015 { 00016 } 00017 00018 00019 public: 00020 00021 void init() { 00022 00023 const FinalState fs(-4.2, 4.2); 00024 addProjection(FastJets(fs, FastJets::CDFJETCLU, 0.7), "Jets"); 00025 00026 _h_m6J = bookHisto1D(1, 1, 1); 00027 _h_X3ppp = bookHisto1D(2, 1, 1); 00028 _h_X4ppp = bookHisto1D(3, 1, 1); 00029 _h_costheta3ppp = bookHisto1D(4, 1, 1); 00030 _h_psi3ppp = bookHisto1D(5, 1, 1); 00031 _h_f3ppp = bookHisto1D(6, 1, 1); 00032 _h_f4ppp = bookHisto1D(6, 1, 2); 00033 _h_f5ppp = bookHisto1D(6, 1, 3); 00034 _h_XApp = bookHisto1D(7, 1, 1); 00035 _h_XCp = bookHisto1D(8, 1, 1); 00036 _h_XE = bookHisto1D(9, 1, 1); 00037 _h_psiAppBpp = bookHisto1D(10, 1, 1); 00038 _h_psiCpDp = bookHisto1D(11, 1, 1); 00039 _h_psiEF = bookHisto1D(12, 1, 1); 00040 _h_fApp = bookHisto1D(13, 1, 1); 00041 _h_fBpp = bookHisto1D(14, 1, 1); 00042 _h_fCp = bookHisto1D(15, 1, 1); 00043 _h_fDp = bookHisto1D(16, 1, 1); 00044 _h_fE = bookHisto1D(17, 1, 1); 00045 _h_fF = bookHisto1D(18, 1, 1); 00046 } 00047 00048 00049 void analyze(const Event& event) { 00050 const double weight = event.weight(); 00051 00052 Jets jets; 00053 double sumEt = 0.0; 00054 FourMomentum jetsystem(0.0, 0.0, 0.0, 0.0); 00055 foreach (const Jet& jet, applyProjection<FastJets>(event, "Jets").jetsByEt()) { 00056 double Et = jet.Et(); 00057 double eta = jet.abseta(); 00058 if (Et > 20.0*GeV && eta < 3.0) { 00059 bool separated=true; 00060 foreach (const Jet& ref, jets) { 00061 if (deltaR(jet.momentum(), ref.momentum()) < 0.9) { 00062 separated=false; 00063 break; 00064 } 00065 } 00066 if (!separated) continue; 00067 jets.push_back(jet); 00068 sumEt += Et; 00069 jetsystem += jet.momentum(); 00070 } 00071 if (jets.size() >= 6) break; 00072 } 00073 00074 if (jets.size() < 6) { 00075 vetoEvent; 00076 } 00077 00078 if (sumEt < 320.0*GeV) { 00079 vetoEvent; 00080 } 00081 00082 double m6J = _safeMass(jetsystem); 00083 if (m6J < 520.0*GeV) { 00084 vetoEvent; 00085 } 00086 00087 LorentzTransform cms_boost(-jetsystem.boostVector()); 00088 vector<FourMomentum> jets6; 00089 foreach (Jet jet, jets) { 00090 jets6.push_back(cms_boost.transform(jet.momentum())); 00091 } 00092 std::sort(jets6.begin(), jets6.end(), FourMomentum::byEDescending()); 00093 00094 FourMomentum pE, pF; 00095 vector<FourMomentum> jets5(_reduce(jets6, pE, pF)); 00096 std::sort(jets5.begin(), jets5.end(), FourMomentum::byEDescending()); 00097 00098 FourMomentum pCp, pDp; 00099 vector<FourMomentum> jets4(_reduce(jets5, pCp, pDp)); 00100 std::sort(jets4.begin(), jets4.end(), FourMomentum::byEDescending()); 00101 00102 FourMomentum pApp, pBpp; 00103 vector<FourMomentum> jets3(_reduce(jets4, pApp, pBpp)); 00104 std::sort(jets3.begin(), jets3.end(), FourMomentum::byEDescending()); 00105 FourMomentum p3ppp(jets3[0]); 00106 FourMomentum p4ppp(jets3[1]); 00107 FourMomentum p5ppp(jets3[2]); 00108 00109 double X3ppp = 2.0*p3ppp.E()/m6J; 00110 if (X3ppp > 0.9) { 00111 vetoEvent; 00112 } 00113 00114 FourMomentum pAV = cms_boost.transform(_avg_beam_in_lab(m6J, jetsystem.rapidity())); 00115 double costheta3ppp = pAV.p3().unit().dot(p3ppp.p3().unit()); 00116 if (fabs(costheta3ppp) > 0.9) { 00117 vetoEvent; 00118 } 00119 00120 // 3-jet-system variables 00121 _h_m6J->fill(m6J, weight); 00122 _h_X3ppp->fill(X3ppp, weight); 00123 _h_X4ppp->fill(2.0*p4ppp.E()/m6J, weight); 00124 _h_costheta3ppp->fill(costheta3ppp, weight); 00125 double psi3ppp = _psi(p3ppp, pAV, p4ppp, p5ppp); 00126 _h_psi3ppp->fill(psi3ppp, weight); 00127 _h_f3ppp->fill(_safeMass(p3ppp)/m6J, weight); 00128 _h_f4ppp->fill(_safeMass(p4ppp)/m6J, weight); 00129 _h_f5ppp->fill(_safeMass(p5ppp)/m6J, weight); 00130 00131 // 4 -> 3 jet variables 00132 _h_fApp->fill(_safeMass(pApp)/m6J, weight); 00133 _h_fBpp->fill(_safeMass(pApp)/m6J, weight); 00134 _h_XApp->fill(pApp.E()/(pApp.E()+pBpp.E()), weight); 00135 double psiAppBpp = _psi(pApp, pBpp, pApp+pBpp, pAV); 00136 _h_psiAppBpp->fill(psiAppBpp, weight); 00137 00138 // 5 -> 4 jet variables 00139 _h_fCp->fill(_safeMass(pCp)/m6J, weight); 00140 _h_fDp->fill(_safeMass(pDp)/m6J, weight); 00141 _h_XCp->fill(pCp.E()/(pCp.E()+pDp.E()), weight); 00142 double psiCpDp = _psi(pCp, pDp, pCp+pDp, pAV); 00143 _h_psiCpDp->fill(psiCpDp, weight); 00144 00145 // 6 -> 5 jet variables 00146 _h_fE->fill(_safeMass(pE)/m6J, weight); 00147 _h_fF->fill(_safeMass(pF)/m6J, weight); 00148 _h_XE->fill(pE.E()/(pE.E()+pF.E()), weight); 00149 double psiEF = _psi(pE, pF, pE+pF, pAV); 00150 _h_psiEF->fill(psiEF, weight); 00151 } 00152 00153 00154 void finalize() { 00155 normalize(_h_m6J); 00156 normalize(_h_X3ppp); 00157 normalize(_h_X4ppp); 00158 normalize(_h_costheta3ppp); 00159 normalize(_h_psi3ppp); 00160 normalize(_h_f3ppp); 00161 normalize(_h_f4ppp); 00162 normalize(_h_f5ppp); 00163 normalize(_h_XApp); 00164 normalize(_h_XCp); 00165 normalize(_h_XE); 00166 normalize(_h_psiAppBpp); 00167 normalize(_h_psiCpDp); 00168 normalize(_h_psiEF); 00169 normalize(_h_fApp); 00170 normalize(_h_fBpp); 00171 normalize(_h_fCp); 00172 normalize(_h_fDp); 00173 normalize(_h_fE); 00174 normalize(_h_fF); 00175 } 00176 00177 00178 private: 00179 00180 vector<FourMomentum> _reduce(const vector<FourMomentum>& jets, 00181 FourMomentum& combined1, 00182 FourMomentum& combined2) { 00183 double minMass2 = 1e9; 00184 size_t idx1(jets.size()), idx2(jets.size()); 00185 for (size_t i = 0; i < jets.size(); ++i) { 00186 for (size_t j = i+1; j < jets.size(); ++j) { 00187 double mass2 = FourMomentum(jets[i] + jets[j]).mass2(); 00188 if (mass2 < minMass2) { 00189 idx1 = i; 00190 idx2 = j; 00191 } 00192 } 00193 } 00194 vector<FourMomentum> newjets; 00195 for (size_t i = 0; i < jets.size(); ++i) { 00196 if (i != idx1 && i != idx2) newjets.push_back(jets[i]); 00197 } 00198 newjets.push_back(jets[idx1] + jets[idx2]); 00199 combined1 = jets[idx1]; 00200 combined2 = jets[idx2]; 00201 return newjets; 00202 } 00203 00204 00205 FourMomentum _avg_beam_in_lab(const double& m, const double& y) { 00206 const double mt = m/2.0; 00207 FourMomentum beam1(mt, 0, 0, mt); 00208 FourMomentum beam2(mt, 0, 0, -mt); 00209 if (fabs(y) > 1e-3) { 00210 FourMomentum boostvec(cosh(y), 0.0, 0.0, sinh(y)); 00211 LorentzTransform cms_boost(-boostvec.boostVector()); 00212 cms_boost = cms_boost.inverse(); 00213 beam1 = cms_boost.transform(beam1); 00214 beam2 = cms_boost.transform(beam2); 00215 } 00216 if (beam1.E() > beam2.E()) { 00217 return beam1 - beam2; 00218 } 00219 else { 00220 return beam2 - beam1; 00221 } 00222 } 00223 00224 00225 double _psi(const FourMomentum& p1, const FourMomentum& p2, 00226 const FourMomentum& p3, const FourMomentum& p4) { 00227 Vector3 p1xp2 = p1.p3().cross(p2.p3()); 00228 Vector3 p3xp4 = p3.p3().cross(p4.p3()); 00229 return mapAngle0ToPi(acos(p1xp2.unit().dot(p3xp4.unit()))); 00230 } 00231 00232 00233 double _safeMass(const FourMomentum& p) { 00234 double mass2 = p.mass2(); 00235 if (mass2 > 0.0) return sqrt(mass2); 00236 if (mass2 < -1e-5) MSG_WARNING("m2 = " << m2 << ". Assuming m2=0."); 00237 return 0.0; 00238 } 00239 00240 00241 private: 00242 00243 Histo1DPtr _h_m6J; 00244 Histo1DPtr _h_X3ppp, _h_X4ppp; 00245 Histo1DPtr _h_costheta3ppp; 00246 Histo1DPtr _h_psi3ppp; 00247 Histo1DPtr _h_f3ppp, _h_f4ppp, _h_f5ppp; 00248 Histo1DPtr _h_XApp, _h_XCp, _h_XE; 00249 Histo1DPtr _h_psiAppBpp, _h_psiCpDp, _h_psiEF; 00250 Histo1DPtr _h_fApp, _h_fBpp, _h_fCp, _h_fDp, _h_fE, _h_fF; 00251 00252 }; 00253 00254 00255 00256 // The hook for the plugin system 00257 DECLARE_RIVET_PLUGIN(CDF_1997_S3541940); 00258 00259 } Generated on Tue Sep 30 2014 19:45:43 for The Rivet MC analysis system by ![]() |