ATLAS_2015_I1345452.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/PromptFinalState.hh" 00007 #include "Rivet/Projections/DressedLeptons.hh" 00008 #include "Rivet/Projections/FastJets.hh" 00009 00010 namespace Rivet { 00011 00012 /// @brief ATLAS 7 TeV pseudo-top analysis 00013 /// 00014 /// @author K .Finelli <kevin.finelli@cern.ch> 00015 /// @author A. Saavedra <a.saavedra@physics.usyd.edu.au> 00016 /// @author L. Lan <llan@physics.usyd.edu.au> 00017 class ATLAS_2015_I1345452 : public Analysis { 00018 public: 00019 00020 /// Constructor 00021 DEFAULT_RIVET_ANALYSIS_CTOR(ATLAS_2015_I1345452); 00022 00023 00024 void init() { 00025 // Eta ranges 00026 Cut eta_full = (Cuts::abseta < 5.0) & (Cuts::pT >= 1.0*MeV); 00027 Cut eta_lep = (Cuts::abseta < 2.5); 00028 00029 // All final state particles 00030 FinalState fs(eta_full); 00031 00032 // Get photons to dress leptons 00033 IdentifiedFinalState photons(fs); 00034 photons.acceptIdPair(PID::PHOTON); 00035 00036 // Projection to find the electrons 00037 IdentifiedFinalState el_id(fs); 00038 el_id.acceptIdPair(PID::ELECTRON); 00039 00040 PromptFinalState electrons(el_id); 00041 electrons.acceptTauDecays(true); 00042 declare(electrons, "electrons"); 00043 00044 DressedLeptons dressedelectrons(photons, electrons, 0.1, eta_lep & (Cuts::pT >= 25.0*GeV), true, true); 00045 declare(dressedelectrons, "dressedelectrons"); 00046 00047 DressedLeptons ewdressedelectrons(photons, electrons, 0.1, eta_full, true, true); 00048 declare(ewdressedelectrons, "ewdressedelectrons"); 00049 00050 DressedLeptons vetodressedelectrons(photons, electrons, 0.1, eta_lep & (Cuts::pT >= 15.0*GeV), true, true); 00051 declare(vetodressedelectrons, "vetodressedelectrons"); 00052 00053 // Projection to find the muons 00054 IdentifiedFinalState mu_id(fs); 00055 mu_id.acceptIdPair(PID::MUON); 00056 PromptFinalState muons(mu_id); 00057 muons.acceptTauDecays(true); 00058 declare(muons, "muons"); 00059 DressedLeptons dressedmuons(photons, muons, 0.1, eta_lep & (Cuts::pT >= 25.0*GeV), true, true); 00060 declare(dressedmuons, "dressedmuons"); 00061 DressedLeptons ewdressedmuons(photons, muons, 0.1, eta_full, true, true); 00062 declare(ewdressedmuons, "ewdressedmuons"); 00063 DressedLeptons vetodressedmuons(photons, muons, 0.1, eta_lep & (Cuts::pT >= 15.0*GeV), true, true); 00064 declare(vetodressedmuons, "vetodressedmuons"); 00065 00066 // Projection to find neutrinos and produce MET 00067 IdentifiedFinalState nu_id; 00068 nu_id.acceptNeutrinos(); 00069 PromptFinalState neutrinos(nu_id); 00070 neutrinos.acceptTauDecays(true); 00071 declare(neutrinos, "neutrinos"); 00072 00073 // Jet clustering. 00074 VetoedFinalState vfs; 00075 vfs.addVetoOnThisFinalState(ewdressedelectrons); 00076 vfs.addVetoOnThisFinalState(ewdressedmuons); 00077 vfs.addVetoOnThisFinalState(neutrinos); 00078 FastJets jets(vfs, FastJets::ANTIKT, 0.4); 00079 jets.useInvisibles(); 00080 declare(jets, "jets"); 00081 00082 //pseudotop leptons and hadrons 00083 _h["ptpseudotophadron_mu"] = bookHisto1D( 1, 1, 2); 00084 _h["ptpseudotophadron_el"] = bookHisto1D( 2, 1, 2); 00085 _h["absrappseudotophadron_mu"] = bookHisto1D( 3, 1, 2); 00086 _h["absrappseudotophadron_el"] = bookHisto1D( 4, 1, 2); 00087 _h["ptpseudotoplepton_mu"] = bookHisto1D( 5, 1, 2); 00088 _h["ptpseudotoplepton_el"] = bookHisto1D( 6, 1, 2); 00089 _h["absrappseudotoplepton_mu"] = bookHisto1D( 7, 1, 2); 00090 _h["absrappseudotoplepton_el"] = bookHisto1D( 8, 1, 2); 00091 _h["ptttbar_mu"] = bookHisto1D( 9, 1, 2); 00092 _h["ptttbar_el"] = bookHisto1D(10, 1, 2); 00093 _h["absrapttbar_mu"] = bookHisto1D(11, 1, 2); 00094 _h["absrapttbar_el"] = bookHisto1D(12, 1, 2); 00095 _h["ttbarmass_mu"] = bookHisto1D(13, 1, 2); 00096 _h["ttbarmass_el"] = bookHisto1D(14, 1, 2); 00097 _h["ptpseudotophadron"] = bookHisto1D(15, 1, 2); 00098 _h["absrappseudotophadron"] = bookHisto1D(16, 1, 2); 00099 _h["ptpseudotoplepton"] = bookHisto1D(17, 1, 2); 00100 _h["absrappseudotoplepton"] = bookHisto1D(18, 1, 2); 00101 _h["ptttbar"] = bookHisto1D(19, 1, 2); 00102 _h["absrapttbar"] = bookHisto1D(20, 1, 2); 00103 _h["ttbarmass"] = bookHisto1D(21, 1, 2); 00104 00105 } 00106 00107 void analyze(const Event& event) { 00108 00109 // Get the selected objects, using the projections. 00110 _dressedelectrons = apply<DressedLeptons>( event, "dressedelectrons").dressedLeptons(); 00111 _vetodressedelectrons = apply<DressedLeptons>( event, "vetodressedelectrons").dressedLeptons(); 00112 _dressedmuons = apply<DressedLeptons>( event, "dressedmuons").dressedLeptons(); 00113 _vetodressedmuons = apply<DressedLeptons>( event, "vetodressedmuons").dressedLeptons(); 00114 _neutrinos = apply<PromptFinalState>(event, "neutrinos").particlesByPt(); 00115 const Jets& all_jets = apply<FastJets>( event, "jets").jetsByPt(Cuts::pT > 25.0*GeV && Cuts::abseta < 2.5); 00116 00117 //get true l+jets events by removing events with more than 1 electron||muon neutrino 00118 unsigned int n_elmu_neutrinos = 0; 00119 foreach (const Particle p, _neutrinos) { 00120 if (p.abspid() == 12 || p.abspid() == 14) ++n_elmu_neutrinos; 00121 } 00122 if (n_elmu_neutrinos != 1) vetoEvent; 00123 00124 DressedLepton *lepton; 00125 if ( _dressedelectrons.size()) lepton = &_dressedelectrons[0]; 00126 else if (_dressedmuons.size()) lepton = &_dressedmuons[0]; 00127 else vetoEvent; 00128 00129 // Calculate the missing ET, using the prompt neutrinos only (really?) 00130 /// @todo Why not use MissingMomentum? 00131 FourMomentum met; 00132 foreach (const Particle& p, _neutrinos) met += p.momentum(); 00133 00134 //remove jets if they are within dR < 0.2 of lepton 00135 Jets jets; 00136 foreach(const Jet& jet, all_jets) { 00137 bool keep = true; 00138 foreach (const DressedLepton& el, _vetodressedelectrons) { 00139 keep &= deltaR(jet, el) >= 0.2; 00140 } 00141 if (keep) jets += jet; 00142 } 00143 00144 bool overlap = false; 00145 Jets bjets, lightjets; 00146 for (unsigned int i = 0; i < jets.size(); ++i) { 00147 const Jet& jet = jets[i]; 00148 foreach (const DressedLepton& el, _dressedelectrons) overlap |= deltaR(jet, el) < 0.4; 00149 foreach (const DressedLepton& mu, _dressedmuons) overlap |= deltaR(jet, mu) < 0.4; 00150 for (unsigned int j = i + 1; j < jets.size(); ++j) { 00151 overlap |= deltaR(jet, jets[j]) < 0.5; 00152 } 00153 //// Count the number of b-tags 00154 bool b_tagged = false; // This is closer to the 00155 Particles bTags = jet.bTags(); // analysis. Something 00156 foreach ( Particle b, bTags ) { // about ghost-associated 00157 b_tagged |= b.pT() > 5*GeV; // B-hadrons 00158 } // 00159 if ( b_tagged ) bjets += jet; 00160 else lightjets += jet; 00161 } 00162 00163 // remove events with object overlap 00164 if (overlap) vetoEvent; 00165 00166 if (bjets.size() < 2 || lightjets.size() < 2) vetoEvent; 00167 00168 FourMomentum pbjet1; //Momentum of bjet1 00169 FourMomentum pbjet2; //Momentum of bjet2 00170 if ( deltaR(bjets[0], *lepton) <= deltaR(bjets[1], *lepton) ) { 00171 pbjet1 = bjets[0].momentum(); 00172 pbjet2 = bjets[1].momentum(); 00173 } else { 00174 pbjet1 = bjets[1].momentum(); 00175 pbjet2 = bjets[0].momentum(); 00176 } 00177 00178 FourMomentum pjet1; // Momentum of jet1 00179 if (lightjets.size()) pjet1 = lightjets[0].momentum(); 00180 00181 FourMomentum pjet2; // Momentum of jet 2 00182 if (lightjets.size() > 1) pjet2 = lightjets[1].momentum(); 00183 00184 double pz = computeneutrinoz(lepton->momentum(), met); 00185 FourMomentum ppseudoneutrino( sqrt(sqr(met.px()) + sqr(met.py()) + sqr(pz)), met.px(), met.py(), pz); 00186 00187 //compute leptonic, hadronic, combined pseudo-top 00188 FourMomentum ppseudotoplepton = lepton->momentum() + ppseudoneutrino + pbjet1; 00189 FourMomentum ppseudotophadron = pbjet2 + pjet1 + pjet2; 00190 FourMomentum pttbar = ppseudotoplepton + ppseudotophadron; 00191 00192 // Evaluate basic event selection 00193 bool pass_eljets = (_dressedelectrons.size() == 1) && 00194 (_vetodressedelectrons.size() < 2) && 00195 (_vetodressedmuons.empty()) && 00196 (met.pT() > 30*GeV) && 00197 (_mT(_dressedelectrons[0].momentum(), met) > 35*GeV) && 00198 (jets.size() >= 4); 00199 bool pass_mujets = (_dressedmuons.size() == 1) && 00200 (_vetodressedmuons.size() < 2) && 00201 (_vetodressedelectrons.empty()) && 00202 (met.pT() > 30*GeV) && 00203 (_mT(_dressedmuons[0].momentum(), met) > 35*GeV) && 00204 (jets.size() >= 4); 00205 00206 // basic event selection requirements 00207 if (!pass_eljets && !pass_mujets) vetoEvent; 00208 00209 // Fill histograms 00210 const double weight = event.weight(); 00211 00212 //pseudotop hadrons and leptons fill histogram 00213 _h["ptpseudotoplepton"]->fill( ppseudotoplepton.pt(), weight); //pT of pseudo top lepton 00214 _h["absrappseudotoplepton"]->fill(ppseudotoplepton.absrap(), weight); //absolute rapidity of pseudo top lepton 00215 _h["ptpseudotophadron"]->fill( ppseudotophadron.pt(), weight); //pT of pseudo top hadron 00216 _h["absrappseudotophadron"]->fill(ppseudotophadron.absrap(), weight); //absolute rapidity of pseudo top hadron 00217 _h["absrapttbar"]->fill( pttbar.absrap(), weight); //absolute rapidity of ttbar 00218 _h["ttbarmass"]->fill( pttbar.mass(), weight); //mass of ttbar 00219 _h["ptttbar"]->fill( pttbar.pt(), weight); //fill pT of ttbar in combined channel 00220 00221 if (pass_eljets) { // electron channel fill histogram 00222 _h["ptpseudotoplepton_el"]->fill( ppseudotoplepton.pt(), weight); //pT of pseudo top lepton 00223 _h["absrappseudotoplepton_el"]->fill(ppseudotoplepton.absrap(), weight); //absolute rapidity of pseudo top lepton 00224 _h["ptpseudotophadron_el"]->fill( ppseudotophadron.pt(), weight); //pT of pseudo top hadron 00225 _h["absrappseudotophadron_el"]->fill(ppseudotophadron.absrap(), weight); //absolute rapidity of pseudo top hadron 00226 _h["absrapttbar_el"]->fill( pttbar.absrap(), weight); //absolute rapidity of ttbar 00227 _h["ttbarmass_el"]->fill( pttbar.mass(), weight); // fill electron channel ttbar mass 00228 _h["ptttbar_el"]->fill( pttbar.pt(), weight); //fill pT of ttbar in electron channel 00229 } 00230 else { // muon channel fill histogram 00231 _h["ptpseudotoplepton_mu"]->fill( ppseudotoplepton.pt(), weight); //pT of pseudo top lepton 00232 _h["absrappseudotoplepton_mu"]->fill(ppseudotoplepton.absrap(), weight); //absolute rapidity of pseudo top lepton 00233 _h["ptpseudotophadron_mu"]->fill( ppseudotophadron.pt(), weight); //pT of pseudo top hadron 00234 _h["absrappseudotophadron_mu"]->fill(ppseudotophadron.absrap(), weight); //absolute rapidity of pseudo top hadron 00235 _h["absrapttbar_mu"]->fill( pttbar.absrap(), weight); //absolute rapidity of ttbar 00236 _h["ttbarmass_mu"]->fill( pttbar.mass(), weight); //fill muon channel histograms 00237 _h["ptttbar_mu"]->fill( pttbar.pt(), weight); //fill pT of ttbar in electron channel 00238 } 00239 } 00240 00241 void finalize() { 00242 // Normalize to cross-section 00243 const double scalefactor(crossSection() / sumOfWeights()); 00244 for (map<string, Histo1DPtr>::iterator hit = _h.begin(); hit != _h.end(); ++hit) { 00245 double sf = scalefactor; 00246 if ( (hit->first).find("_") == std::string::npos ) sf *= 0.5; 00247 scale(hit->second, sf); 00248 } 00249 } 00250 00251 private: 00252 00253 00254 double computeneutrinoz(const FourMomentum& lepton, FourMomentum& met) const { 00255 //computing z component of neutrino momentum given lepton and met 00256 double pzneutrino; 00257 double m_W = 80.399; // in GeV, given in the paper 00258 double k = (( sqr( m_W ) - sqr( lepton.mass() ) ) / 2 ) + (lepton.px() * met.px() + lepton.py() * met.py()); 00259 double a = sqr ( lepton.E() )- sqr ( lepton.pz() ); 00260 double b = -2*k*lepton.pz(); 00261 double c = sqr( lepton.E() ) * sqr( met.pT() ) - sqr( k ); 00262 double discriminant = sqr(b) - 4 * a * c; 00263 double quad[2] = { (- b - sqrt(discriminant)) / (2 * a), (- b + sqrt(discriminant)) / (2 * a) }; //two possible quadratic solns 00264 if (discriminant < 0) pzneutrino = - b / (2 * a); //if the discriminant is negative 00265 else { //if the discriminant is greater than or equal to zero, take the soln with smallest absolute value 00266 double absquad[2]; 00267 for (int n=0; n<2; ++n) absquad[n] = fabs(quad[n]); 00268 if (absquad[0] < absquad[1]) pzneutrino = quad[0]; 00269 else pzneutrino = quad[1]; 00270 } 00271 if ( !std::isfinite(pzneutrino) ) std::cout << "Found non-finite value" << std::endl; 00272 return pzneutrino; 00273 } 00274 00275 double _mT(const FourMomentum &l, FourMomentum &nu) const { 00276 return sqrt( 2 * l.pT() * nu.pT() * (1 - cos(deltaPhi(l, nu))) ); 00277 } 00278 00279 /// @name Objects that are used by the event selection decisions 00280 vector<DressedLepton> _dressedelectrons, _vetodressedelectrons, _dressedmuons, _vetodressedmuons; 00281 Particles _neutrinos; 00282 map<string, Histo1DPtr> _h; 00283 }; 00284 00285 // The hook for the plugin system 00286 DECLARE_RIVET_PLUGIN(ATLAS_2015_I1345452); 00287 00288 } Generated on Tue Dec 13 2016 16:32:35 for The Rivet MC analysis system by ![]() |