ATLAS_2013_I1243871.cc
Go to the documentation of this file.
00001 // -*- C++ -*- 00002 #include "Rivet/Analysis.hh" 00003 #include "Rivet/Tools/Logging.hh" 00004 #include "Rivet/Projections/FinalState.hh" 00005 #include "Rivet/Projections/IdentifiedFinalState.hh" 00006 #include "Rivet/Projections/VetoedFinalState.hh" 00007 #include "Rivet/Projections/FastJets.hh" 00008 #include "Rivet/Tools/ParticleIdUtils.hh" 00009 #include "Rivet/Particle.hh" 00010 00011 namespace Rivet { 00012 00013 00014 class ATLAS_2013_I1243871 : public Analysis { 00015 public: 00016 00017 /// Constructor 00018 ATLAS_2013_I1243871() 00019 : Analysis("ATLAS_2013_I1243871") 00020 { } 00021 00022 00023 /// Book histograms and initialise projections before the run 00024 void init() { 00025 // Set up projections 00026 const FinalState fs(-4.5, 4.5); 00027 addProjection(fs, "ALL_FS"); 00028 00029 /// Get electrons from truth record 00030 IdentifiedFinalState elec_fs(Cuts::abseta < 2.47 && Cuts::pT > 25*GeV); 00031 elec_fs.acceptIdPair(PID::ELECTRON); 00032 addProjection(elec_fs, "ELEC_FS"); 00033 00034 /// Get muons which pass the initial kinematic cuts: 00035 IdentifiedFinalState muon_fs(Cuts::abseta < 2.5 && Cuts::pT > 20*GeV); 00036 muon_fs.acceptIdPair(PID::MUON); 00037 addProjection(muon_fs, "MUON_FS"); 00038 00039 // Final state used as input for jet-finding. 00040 // We include everything except the muons and neutrinos 00041 VetoedFinalState jet_input(fs); 00042 jet_input.vetoNeutrinos(); 00043 jet_input.addVetoPairId(PID::MUON); 00044 addProjection(jet_input, "JET_INPUT"); 00045 00046 // Get the jets 00047 FastJets jets(jet_input, FastJets::ANTIKT, 0.4); 00048 addProjection(jets, "JETS"); 00049 00050 // Book histograms 00051 for (size_t d = 0; d < 5; ++d) { 00052 _p_b_rho[d] = bookProfile1D(d+1, 1, 1); 00053 _p_l_rho[d] = bookProfile1D(d+1, 2, 1); 00054 _p_b_Psi[d] = bookProfile1D(d+1, 1, 2); 00055 _p_l_Psi[d] = bookProfile1D(d+1, 2, 2); 00056 } 00057 } 00058 00059 00060 /// Perform the per-event analysis 00061 void analyze(const Event& event) { 00062 const double weight = event.weight(); 00063 00064 /// Get the various sets of final state particles 00065 const ParticleVector& elecFS = applyProjection<IdentifiedFinalState>(event, "ELEC_FS").particlesByPt(); 00066 const ParticleVector& muonFS = applyProjection<IdentifiedFinalState>(event, "MUON_FS").particlesByPt(); 00067 00068 // Get all jets with pT > 7 GeV (ATLAS standard jet collection) 00069 /// @todo Why rewrite the jets collection as a vector of pointers? 00070 const Jets& jets = applyProjection<FastJets>(event, "JETS").jetsByPt(7*GeV); 00071 vector<const Jet*> allJets; 00072 foreach(const Jet& j, jets) allJets.push_back(&j); 00073 00074 // Keep any jets that pass the pt cut 00075 vector<const Jet*> pt_jets; 00076 foreach (const Jet* j, allJets) { 00077 /// @todo Use direct kinematics access 00078 const double pt = j->momentum().pT(); 00079 const double eta = j->momentum().eta(); 00080 if (pt > 25*GeV && fabs(eta) < 2.5) pt_jets.push_back(j); 00081 } 00082 00083 // Remove jets too close to an electron 00084 vector<const Jet*> good_jets; 00085 foreach (const Jet* j, pt_jets) { 00086 bool isElectron = 0; 00087 foreach (const Particle& e, elecFS) { 00088 const double elec_jet_dR = deltaR(e.momentum(), j->momentum()); 00089 if (elec_jet_dR < 0.2) { isElectron = true; break; } 00090 } 00091 if (!isElectron) good_jets.push_back(j); 00092 } 00093 00094 // Classify the event type 00095 const size_t nElec = elecFS.size(); 00096 const size_t nMuon = muonFS.size(); 00097 bool isSemilepton = false, isDilepton = false; 00098 if (nElec == 1 && nMuon == 0) { 00099 isSemilepton = true; 00100 } else if (nElec == 0 && nMuon == 1) { 00101 isSemilepton = true; 00102 } else if (nElec == 2 && nMuon == 0) { 00103 if (charge(elecFS[0]) != charge(elecFS[1])) isDilepton = true; 00104 } else if (nElec == 1 && nMuon == 1) { 00105 if (charge(elecFS[0]) != charge(muonFS[0])) isDilepton = true; 00106 } else if (nElec == 0 && nMuon == 2) { 00107 if (charge(muonFS[0]) != charge(muonFS[1])) isDilepton = true; 00108 } 00109 const bool isGoodEvent = (isSemilepton && good_jets.size() >= 4) || (isDilepton && good_jets.size() >= 2); 00110 if (!isGoodEvent) vetoEvent; 00111 00112 // Select b-hadrons 00113 /// @todo Use built-in identification on Particle, avoid HepMC 00114 vector<const GenParticle *> b_hadrons; 00115 vector<const GenParticle *> allParticles = particles(event.genEvent()); 00116 for (size_t i = 0; i < allParticles.size(); i++) { 00117 const GenParticle* p = allParticles.at(i); 00118 if ( !(PID::isHadron( p->pdg_id() ) && PID::hasBottom( p->pdg_id() )) ) continue; 00119 if (p->momentum().perp() < 5*GeV) continue; 00120 b_hadrons.push_back(p); 00121 } 00122 00123 // Select b-jets as those containing a b-hadron 00124 /// @todo Use built-in dR < 0.3 Jet tagging, avoid HepMC 00125 vector<const Jet*> b_jets; 00126 foreach (const Jet* j, good_jets) { 00127 bool isbJet = false; 00128 foreach (const GenParticle* b, b_hadrons) { 00129 /// @todo Use direct momentum accessor / delta functions 00130 const FourMomentum hadron = b->momentum(); 00131 const double hadron_jet_dR = deltaR(j->momentum(), hadron); 00132 if (hadron_jet_dR < 0.3) { isbJet = true; break; } 00133 } 00134 // Check if it is overlapped to any other jet 00135 bool isOverlapped = false; 00136 foreach (const Jet* k, allJets) { 00137 if (j == k) continue; 00138 double dRjj = deltaR(j->momentum(), k->momentum()); 00139 if (dRjj < 0.8) { isOverlapped = true; break; } 00140 } 00141 if (isbJet && !isOverlapped) b_jets.push_back(j); 00142 } 00143 MSG_DEBUG(b_jets.size() << " b-jets selected"); 00144 00145 00146 // Select light-jets as the pair of non-b-jets with invariant mass closest to the W mass 00147 /// @todo Use built-in b-tagging (dR < 0.3 defn), avoid HepMC 00148 const double nominalW = 80.4*GeV; 00149 double deltaM = 500*GeV; 00150 const Jet* light1 = NULL; const Jet* light2 = NULL; // NB: const Jets, not const pointers! 00151 foreach (const Jet* i, good_jets) { 00152 bool isbJet1 = false; 00153 foreach (const GenParticle* b, b_hadrons) { 00154 /// @todo Use direct momentum accessor / delta functions 00155 const FourMomentum hadron = b->momentum(); 00156 const double hadron_jet_dR = deltaR(i->momentum(), hadron); 00157 if (hadron_jet_dR < 0.3) { isbJet1 = true; break; } 00158 } 00159 if (isbJet1) continue; 00160 foreach (const Jet* j, good_jets) { 00161 bool isbJet2 = false; 00162 foreach (const GenParticle* b, b_hadrons) { 00163 FourMomentum hadron = b->momentum(); 00164 double hadron_jet_dR = deltaR(j->momentum(), hadron); 00165 if (hadron_jet_dR < 0.3) { isbJet2 = true; break; } 00166 } 00167 if (isbJet2) continue; 00168 double invMass = (i->momentum()+j->momentum()).mass(); 00169 if (fabs(invMass-nominalW) < deltaM){ 00170 deltaM = fabs(invMass - nominalW); 00171 light1 = i; 00172 light2 = j; 00173 } 00174 } 00175 } 00176 00177 // Check that both jets are not overlapped, and populate the light jets list 00178 vector<const Jet*> light_jets; 00179 const bool hasGoodLight = light1 != NULL && light2 != NULL && light1 != light2; 00180 if (hasGoodLight) { 00181 bool isOverlap1 = false, isOverlap2 = false; 00182 foreach (const Jet* j, allJets) { 00183 if (light1 == j) continue; 00184 const double dR1j = deltaR(light1->momentum(), j->momentum()); 00185 if (dR1j < 0.8) { isOverlap1 = true; break; } 00186 } 00187 foreach (const Jet* j, allJets) { 00188 if (light2 == j) continue; 00189 const double dR2j = deltaR(light2->momentum(), j->momentum()); 00190 if (dR2j < 0.8) { isOverlap2 = true; break; } 00191 } 00192 if (!isOverlap1 && !isOverlap2) { 00193 light_jets.push_back(light1); 00194 light_jets.push_back(light2); 00195 } 00196 } 00197 MSG_DEBUG(light_jets.size() << " light jets selected"); 00198 00199 00200 // Calculate the jet shapes 00201 /// @todo Use C++11 vector/array initialization 00202 const double binWidth = 0.04; // -> 10 bins from 0.0-0.4 00203 vector<double> ptEdges; ptEdges += 30, 40, 50, 70, 100, 150; 00204 00205 // b-jet shapes 00206 MSG_DEBUG("Filling b-jet shapes"); 00207 foreach (const Jet* bJet, b_jets) { 00208 // Work out jet pT bin and skip this jet if out of range 00209 const double jetPt = bJet->momentum().pT(); 00210 MSG_DEBUG("Jet pT = " << jetPt/GeV << " GeV"); 00211 if (!inRange(jetPt/GeV, 30., 150.)) continue; 00212 /// @todo Use YODA bin index lookup tools 00213 size_t ipt; for (ipt = 0; ipt < 5; ++ipt) if (inRange(jetPt/GeV, ptEdges[ipt], ptEdges[ipt+1])) break; 00214 MSG_DEBUG("Jet pT index = " << ipt); 00215 00216 // Calculate jet shape 00217 vector<double> rings(10, 0); 00218 foreach (const Particle& p, bJet->particles()) { 00219 const double dR = deltaR(bJet->momentum(), p.momentum()); 00220 const size_t idR = (size_t) floor(dR/binWidth); 00221 for (size_t i = idR; i < 10; ++i) rings[i] += p.pT(); 00222 } 00223 00224 // Fill each dR bin of the histos for this jet pT 00225 for (int iBin = 0; iBin < 10; ++iBin) { 00226 const double rcenter = 0.02 + iBin*binWidth; 00227 const double rhoval = (iBin != 0 ? (rings[iBin]-rings[iBin-1]) : rings[iBin]) / binWidth / rings[9]; 00228 const double psival = rings[iBin] / rings[9]; 00229 MSG_DEBUG(rcenter << ", " << rhoval << ", " << psival); 00230 _p_b_rho[ipt]->fill(rcenter, rhoval, weight); 00231 _p_b_Psi[ipt]->fill(rcenter, psival, weight); 00232 } 00233 } 00234 00235 // Light jet shapes 00236 MSG_DEBUG("Filling light jet shapes"); 00237 foreach (const Jet* lJet, light_jets) { 00238 // Work out jet pT bin and skip this jet if out of range 00239 const double jetPt = lJet->momentum().pT(); 00240 MSG_DEBUG("Jet pT = " << jetPt/GeV << " GeV"); 00241 if (!inRange(jetPt/GeV, 30., 150.)) continue; 00242 /// @todo Use YODA bin index lookup tools 00243 size_t ipt; for (ipt = 0; ipt < 5; ++ipt) if (inRange(jetPt/GeV, ptEdges[ipt], ptEdges[ipt+1])) break; 00244 MSG_DEBUG("Jet pT index = " << ipt); 00245 00246 // Calculate jet shape 00247 vector<double> rings(10, 0); 00248 foreach (const Particle& p, lJet->particles()) { 00249 const double dR = deltaR(lJet->momentum(), p.momentum()); 00250 const size_t idR = (size_t) floor(dR/binWidth); 00251 for (size_t i = idR; i < 10; ++i) rings[i] += p.pT(); 00252 } 00253 00254 // Fill each dR bin of the histos for this jet pT 00255 for (int iBin = 0; iBin < 10; ++iBin) { 00256 const double rcenter = 0.02 + iBin*binWidth; 00257 const double rhoval = (iBin != 0 ? (rings[iBin]-rings[iBin-1]) : rings[iBin]) / binWidth / rings[9]; 00258 const double psival = rings[iBin] / rings[9]; 00259 _p_l_rho[ipt]->fill(rcenter, rhoval, weight); 00260 _p_l_Psi[ipt]->fill(rcenter, psival, weight); 00261 } 00262 } 00263 00264 } 00265 00266 00267 private: 00268 00269 Profile1DPtr _p_b_rho[5]; 00270 Profile1DPtr _p_l_rho[5]; 00271 Profile1DPtr _p_b_Psi[5]; 00272 Profile1DPtr _p_l_Psi[5]; 00273 00274 }; 00275 00276 00277 // The hook for the plugin system 00278 DECLARE_RIVET_PLUGIN(ATLAS_2013_I1243871); 00279 00280 } Generated on Wed Oct 7 2015 12:09:11 for The Rivet MC analysis system by ![]() |