ATLAS_2013_I1217867.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/VetoedFinalState.hh" 00005 #include "Rivet/Projections/IdentifiedFinalState.hh" 00006 #include "Rivet/Projections/DressedLeptons.hh" 00007 #include "Rivet/Projections/FastJets.hh" 00008 00009 namespace Rivet { 00010 00011 class ATLAS_2013_I1217867 : public Analysis { 00012 public: 00013 00014 /// @name Constructors etc. 00015 //@{ 00016 00017 /// Constructor 00018 ATLAS_2013_I1217867() 00019 : Analysis("ATLAS_2013_I1217867") 00020 { 00021 m_njet=4; 00022 _h_dI.resize(2, std::vector<Histo1DPtr>(m_njet)); 00023 _h_dI_ratio.resize(2, std::vector<Histo1DPtr>(m_njet-1)); 00024 } 00025 00026 //@} 00027 00028 00029 public: 00030 00031 /// @name Analysis methods 00032 //@{ 00033 00034 /// Book histograms and initialise projections before the run 00035 void init() { 00036 // Initialise projections 00037 00038 FinalState fs(-5.0, 5.0, 0.0*GeV); 00039 00040 IdentifiedFinalState bareElectrons(fs); 00041 bareElectrons.acceptIdPair(PID::ELECTRON); 00042 00043 Cut cuts = ( EtaIn(-2.47, -1.52) 00044 | EtaIn(-1.37, 1.37) 00045 | EtaIn( 1.52, 2.47) ) & (Cuts::pT >= 20.0*GeV); 00046 00047 DressedLeptons electronClusters(fs, bareElectrons, 0.1, true, cuts); 00048 addProjection(electronClusters, "electronClusters"); 00049 00050 IdentifiedFinalState bareMuons(fs); 00051 bareMuons.acceptIdPair(PID::MUON); 00052 Cut mucuts = EtaIn(-2.4,2.4) & (Cuts::pT >= 20.0*GeV); 00053 DressedLeptons muonClusters(fs, bareMuons, 0.1, true, mucuts); 00054 addProjection(muonClusters, "muonClusters"); 00055 00056 IdentifiedFinalState neutrinos(-MAXRAPIDITY, MAXRAPIDITY, 25.0*GeV); 00057 neutrinos.acceptNeutrinos(); 00058 addProjection(neutrinos, "neutrinos"); 00059 00060 VetoedFinalState jetFS(fs); 00061 jetFS.addVetoOnThisFinalState(electronClusters); 00062 jetFS.addVetoOnThisFinalState(muonClusters); 00063 jetFS.addVetoOnThisFinalState(neutrinos); 00064 FastJets jetpro(jetFS, FastJets::KT, 0.6); 00065 jetpro.useInvisibles(true); 00066 addProjection(jetpro, "jets"); 00067 00068 // Book histograms 00069 for (size_t flav=0; flav < 2; ++flav) { 00070 for (size_t i=0; i < m_njet; ++i) _h_dI[flav][i] = bookHisto1D(i+1, 1, flav+1); 00071 for (size_t i=0; i < m_njet-1; ++i) _h_dI_ratio[flav][i] = bookHisto1D(4+i+1, 1, flav+1); 00072 } 00073 } 00074 00075 00076 /// Perform the per-event analysis 00077 void analyze(const Event& e) { 00078 const double weight = e.weight(); 00079 00080 const DressedLeptons& electronClusters = applyProjection<DressedLeptons>(e, "electronClusters"); 00081 const DressedLeptons& muonClusters = applyProjection<DressedLeptons>(e, "muonClusters"); 00082 int ne = electronClusters.clusteredLeptons().size(); 00083 int nmu = muonClusters.clusteredLeptons().size(); 00084 00085 FourMomentum lepton; 00086 size_t flav = 2; 00087 if (ne==1) { 00088 lepton=electronClusters.clusteredLeptons()[0].momentum(); 00089 flav = 0; 00090 if (nmu > 0) vetoEvent; 00091 } 00092 else if (nmu == 1) { 00093 lepton=muonClusters.clusteredLeptons()[0].momentum(); 00094 flav = 1; 00095 if (ne > 0) vetoEvent; 00096 } 00097 else { 00098 vetoEvent; 00099 } 00100 00101 const Particles& neutrinos = applyProjection<FinalState>(e, "neutrinos").particlesByPt(); 00102 if (neutrinos.size() < 1) vetoEvent; 00103 FourMomentum neutrino = neutrinos[0].momentum(); 00104 00105 double mtW=sqrt(2.0*lepton.pT()*neutrino.pT()*(1-cos(lepton.phi()-neutrino.phi()))); 00106 if (mtW<40.0*GeV) vetoEvent; 00107 00108 const fastjet::ClusterSequence* seq = applyProjection<FastJets>(e, "jets").clusterSeq(); 00109 if (seq != NULL) { 00110 for (size_t i = 0; i < min(m_njet,(size_t)seq->n_particles()); ++i) { 00111 double d_ij = sqrt(seq->exclusive_dmerge_max(i)); 00112 _h_dI[flav][i]->fill(d_ij, weight); 00113 00114 if (i<m_njet-1) { 00115 if (d_ij>20.0*GeV) { 00116 double d_ijplus1 = sqrt(seq->exclusive_dmerge_max(i+1)); 00117 _h_dI_ratio[flav][i]->fill(d_ijplus1/d_ij, weight); 00118 } 00119 } 00120 } 00121 } 00122 00123 } 00124 00125 00126 /// Normalise histograms etc., after the run 00127 void finalize() { 00128 00129 for (size_t flav = 0; flav < 2; ++flav) { 00130 for (size_t i = 0; i < m_njet; ++i) { 00131 normalize(_h_dI[flav][i], 1.0, false); 00132 if (i < m_njet-1) normalize(_h_dI_ratio[flav][i], 1.0, false); 00133 } 00134 } 00135 } 00136 00137 //@} 00138 00139 00140 private: 00141 00142 /// @name Histograms 00143 //@{ 00144 std::vector<std::vector<Histo1DPtr> > _h_dI; 00145 std::vector<std::vector<Histo1DPtr> > _h_dI_ratio; 00146 00147 //@} 00148 00149 size_t m_njet; 00150 }; 00151 00152 00153 00154 // The hook for the plugin system 00155 DECLARE_RIVET_PLUGIN(ATLAS_2013_I1217867); 00156 00157 00158 } Generated on Thu Feb 6 2014 17:38:42 for The Rivet MC analysis system by ![]() |