Go to the documentation of this file.00001 #include "PhysicsTools/JetMCUtils/interface/combination.h"
00002
00003 #include "TopQuarkAnalysis/TopEventSelection/plugins/TtSemiLepSignalSelMVAComputer.h"
00004 #include "TopQuarkAnalysis/TopEventSelection/interface/TtSemiLepSignalSelEval.h"
00005
00006 #include "DataFormats/RecoCandidate/interface/RecoCandidate.h"
00007 #include "DataFormats/PatCandidates/interface/Jet.h"
00008 #include "DataFormats/PatCandidates/interface/MET.h"
00009 #include "DataFormats/PatCandidates/interface/Muon.h"
00010 #include "DataFormats/PatCandidates/interface/Electron.h"
00011
00012 #include "DataFormats/Common/interface/TriggerResults.h"
00013 #include "FWCore/Framework/interface/TriggerNamesService.h"
00014 #include "FWCore/ServiceRegistry/interface/Service.h"
00015 #include "DataFormats/PatCandidates/interface/Flags.h"
00016
00017
00018 TtSemiLepSignalSelMVAComputer::TtSemiLepSignalSelMVAComputer(const edm::ParameterSet& cfg):
00019 muons_ (cfg.getParameter<edm::InputTag>("muons")),
00020 jets_ (cfg.getParameter<edm::InputTag>("jets")),
00021 METs_ (cfg.getParameter<edm::InputTag>("mets")),
00022 electrons_ (cfg.getParameter<edm::InputTag>("elecs"))
00023 {
00024 produces< double >("DiscSel");
00025 }
00026
00027
00028
00029 TtSemiLepSignalSelMVAComputer::~TtSemiLepSignalSelMVAComputer()
00030 {
00031 }
00032
00033 void
00034 TtSemiLepSignalSelMVAComputer::produce(edm::Event& evt, const edm::EventSetup& setup)
00035 {
00036 std::auto_ptr< double > pOutDisc (new double);
00037
00038 mvaComputer.update<TtSemiLepSignalSelMVARcd>(setup, "ttSemiLepSignalSelMVA");
00039
00040
00041
00042 edm::ESHandle<PhysicsTools::Calibration::MVAComputerContainer> calibContainer;
00043 setup.get<TtSemiLepSignalSelMVARcd>().get( calibContainer );
00044 std::vector<PhysicsTools::Calibration::VarProcessor*> processors
00045 = (calibContainer->find("ttSemiLepSignalSelMVA")).getProcessors();
00046
00047
00048 edm::Handle<edm::View<pat::MET> > MET_handle;
00049 evt.getByLabel(METs_,MET_handle);
00050 if(!MET_handle.isValid()) return;
00051 const edm::View<pat::MET> MET = *MET_handle;
00052
00053 edm::Handle< std::vector<pat::Jet> > jet_handle;
00054 evt.getByLabel(jets_, jet_handle);
00055 if(!jet_handle.isValid()) return;
00056 const std::vector<pat::Jet> jets = *jet_handle;
00057 unsigned int nJets = 0;
00058 std::vector<pat::Jet> seljets;
00059 for(std::vector<pat::Jet>::const_iterator it = jets.begin(); it != jets.end(); it++) {
00060 if(!(pat::Flags::test(*it, pat::Flags::Overlap::Electrons))) continue;
00061 if(it->pt()>30. && fabs(it->eta())<2.4) {
00062 seljets.push_back(*it);
00063 nJets++;
00064 }
00065 }
00066
00067 edm::Handle< edm::View<pat::Muon> > muon_handle;
00068 evt.getByLabel(muons_, muon_handle);
00069 if(!muon_handle.isValid()) return;
00070 const edm::View<pat::Muon> muons = *muon_handle;
00071 int nmuons = 0;
00072 std::vector<pat::Muon> selMuons;
00073 for(edm::View<pat::Muon>::const_iterator it = muons.begin(); it!=muons.end(); it++) {
00074 reco::TrackRef gltr = it->track();
00075 reco::TrackRef trtr = it->innerTrack();
00076 if(it->pt()>30 && fabs(it->eta())<2.1 && (it->pt()/(it->pt()+it->trackIso()+it->caloIso()))>0.95 && it->isGlobalMuon()){
00077 if(gltr.isNull()) continue;
00078 if((gltr->chi2()/gltr->ndof())<10 && trtr->numberOfValidHits()>11) {
00079 double dRmin = 9999.;
00080 for(std::vector<pat::Jet>::const_iterator ajet = seljets.begin(); ajet != seljets.end(); ajet++) {
00081 math::XYZTLorentzVector jet = ajet->p4();
00082 math::XYZTLorentzVector muon = it->p4();
00083 double tmpdR = DeltaR(muon,jet);
00084 if(tmpdR<dRmin) dRmin = tmpdR;
00085 }
00086 if(dRmin>0.3) {
00087 nmuons++;
00088 selMuons.push_back(*it);
00089 }
00090 }
00091 }
00092 }
00093
00094 edm::Handle< edm::View<pat::Electron> > electron_handle;
00095 evt.getByLabel(electrons_, electron_handle);
00096 if(!electron_handle.isValid()) return;
00097 const edm::View<pat::Electron> electrons = *electron_handle;
00098 int nelectrons = 0;
00099 for(edm::View<pat::Electron>::const_iterator it = electrons.begin(); it!=electrons.end(); it++) {
00100 if(it->pt()>30 && fabs(it->eta())<2.4 && (it->pt()/(it->pt()+it->trackIso()+it->caloIso()))>0.95 && it->isElectronIDAvailable("eidTight"))
00101 {
00102 if(it->electronID("eidTight")==1) nelectrons++;
00103 }
00104 }
00105
00106
00107 double discrim;
00108
00109 if( nmuons!=1 ||
00110 nJets < 4 ||
00111 nelectrons > 0 ) discrim = -1.;
00112 else {
00113
00114 math::XYZTLorentzVector muon = selMuons.begin()->p4();
00115
00116 TtSemiLepSignalSel selection(seljets,muon,MET);
00117
00118 discrim = evaluateTtSemiLepSignalSel(mvaComputer, selection);
00119 }
00120
00121 *pOutDisc = discrim;
00122
00123 evt.put(pOutDisc, "DiscSel");
00124
00125 DiscSel = discrim;
00126
00127 }
00128
00129 void
00130 TtSemiLepSignalSelMVAComputer::beginJob()
00131 {
00132 }
00133
00134 void
00135 TtSemiLepSignalSelMVAComputer::endJob()
00136 {
00137 }
00138
00139 double TtSemiLepSignalSelMVAComputer::DeltaPhi(const math::XYZTLorentzVector& v1, const math::XYZTLorentzVector& v2)
00140 {
00141 double dPhi = fabs(v1.Phi() - v2.Phi());
00142 if (dPhi > TMath::Pi()) dPhi = 2*TMath::Pi() - dPhi;
00143 return dPhi;
00144 }
00145
00146 double TtSemiLepSignalSelMVAComputer::DeltaR(const math::XYZTLorentzVector& v1, const math::XYZTLorentzVector& v2)
00147 {
00148 double dPhi = DeltaPhi(v1,v2);
00149 double dR = TMath::Sqrt((v1.Eta()-v2.Eta())*(v1.Eta()-v2.Eta())+dPhi*dPhi);
00150 return dR;
00151 }
00152
00153
00154
00155
00156 MVA_COMPUTER_CONTAINER_IMPLEMENT(TtSemiLepSignalSelMVA);