CMS_2013_I1273574.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 /// CMS 4-jet production at 7 TeV 00010 class CMS_2013_I1273574 : public Analysis { 00011 public: 00012 00013 /// Constructor 00014 CMS_2013_I1273574() 00015 : Analysis("CMS_2013_I1273574") 00016 { } 00017 00018 00019 /// Book histograms and initialise projections before the run 00020 void init() { 00021 const FinalState cnfs(-4.7, 4.7); 00022 addProjection(FastJets(cnfs, FastJets::ANTIKT, 0.5), "Jets"); 00023 00024 _h_jetetas[0] = bookHisto1D(1,1,1); 00025 _h_jetpts[0] = bookHisto1D(2,1,1); 00026 _h_DeltaS = bookHisto1D(3,1,1); 00027 _h_DeltaPhiSoft = bookHisto1D(4,1,1); 00028 _h_DeltaPtRelSoft = bookHisto1D(5,1,1); 00029 _h_jetetas[2] = bookHisto1D(6,1,1); 00030 _h_jetpts[2] = bookHisto1D(7,1,1); 00031 _h_jetetas[3] = bookHisto1D(8,1,1); 00032 _h_jetpts[3] = bookHisto1D(9,1,1); 00033 _h_jetetas[1] = bookHisto1D(10,1,1); 00034 _h_jetpts[1] = bookHisto1D(11,1,1); 00035 } 00036 00037 00038 /// Perform the per-event analysis 00039 void analyze(const Event& event) { 00040 /// @todo Use jetsByPt(ptGtr(20*GeV) & absetaIn(4.7)), then no need for the lower loop; 00041 const Jets jets = applyProjection<FastJets>(event, "Jets").jetsByPt(20*GeV); 00042 if (jets.size() < 4) vetoEvent; 00043 00044 // Ensure that there are exactly 4 jets > 20 GeV, with two above 50 GeV 00045 Jets hardjets, alljets; 00046 foreach (const Jet& j, jets) { 00047 if (j.abseta() > 4.7) continue; 00048 if (j.pT() > 50*GeV) hardjets.push_back(j); 00049 if (j.pT() > 20*GeV) alljets.push_back(j); 00050 } 00051 if (hardjets.size() < 2 || alljets.size() != 4) vetoEvent; 00052 const double weight = event.weight(); 00053 00054 // Histogram pT and eta of all 4 jets 00055 for (size_t i = 0; i < 4; ++i) { 00056 _h_jetpts[i]->fill(alljets[i].pT()/GeV, weight); 00057 _h_jetetas[i]->fill(alljets[i].eta(), weight); 00058 } 00059 00060 // Create vector sums of the hard and soft pairs of jets 00061 const FourMomentum p12 = alljets[0].momentum() + alljets[1].momentum(); 00062 const FourMomentum p34 = alljets[2].momentum() + alljets[3].momentum(); 00063 00064 // Fill the delta(phi) between the soft jets 00065 const double dphisoft = deltaPhi(alljets[2], alljets[3]); 00066 _h_DeltaPhiSoft->fill(dphisoft, weight); 00067 00068 // Fill the pT balance between the soft jets 00069 const double ptbalanceSoft = p34.pT() / (alljets[2].pT() + alljets[3].pT()); 00070 _h_DeltaPtRelSoft->fill(ptbalanceSoft, weight); 00071 00072 // Fill the azimuthal angle difference between the two jet pairs 00073 const double p12p34_trans = p12.px()*p34.px() + p12.py()*p34.py(); 00074 const double DeltaS = acos( p12p34_trans / p12.pT() / p34.pT() ); 00075 _h_DeltaS->fill(DeltaS, weight); 00076 } 00077 00078 00079 /// Normalise histograms (mostly to cross-section) 00080 void finalize() { 00081 const double invlumi = crossSection()/picobarn/sumOfWeights(); 00082 for (size_t i = 0; i < 4; ++i) { 00083 scale(_h_jetpts[i], invlumi); 00084 scale(_h_jetetas[i], invlumi); 00085 } 00086 normalize(_h_DeltaPtRelSoft); 00087 normalize(_h_DeltaPhiSoft); 00088 normalize(_h_DeltaS); 00089 } 00090 00091 00092 private: 00093 00094 Histo1DPtr _h_jetpts[4], _h_jetetas[4]; 00095 Histo1DPtr _h_DeltaS, _h_DeltaPhiSoft, _h_DeltaPtRelSoft; 00096 00097 }; 00098 00099 00100 DECLARE_RIVET_PLUGIN(CMS_2013_I1273574); 00101 00102 } Generated on Tue Sep 30 2014 19:45:43 for The Rivet MC analysis system by ![]() |