#include <TtSemiLepSignalSelMVAComputer.h>
Public Member Functions | |
TtSemiLepSignalSelMVAComputer (const edm::ParameterSet &) | |
~TtSemiLepSignalSelMVAComputer () | |
Private Member Functions | |
virtual void | beginJob () |
double | DeltaPhi (math::XYZTLorentzVector v1, math::XYZTLorentzVector v2) |
double | DeltaR (math::XYZTLorentzVector v1, math::XYZTLorentzVector v2) |
virtual void | endJob () |
virtual void | produce (edm::Event &evt, const edm::EventSetup &setup) |
Private Attributes | |
double | DiscSel |
edm::InputTag | electrons_ |
edm::InputTag | jets_ |
edm::InputTag | METs_ |
edm::InputTag | muons_ |
PhysicsTools::MVAComputerCache | mvaComputer |
Definition at line 17 of file TtSemiLepSignalSelMVAComputer.h.
TtSemiLepSignalSelMVAComputer::TtSemiLepSignalSelMVAComputer | ( | const edm::ParameterSet & | cfg | ) | [explicit] |
Definition at line 18 of file TtSemiLepSignalSelMVAComputer.cc.
: muons_ (cfg.getParameter<edm::InputTag>("muons")), jets_ (cfg.getParameter<edm::InputTag>("jets")), METs_ (cfg.getParameter<edm::InputTag>("mets")), electrons_ (cfg.getParameter<edm::InputTag>("elecs")) { produces< double >("DiscSel"); }
TtSemiLepSignalSelMVAComputer::~TtSemiLepSignalSelMVAComputer | ( | ) |
Definition at line 29 of file TtSemiLepSignalSelMVAComputer.cc.
{ }
void TtSemiLepSignalSelMVAComputer::beginJob | ( | void | ) | [private, virtual] |
Reimplemented from edm::EDProducer.
Definition at line 130 of file TtSemiLepSignalSelMVAComputer.cc.
{ }
double TtSemiLepSignalSelMVAComputer::DeltaPhi | ( | math::XYZTLorentzVector | v1, |
math::XYZTLorentzVector | v2 | ||
) | [private] |
double TtSemiLepSignalSelMVAComputer::DeltaR | ( | math::XYZTLorentzVector | v1, |
math::XYZTLorentzVector | v2 | ||
) | [private] |
Definition at line 146 of file TtSemiLepSignalSelMVAComputer.cc.
References DeltaPhi(), and dPhi().
Referenced by produce().
void TtSemiLepSignalSelMVAComputer::endJob | ( | void | ) | [private, virtual] |
Reimplemented from edm::EDProducer.
Definition at line 135 of file TtSemiLepSignalSelMVAComputer.cc.
{ }
void TtSemiLepSignalSelMVAComputer::produce | ( | edm::Event & | evt, |
const edm::EventSetup & | setup | ||
) | [private, virtual] |
Implements edm::EDProducer.
Definition at line 34 of file TtSemiLepSignalSelMVAComputer.cc.
References edm::View< T >::begin(), DeltaR(), DiscSel, pat::Flags::Overlap::Electrons, electrons_, edm::View< T >::end(), evaluateTtSemiLepSignalSel(), edm::Event::getByLabel(), edm::Ref< C, T, F >::isNull(), edm::HandleBase::isValid(), metsig::jet, analyzePatCleaning_cfg::jets, jets_, METs_, metsig::muon, ExpressReco_HICollisions_FallBack::muons, muons_, mvaComputer, edm::Event::put(), elec_selection::selection(), pat::Flags::test(), and PhysicsTools::MVAComputerCache::update().
{ std::auto_ptr< double > pOutDisc (new double); mvaComputer.update<TtSemiLepSignalSelMVARcd>(setup, "ttSemiLepSignalSelMVA"); // read name of the last processor in the MVA calibration // (to be used as meta information) edm::ESHandle<PhysicsTools::Calibration::MVAComputerContainer> calibContainer; setup.get<TtSemiLepSignalSelMVARcd>().get( calibContainer ); std::vector<PhysicsTools::Calibration::VarProcessor*> processors = (calibContainer->find("ttSemiLepSignalSelMVA")).getProcessors(); //make your preselection! This must!! be the same one as in TraintreeSaver.cc edm::Handle<edm::View<pat::MET> > MET_handle; evt.getByLabel(METs_,MET_handle); if(!MET_handle.isValid()) return; const edm::View<pat::MET> MET = *MET_handle; edm::Handle< std::vector<pat::Jet> > jet_handle; evt.getByLabel(jets_, jet_handle); if(!jet_handle.isValid()) return; const std::vector<pat::Jet> jets = *jet_handle; unsigned int nJets = 0; std::vector<pat::Jet> seljets; for(std::vector<pat::Jet>::const_iterator it = jets.begin(); it != jets.end(); it++) { if(!(pat::Flags::test(*it, pat::Flags::Overlap::Electrons))) continue; if(it->pt()>30. && fabs(it->eta())<2.4) { seljets.push_back(*it); nJets++; } } edm::Handle< edm::View<pat::Muon> > muon_handle; evt.getByLabel(muons_, muon_handle); if(!muon_handle.isValid()) return; const edm::View<pat::Muon> muons = *muon_handle; int nmuons = 0; std::vector<pat::Muon> selMuons; for(edm::View<pat::Muon>::const_iterator it = muons.begin(); it!=muons.end(); it++) { reco::TrackRef gltr = it->track(); // global track reco::TrackRef trtr = it->innerTrack(); // tracker track if(it->pt()>30 && fabs(it->eta())<2.1 && (it->pt()/(it->pt()+it->trackIso()+it->caloIso()))>0.95 && it->isGlobalMuon()){ if(gltr.isNull()) continue; //temporary problems with dead trackrefs if((gltr->chi2()/gltr->ndof())<10 && trtr->numberOfValidHits()>11) { double dRmin = 9999.; for(std::vector<pat::Jet>::const_iterator ajet = seljets.begin(); ajet != seljets.end(); ajet++) { math::XYZTLorentzVector jet = ajet->p4(); math::XYZTLorentzVector muon = it->p4(); double tmpdR = DeltaR(muon,jet); if(tmpdR<dRmin) dRmin = tmpdR; } if(dRmin>0.3) { //temporary problems with muon isolation nmuons++; selMuons.push_back(*it); } } } } edm::Handle< edm::View<pat::Electron> > electron_handle; evt.getByLabel(electrons_, electron_handle); if(!electron_handle.isValid()) return; const edm::View<pat::Electron> electrons = *electron_handle; int nelectrons = 0; for(edm::View<pat::Electron>::const_iterator it = electrons.begin(); it!=electrons.end(); it++) { if(it->pt()>30 && fabs(it->eta())<2.4 && (it->pt()/(it->pt()+it->trackIso()+it->caloIso()))>0.95 && it->isElectronIDAvailable("eidTight")) { if(it->electronID("eidTight")==1) nelectrons++; } } double discrim; // discriminator output for events which do not pass the preselection is set to -1 if( nmuons!=1 || nJets < 4 || nelectrons > 0 ) discrim = -1.; //std::cout<<"nJets: "<<seljets.size()<<" numLeptons: "<<nleptons<<std::endl;} else { //check wheter a event was already selected (problem with duplicated events) math::XYZTLorentzVector muon = selMuons.begin()->p4(); TtSemiLepSignalSel selection(seljets,muon,MET); discrim = evaluateTtSemiLepSignalSel(mvaComputer, selection); } *pOutDisc = discrim; evt.put(pOutDisc, "DiscSel"); DiscSel = discrim; }
double TtSemiLepSignalSelMVAComputer::DiscSel [private] |
Definition at line 40 of file TtSemiLepSignalSelMVAComputer.h.
Referenced by produce().
Definition at line 36 of file TtSemiLepSignalSelMVAComputer.h.
Referenced by produce().
Definition at line 34 of file TtSemiLepSignalSelMVAComputer.h.
Referenced by produce().
Definition at line 35 of file TtSemiLepSignalSelMVAComputer.h.
Referenced by produce().
Definition at line 33 of file TtSemiLepSignalSelMVAComputer.h.
Referenced by produce().
Definition at line 38 of file TtSemiLepSignalSelMVAComputer.h.
Referenced by produce().