rivet is hosted by Hepforge, IPPP Durham
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").jets(cmpMomByEt)) {
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) vetoEvent;
00075       if (sumEt < 320.0*GeV) vetoEvent;
00076 
00077       double m6J = _safeMass(jetsystem);
00078       if (m6J < 520.0*GeV) vetoEvent;
00079 
00080       LorentzTransform cms_boost(-jetsystem.boostVector());
00081       vector<FourMomentum> jets6;
00082       foreach (Jet jet, jets) {
00083         jets6.push_back(cms_boost.transform(jet.momentum()));
00084       }
00085       std::sort(jets6.begin(), jets6.end(), FourMomentum::byEDescending());
00086 
00087       FourMomentum pE, pF;
00088       vector<FourMomentum> jets5(_reduce(jets6, pE, pF));
00089       std::sort(jets5.begin(), jets5.end(), FourMomentum::byEDescending());
00090 
00091       FourMomentum pCp, pDp;
00092       vector<FourMomentum> jets4(_reduce(jets5, pCp, pDp));
00093       std::sort(jets4.begin(), jets4.end(), FourMomentum::byEDescending());
00094 
00095       FourMomentum pApp, pBpp;
00096       vector<FourMomentum> jets3(_reduce(jets4, pApp, pBpp));
00097       std::sort(jets3.begin(), jets3.end(), FourMomentum::byEDescending());
00098       FourMomentum p3ppp(jets3[0]);
00099       FourMomentum p4ppp(jets3[1]);
00100       FourMomentum p5ppp(jets3[2]);
00101 
00102       double X3ppp = 2.0*p3ppp.E()/m6J;
00103       if (X3ppp > 0.9) {
00104         vetoEvent;
00105       }
00106 
00107       FourMomentum pAV = cms_boost.transform(_avg_beam_in_lab(m6J, jetsystem.rapidity()));
00108       double costheta3ppp = pAV.p3().unit().dot(p3ppp.p3().unit());
00109       if (fabs(costheta3ppp) > 0.9) {
00110         vetoEvent;
00111       }
00112 
00113       // 3-jet-system variables
00114       _h_m6J->fill(m6J, weight);
00115       _h_X3ppp->fill(X3ppp, weight);
00116       _h_X4ppp->fill(2.0*p4ppp.E()/m6J, weight);
00117       _h_costheta3ppp->fill(costheta3ppp, weight);
00118       double psi3ppp = _psi(p3ppp, pAV, p4ppp, p5ppp);
00119       _h_psi3ppp->fill(psi3ppp, weight);
00120       _h_f3ppp->fill(_safeMass(p3ppp)/m6J, weight);
00121       _h_f4ppp->fill(_safeMass(p4ppp)/m6J, weight);
00122       _h_f5ppp->fill(_safeMass(p5ppp)/m6J, weight);
00123 
00124       // 4 -> 3 jet variables
00125       _h_fApp->fill(_safeMass(pApp)/m6J, weight);
00126       _h_fBpp->fill(_safeMass(pApp)/m6J, weight);
00127       _h_XApp->fill(pApp.E()/(pApp.E()+pBpp.E()), weight);
00128       double psiAppBpp = _psi(pApp, pBpp, pApp+pBpp, pAV);
00129       _h_psiAppBpp->fill(psiAppBpp, weight);
00130 
00131       // 5 -> 4 jet variables
00132       _h_fCp->fill(_safeMass(pCp)/m6J, weight);
00133       _h_fDp->fill(_safeMass(pDp)/m6J, weight);
00134       _h_XCp->fill(pCp.E()/(pCp.E()+pDp.E()), weight);
00135       double psiCpDp = _psi(pCp, pDp, pCp+pDp, pAV);
00136       _h_psiCpDp->fill(psiCpDp, weight);
00137 
00138       // 6 -> 5 jet variables
00139       _h_fE->fill(_safeMass(pE)/m6J, weight);
00140       _h_fF->fill(_safeMass(pF)/m6J, weight);
00141       _h_XE->fill(pE.E()/(pE.E()+pF.E()), weight);
00142       double psiEF = _psi(pE, pF, pE+pF, pAV);
00143       _h_psiEF->fill(psiEF, weight);
00144     }
00145 
00146 
00147     void finalize() {
00148       normalize(_h_m6J);
00149       normalize(_h_X3ppp);
00150       normalize(_h_X4ppp);
00151       normalize(_h_costheta3ppp);
00152       normalize(_h_psi3ppp);
00153       normalize(_h_f3ppp);
00154       normalize(_h_f4ppp);
00155       normalize(_h_f5ppp);
00156       normalize(_h_XApp);
00157       normalize(_h_XCp);
00158       normalize(_h_XE);
00159       normalize(_h_psiAppBpp);
00160       normalize(_h_psiCpDp);
00161       normalize(_h_psiEF);
00162       normalize(_h_fApp);
00163       normalize(_h_fBpp);
00164       normalize(_h_fCp);
00165       normalize(_h_fDp);
00166       normalize(_h_fE);
00167       normalize(_h_fF);
00168     }
00169 
00170 
00171   private:
00172 
00173     vector<FourMomentum> _reduce(const vector<FourMomentum>& jets,
00174                                  FourMomentum& combined1,
00175                                  FourMomentum& combined2) {
00176       double minMass2 = 1e9;
00177       size_t idx1(jets.size()), idx2(jets.size());
00178       for (size_t i = 0; i < jets.size(); ++i) {
00179         for (size_t j = i+1; j < jets.size(); ++j) {
00180           double mass2 = FourMomentum(jets[i] + jets[j]).mass2();
00181           if (mass2 < minMass2) {
00182             idx1 = i;
00183             idx2 = j;
00184           }
00185         }
00186       }
00187       vector<FourMomentum> newjets;
00188       for (size_t i = 0; i < jets.size(); ++i) {
00189         if (i != idx1 && i != idx2) newjets.push_back(jets[i]);
00190       }
00191       newjets.push_back(jets[idx1] + jets[idx2]);
00192       combined1 = jets[idx1];
00193       combined2 = jets[idx2];
00194       return newjets;
00195     }
00196 
00197 
00198     FourMomentum _avg_beam_in_lab(const double& m, const double& y) {
00199       const double mt = m/2.0;
00200       FourMomentum beam1(mt, 0, 0,  mt);
00201       FourMomentum beam2(mt, 0, 0, -mt);
00202       if (fabs(y) > 1e-3) {
00203         FourMomentum boostvec(cosh(y), 0.0, 0.0, sinh(y));
00204         LorentzTransform cms_boost(-boostvec.boostVector());
00205         cms_boost = cms_boost.inverse();
00206         beam1 = cms_boost.transform(beam1);
00207         beam2 = cms_boost.transform(beam2);
00208       }
00209       if (beam1.E() > beam2.E()) {
00210         return beam1 - beam2;
00211       }
00212       else {
00213         return beam2 - beam1;
00214       }
00215     }
00216 
00217 
00218     double _psi(const FourMomentum& p1, const FourMomentum& p2,
00219                 const FourMomentum& p3, const FourMomentum& p4) {
00220       Vector3 p1xp2 = p1.p3().cross(p2.p3());
00221       Vector3 p3xp4 = p3.p3().cross(p4.p3());
00222       return mapAngle0ToPi(acos(p1xp2.unit().dot(p3xp4.unit())));
00223     }
00224 
00225 
00226     double _safeMass(const FourMomentum& p) {
00227       double mass2 = p.mass2();
00228       if (mass2 > 0.0) return sqrt(mass2);
00229       if (mass2 < -1e-5) MSG_WARNING("m2 = " << m2 << ". Assuming m2=0.");
00230       return 0.0;
00231     }
00232 
00233 
00234   private:
00235 
00236     Histo1DPtr _h_m6J;
00237     Histo1DPtr _h_X3ppp, _h_X4ppp;
00238     Histo1DPtr _h_costheta3ppp;
00239     Histo1DPtr _h_psi3ppp;
00240     Histo1DPtr _h_f3ppp, _h_f4ppp, _h_f5ppp;
00241     Histo1DPtr _h_XApp, _h_XCp, _h_XE;
00242     Histo1DPtr _h_psiAppBpp, _h_psiCpDp, _h_psiEF;
00243     Histo1DPtr _h_fApp, _h_fBpp, _h_fCp, _h_fDp, _h_fE, _h_fF;
00244 
00245   };
00246 
00247 
00248 
00249   // The hook for the plugin system
00250   DECLARE_RIVET_PLUGIN(CDF_1997_S3541940);
00251 
00252 }