MC_HHJETS.cc
Go to the documentation of this file.
00001 // -*- C++ -*- 00002 #include "Rivet/Analyses/MC_JetAnalysis.hh" 00003 #include "Rivet/Projections/ZFinder.hh" 00004 #include "Rivet/Projections/VetoedFinalState.hh" 00005 #include "Rivet/Projections/IdentifiedFinalState.hh" 00006 #include "Rivet/Projections/FastJets.hh" 00007 00008 namespace Rivet { 00009 00010 /// @brief MC validation analysis for higgs pairs events (stable Higgses) 00011 class MC_HHJETS : public MC_JetAnalysis { 00012 public: 00013 00014 /// Default constructor 00015 MC_HHJETS() 00016 : MC_JetAnalysis("MC_HHJETS", 4, "Jets") 00017 { } 00018 00019 00020 /// @name Analysis methods 00021 //@{ 00022 00023 /// Book histograms 00024 void init() { 00025 IdentifiedFinalState ifs(Cuts::abseta < 10.0 && Cuts::pT > 0*GeV); 00026 ifs.acceptId(25); 00027 declare(ifs,"IFS"); 00028 00029 VetoedFinalState vfs; 00030 vfs.addVetoPairId(25); 00031 declare(FastJets(vfs, FastJets::ANTIKT, 0.4), "Jets"); 00032 00033 _h_HH_mass = bookHisto1D("HH_mass", 250, 240, 4000.0); 00034 _h_HH_dR = bookHisto1D("HH_dR", 25, 0.5, 10.0); 00035 _h_HH_dPhi = bookHisto1D("HH_dPhi", 64, 0, 3.2); 00036 _h_HH_deta= bookHisto1D("HH_deta", 50, -5, 5); 00037 _h_H_pT = bookHisto1D("H_pT", 50, 0, 2000.0); 00038 _h_HH_pT = bookHisto1D("HH_pT", 200, 0, 2000.0); 00039 _h_H_pT1 = bookHisto1D("H_pT1", 200, 0, 2000.0); 00040 _h_H_pT2 = bookHisto1D("H_pT2", 200, 0, 2000.0); 00041 _h_H_eta = bookHisto1D("H_eta", 50, -5.0, 5.0); 00042 _h_H_eta1 = bookHisto1D("H_eta1", 50, -5.0, 5.0); 00043 _h_H_eta2 = bookHisto1D("H_eta2", 50, -5.0, 5.0); 00044 _h_H_phi = bookHisto1D("H_phi", 25, 0.0, TWOPI); 00045 _h_H_jet1_deta = bookHisto1D("H_jet1_deta", 50, -5.0, 5.0); 00046 _h_H_jet1_dR = bookHisto1D("H_jet1_dR", 25, 0.5, 7.0); 00047 00048 MC_JetAnalysis::init(); 00049 } 00050 00051 00052 00053 /// Do the analysis 00054 void analyze(const Event & e) { 00055 00056 const IdentifiedFinalState& ifs = apply<IdentifiedFinalState>(e, "IFS"); 00057 Particles allp = ifs.particlesByPt(); 00058 if (allp.empty()) vetoEvent; 00059 00060 const double weight = e.weight(); 00061 00062 FourMomentum hmom = allp[0].momentum(); 00063 if (allp.size() > 1) { 00064 FourMomentum hmom2(allp[1].momentum()); 00065 _h_HH_dR->fill(deltaR(hmom, hmom2), weight); 00066 _h_HH_dPhi->fill(deltaPhi(hmom, hmom2), weight); 00067 _h_HH_deta->fill(hmom.eta()-hmom2.eta(), weight); 00068 _h_HH_pT->fill((hmom+hmom2).pT(), weight); 00069 _h_HH_mass->fill((hmom+hmom2).mass(), weight); 00070 00071 if (hmom.pT() > hmom2.pT()) { 00072 _h_H_pT1->fill(hmom.pT(), weight); 00073 _h_H_eta1->fill(hmom.eta(), weight); 00074 _h_H_pT2->fill(hmom2.pT(), weight); 00075 _h_H_eta2->fill(hmom2.eta(), weight); 00076 } else { 00077 _h_H_pT1->fill(hmom2.pT(), weight); 00078 _h_H_eta1->fill(hmom2.eta(), weight); 00079 _h_H_pT2->fill(hmom.pT(), weight); 00080 _h_H_eta2->fill(hmom.eta(), weight); 00081 } 00082 } 00083 _h_H_pT->fill(hmom.pT(), weight); 00084 _h_H_eta->fill(hmom.eta(), weight); 00085 _h_H_phi->fill(hmom.azimuthalAngle(), weight); 00086 00087 00088 // Get the jet candidates 00089 Jets jets = apply<FastJets>(e, "Jets").jetsByPt(20.0*GeV); 00090 if (!jets.empty()) { 00091 _h_H_jet1_deta->fill(deltaEta(hmom, jets[0]), weight); 00092 _h_H_jet1_dR->fill(deltaR(hmom, jets[0]), weight); 00093 } 00094 00095 MC_JetAnalysis::analyze(e); 00096 } 00097 00098 00099 /// Finalize 00100 void finalize() { 00101 scale(_h_HH_mass, crossSection()/sumOfWeights()); 00102 scale(_h_HH_dR, crossSection()/sumOfWeights()); 00103 scale(_h_HH_deta, crossSection()/sumOfWeights()); 00104 scale(_h_HH_dPhi, crossSection()/sumOfWeights()); 00105 scale(_h_H_pT, crossSection()/sumOfWeights()); 00106 scale(_h_H_pT1, crossSection()/sumOfWeights()); 00107 scale(_h_H_pT2, crossSection()/sumOfWeights()); 00108 scale(_h_HH_pT, crossSection()/sumOfWeights()); 00109 scale(_h_H_eta, crossSection()/sumOfWeights()); 00110 scale(_h_H_eta1, crossSection()/sumOfWeights()); 00111 scale(_h_H_eta2, crossSection()/sumOfWeights()); 00112 scale(_h_H_phi, crossSection()/sumOfWeights()); 00113 scale(_h_H_jet1_deta, crossSection()/sumOfWeights()); 00114 scale(_h_H_jet1_dR, crossSection()/sumOfWeights()); 00115 00116 MC_JetAnalysis::finalize(); 00117 } 00118 00119 //@} 00120 00121 00122 private: 00123 00124 /// @name Histograms 00125 //@{ 00126 Histo1DPtr _h_HH_mass, _h_HH_pT, _h_HH_dR, _h_HH_deta, _h_HH_dPhi; 00127 Histo1DPtr _h_H_pT, _h_H_pT1, _h_H_pT2, _h_H_eta, _h_H_eta1, _h_H_eta2, _h_H_phi; 00128 Histo1DPtr _h_H_jet1_deta, _h_H_jet1_dR; 00129 //@} 00130 00131 }; 00132 00133 00134 // The hook for the plugin system 00135 DECLARE_RIVET_PLUGIN(MC_HHJETS); 00136 00137 } Generated on Tue Dec 13 2016 16:32:38 for The Rivet MC analysis system by ![]() |