ATLAS_2011_I930220.cc
Go to the documentation of this file.
00001 // -*- C++ -*- 00002 #include "Rivet/Analysis.hh" 00003 #include "Rivet/Projections/FastJets.hh" 00004 #include "Rivet/Projections/HeavyHadrons.hh" 00005 #include "Rivet/Tools/BinnedHistogram.hh" 00006 00007 namespace Rivet { 00008 00009 00010 /// @brief ATLAS inclusive b-jet pT spectrum, di-jet mass and di-jet chi 00011 class ATLAS_2011_I930220: public Analysis { 00012 public: 00013 00014 ATLAS_2011_I930220() 00015 : Analysis("ATLAS_2011_I930220") 00016 { } 00017 00018 00019 void init() { 00020 FinalState fs(-3.5, 3.5); 00021 addProjection(fs, "FinalState"); 00022 FastJets fj(fs, FastJets::ANTIKT, 0.4); 00023 fj.useInvisibles(); 00024 addProjection(fj, "Jets"); 00025 addProjection(HeavyHadrons(Cuts::abseta < 3.5 && Cuts::pT > 5*GeV), "BHadrons"); 00026 00027 double ybins[] = { 0.0, 0.3, 0.8, 1.2, 2.1 }; 00028 for (size_t i = 0; i < 4; ++i) 00029 _bjetpT_SV0.addHistogram(ybins[i], ybins[i+1], bookHisto1D(i+1, 1, 1)); 00030 00031 _bjetpT_SV0_All = bookHisto1D(5, 1, 1); 00032 _bjetpT_pTRel = bookHisto1D(6, 1, 1); 00033 _dijet_mass = bookHisto1D(7, 1, 1); 00034 _dijet_phi = bookHisto1D(8, 1, 1); 00035 _dijet_chi_110_370 = bookHisto1D(9, 1, 1); 00036 _dijet_chi_370_850 = bookHisto1D(10, 1, 1); 00037 00038 _chiCounter1 = 0.0; 00039 _chiCounter2 = 0.0; 00040 _phiCounter = 0.0; 00041 } 00042 00043 00044 void analyze(const Event& evt) { 00045 const double weight = evt.weight(); 00046 00047 const Particles& bHadrons = applyProjection<HeavyHadrons>(evt, "BHadrons").bHadrons(); 00048 const Jets& jets = applyProjection<JetAlg>(evt, "Jets").jetsByPt(15*GeV); 00049 00050 FourMomentum leadingJet, subleadingJet; 00051 int leadJet = 0, subJet = 0; 00052 foreach (const Jet& j, jets) { 00053 bool hasB = false; 00054 foreach (const Particle& b, bHadrons) 00055 if (deltaR(j, b) < 0.3) { hasB = true; break; } 00056 00057 // Identify and classify the leading and subleading jets 00058 if (j.absrap() < 2.1) { ///< Move this into the jets defn 00059 if (!leadJet) { 00060 leadingJet = j.momentum(); 00061 leadJet = (hasB && j.pT() > 40*GeV) ? 2 : 1; 00062 } 00063 else if (leadJet && !subJet) { 00064 subleadingJet = j.momentum(); 00065 subJet = (hasB && j.pT() > 40*GeV) ? 2 : 1; 00066 } 00067 if (hasB) { 00068 _bjetpT_SV0.fill(j.absrap(), j.pT()/GeV, weight); 00069 _bjetpT_SV0_All->fill(j.pT()/GeV, weight); 00070 _bjetpT_pTRel->fill(j.pT()/GeV, weight); 00071 } 00072 } 00073 } 00074 00075 // Di-b-jet plots require both the leading and subleading jets to be b-tagged and have pT > 40 GeV 00076 if (leadJet == 2 && subJet == 2) { 00077 const double mass = FourMomentum( leadingJet + subleadingJet ).mass(); 00078 _dijet_mass->fill(mass/GeV, weight); 00079 00080 // Plot dphi for high-mass di-b-jets 00081 if (mass > 110*GeV) { 00082 _phiCounter += weight; 00083 const double d_phi = deltaPhi( leadingJet.phi(), subleadingJet.phi() ); 00084 _dijet_phi->fill(fabs(d_phi), weight); 00085 } 00086 00087 // Plot chi for low y_boost di-b-jets (in two high-mass bins) 00088 const double y_boost = 0.5 * (leadingJet.rapidity() + subleadingJet.rapidity()); 00089 const double chi = exp( fabs( leadingJet.rapidity() - subleadingJet.rapidity() ) ); 00090 if ( fabs(y_boost) < 1.1 ) { 00091 if (inRange(mass/GeV, 110, 370)) { 00092 _chiCounter1 += weight; 00093 _dijet_chi_110_370->fill(chi, weight); 00094 } else if (inRange(mass/GeV, 370, 850)) { 00095 _chiCounter2 += weight; 00096 _dijet_chi_370_850->fill(chi, weight); 00097 } 00098 } 00099 } 00100 } 00101 00102 00103 void finalize() { 00104 // Normalizing to cross-section and mass 00105 // Additional factors represent the division by rapidity 00106 const double xsec = crossSectionPerEvent()/(picobarn); 00107 const double chiScale1 = 1 / _chiCounter1 / 260.0; 00108 const double chiScale2 = 1 / _chiCounter2 / 480.0; 00109 const double phiScale = 1 / _phiCounter; 00110 00111 _bjetpT_SV0.scale(xsec/2, this); 00112 scale(_bjetpT_SV0_All, xsec); 00113 scale(_bjetpT_pTRel, xsec); 00114 scale(_dijet_mass, xsec); 00115 scale(_dijet_phi, phiScale ); 00116 scale(_dijet_chi_110_370, chiScale1); 00117 scale(_dijet_chi_370_850, chiScale2); 00118 } 00119 00120 00121 private: 00122 00123 BinnedHistogram<double> _bjetpT_SV0; 00124 00125 Histo1DPtr _bjetpT_SV0_All; 00126 Histo1DPtr _bjetpT_pTRel; 00127 Histo1DPtr _dijet_mass; 00128 Histo1DPtr _dijet_phi; 00129 Histo1DPtr _dijet_chi_110_370; 00130 Histo1DPtr _dijet_chi_370_850; 00131 00132 double _chiCounter1; 00133 double _chiCounter2; 00134 double _phiCounter; 00135 }; 00136 00137 00138 // The hook for the plugin system 00139 DECLARE_RIVET_PLUGIN(ATLAS_2011_I930220); 00140 00141 } Generated on Tue Mar 24 2015 17:35:24 for The Rivet MC analysis system by ![]() |