#include <TtSemiLepSignalSelMVATrainer.h>
Classes | |
struct | JetwithHigherPt |
Public Member Functions | |
TtSemiLepSignalSelMVATrainer (const edm::ParameterSet &) | |
~TtSemiLepSignalSelMVATrainer () | |
Private Member Functions | |
virtual void | analyze (const edm::Event &evt, const edm::EventSetup &setup) |
virtual void | beginJob () |
double | DeltaPhi (math::XYZTLorentzVector v1, math::XYZTLorentzVector v2) |
double | DeltaR (math::XYZTLorentzVector v1, math::XYZTLorentzVector v2) |
Private Attributes | |
edm::InputTag | electrons_ |
edm::InputTag | jets_ |
int | lepChannel_ |
int | maxEv_ |
edm::InputTag | METs_ |
edm::InputTag | muons_ |
PhysicsTools::MVAComputerCache | mvaComputer |
int | selEv |
int | whatData_ |
Definition at line 24 of file TtSemiLepSignalSelMVATrainer.h.
TtSemiLepSignalSelMVATrainer::TtSemiLepSignalSelMVATrainer | ( | const edm::ParameterSet & | cfg | ) | [explicit] |
Definition at line 19 of file TtSemiLepSignalSelMVATrainer.cc.
: muons_ (cfg.getParameter<edm::InputTag>("muons")), electrons_ (cfg.getParameter<edm::InputTag>("elecs")), jets_ (cfg.getParameter<edm::InputTag>("jets")), METs_ (cfg.getParameter<edm::InputTag>("mets")), lepChannel_(cfg.getParameter<int>("lepChannel")), whatData_ (cfg.getParameter<int>("whatData")), maxEv_ (cfg.getParameter<int>("maxEv")) { }
TtSemiLepSignalSelMVATrainer::~TtSemiLepSignalSelMVATrainer | ( | ) |
Definition at line 30 of file TtSemiLepSignalSelMVATrainer.cc.
{ }
void TtSemiLepSignalSelMVATrainer::analyze | ( | const edm::Event & | evt, |
const edm::EventSetup & | setup | ||
) | [private, virtual] |
Implements edm::EDAnalyzer.
Definition at line 35 of file TtSemiLepSignalSelMVATrainer.cc.
References edm::View< T >::begin(), gather_cfg::cout, DeltaR(), HI_PhotonSkim_cff::electrons, pat::Flags::Overlap::Electrons, electrons_, edm::View< T >::end(), evaluateTtSemiLepSignalSel(), TtGenEvtProducer_cfi::genEvt, edm::Event::getByLabel(), edm::Ref< C, T, F >::isNull(), edm::HandleBase::isValid(), metsig::jet, analyzePatCleaning_cfg::jets, jets_, lepChannel_, maxEv_, METs_, metsig::muon, patZpeak::muons, muons_, mvaComputer, corrVsCorr::selection, selEv, HcalObjRepresent::setup(), python::multivaluedict::sort(), pat::Flags::test(), PhysicsTools::MVAComputerCache::update(), CommonMethods::weight(), and whatData_.
{ //communication with CMSSW CondDB mvaComputer.update<TtSemiLepSignalSelMVARcd>("trainer", setup, "traintreeSaver"); // can occur in the last iteration when the // MVATrainer is about to save the result if(!mvaComputer) return; //make your preselection here!! //the following code is for the default example 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; //std::cout<<"number of jets: "<<jets.size()<<std::endl; for(std::vector<pat::Jet>::const_iterator it = jets.begin(); it != jets.end(); it++) { //std::cout<<"Jet Et: "<<it->et()<<" Eta: "<<fabs(it->eta())<<std::endl; if(!(pat::Flags::test(*it, pat::Flags::Overlap::Electrons))) continue; if(it->et()>30. && fabs(it->eta())<2.4) { seljets.push_back(*it); nJets++; } } //std::cout<<"selected Jets: "<<nJets<<std::endl; if(nJets<4) return; //sort by Pt sort(seljets.begin(),seljets.end(),JetwithHigherPt()); 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; } reco::TrackRef trtr = it->track(); // tracker track if(dRmin>0.3) { //temporary problems with muon isolation nmuons++; selMuons.push_back(*it); } } } } //std::cout<<"selected Muons: "<<nleptons<<std::endl; if(nmuons!=1) return; 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++; } } if(nelectrons>0) return; //end of the preselection math::XYZTLorentzVector muon = selMuons.begin()->p4(); //count the number of selected events selEv++; //skip event if enough events are already selected if(selEv>maxEv_ && maxEv_!=-1) return; //calculation of InputVariables //see TopQuarkAnalysis/TopTools/interface/TtSemiLepSignalSel.h // /src/TtSemiLepSignalSel.cc //all objects i.e. jets, muons, electrons... which are needed for the calculation //of the input-variables have to be passed to this class TtSemiLepSignalSel selection(seljets,muon,MET); //this is only needed for the default example edm::Handle<TtGenEvent> genEvt; evt.getByLabel("genEvt", genEvt); double weight = 1.0; //standard no weight, i.e. weight=1.0, set this to the corresponding weight if //different weights for different events are available if(whatData_==-1) { //your training-file contains both, signal and background events bool isSignal; isSignal = true;//true for signal, false for background this has to be derived in some way evaluateTtSemiLepSignalSel(mvaComputer, selection, weight, true, isSignal); } else { if(whatData_==1){ //your tree contains only signal events //if needed do a special signal selection here //the following code is for the default example if(genEvt->isSemiLeptonic() && genEvt->semiLeptonicChannel() == lepChannel_) { //std::cout<<"a tt_semlep muon event"<<std::endl; evaluateTtSemiLepSignalSel(mvaComputer, selection, weight, true, true); } else selEv--; } else if(whatData_==0){ //std::cout<<"a Wjets event"<<std::endl; evaluateTtSemiLepSignalSel(mvaComputer, selection, weight, true, false); } else std::cout<<"Config File Error!! Please check <whatData> in TtSemiLepSignalSelMVATrainer.cfi"; } }
void TtSemiLepSignalSelMVATrainer::beginJob | ( | void | ) | [private, virtual] |
Reimplemented from edm::EDAnalyzer.
Definition at line 164 of file TtSemiLepSignalSelMVATrainer.cc.
References gather_cfg::cout, selEv, and whatData_.
double TtSemiLepSignalSelMVATrainer::DeltaPhi | ( | math::XYZTLorentzVector | v1, |
math::XYZTLorentzVector | v2 | ||
) | [private] |
double TtSemiLepSignalSelMVATrainer::DeltaR | ( | math::XYZTLorentzVector | v1, |
math::XYZTLorentzVector | v2 | ||
) | [private] |
Definition at line 179 of file TtSemiLepSignalSelMVATrainer.cc.
References DeltaPhi(), dPhi(), and PFRecoTauDiscriminationAgainstElectronDeadECAL_cfi::dR.
Referenced by analyze().
Definition at line 47 of file TtSemiLepSignalSelMVATrainer.h.
Referenced by analyze().
Definition at line 48 of file TtSemiLepSignalSelMVATrainer.h.
Referenced by analyze().
int TtSemiLepSignalSelMVATrainer::lepChannel_ [private] |
Definition at line 51 of file TtSemiLepSignalSelMVATrainer.h.
Referenced by analyze().
int TtSemiLepSignalSelMVATrainer::maxEv_ [private] |
Definition at line 53 of file TtSemiLepSignalSelMVATrainer.h.
Referenced by analyze().
Definition at line 49 of file TtSemiLepSignalSelMVATrainer.h.
Referenced by analyze().
Definition at line 46 of file TtSemiLepSignalSelMVATrainer.h.
Referenced by analyze().
Definition at line 56 of file TtSemiLepSignalSelMVATrainer.h.
Referenced by analyze().
int TtSemiLepSignalSelMVATrainer::selEv [private] |
Definition at line 54 of file TtSemiLepSignalSelMVATrainer.h.
Referenced by analyze(), and beginJob().
int TtSemiLepSignalSelMVATrainer::whatData_ [private] |
Definition at line 52 of file TtSemiLepSignalSelMVATrainer.h.
Referenced by analyze(), and beginJob().