CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_1/src/TopQuarkAnalysis/TopEventSelection/plugins/TtSemiLepSignalSelMVATrainer.cc

Go to the documentation of this file.
00001 #include "TMath.h"
00002 #include <algorithm>
00003 
00004 #include "PhysicsTools/JetMCUtils/interface/combination.h"
00005 #include "PhysicsTools/MVATrainer/interface/HelperMacros.h"
00006 
00007 #include "AnalysisDataFormats/TopObjects/interface/TtEvent.h"
00008 #include "TopQuarkAnalysis/TopEventSelection/plugins/TtSemiLepSignalSelMVATrainer.h"
00009 #include "TopQuarkAnalysis/TopEventSelection/interface/TtSemiLepSignalSelEval.h"
00010 
00011 #include "DataFormats/PatCandidates/interface/MET.h"
00012 #include "DataFormats/PatCandidates/interface/Jet.h"
00013 #include "DataFormats/PatCandidates/interface/Muon.h"
00014 #include "DataFormats/PatCandidates/interface/Flags.h"
00015 #include "DataFormats/PatCandidates/interface/Electron.h"
00016 #include "DataFormats/RecoCandidate/interface/RecoCandidate.h"
00017 
00018 
00019 TtSemiLepSignalSelMVATrainer::TtSemiLepSignalSelMVATrainer(const edm::ParameterSet& cfg):
00020   muons_     (cfg.getParameter<edm::InputTag>("muons")),
00021   electrons_ (cfg.getParameter<edm::InputTag>("elecs")),
00022   jets_      (cfg.getParameter<edm::InputTag>("jets")),
00023   METs_      (cfg.getParameter<edm::InputTag>("mets")),
00024   lepChannel_(cfg.getParameter<int>("lepChannel")),
00025   whatData_  (cfg.getParameter<int>("whatData")),
00026   maxEv_     (cfg.getParameter<int>("maxEv"))
00027 {
00028 }
00029 
00030 TtSemiLepSignalSelMVATrainer::~TtSemiLepSignalSelMVATrainer()
00031 {
00032 }
00033 
00034 void
00035 TtSemiLepSignalSelMVATrainer::analyze(const edm::Event& evt, const edm::EventSetup& setup)
00036 {
00037   //communication with CMSSW CondDB
00038   mvaComputer.update<TtSemiLepSignalSelMVARcd>("trainer", setup, "traintreeSaver");
00039 
00040   // can occur in the last iteration when the 
00041   // MVATrainer is about to save the result
00042   if(!mvaComputer) return;
00043 
00044 
00045   //make your preselection here!!
00046   //the following code is for the default example  
00047   edm::Handle<edm::View<pat::MET> > MET_handle;
00048   evt.getByLabel(METs_,MET_handle);
00049   if(!MET_handle.isValid()) return;
00050   const edm::View<pat::MET> MET = *MET_handle;
00051 
00052   edm::Handle< std::vector<pat::Jet> > jet_handle;
00053   evt.getByLabel(jets_, jet_handle);
00054   if(!jet_handle.isValid()) return;
00055   const std::vector<pat::Jet> jets = *jet_handle;
00056   unsigned int nJets = 0;
00057   std::vector<pat::Jet> seljets;
00058   //std::cout<<"number of jets: "<<jets.size()<<std::endl;
00059   for(std::vector<pat::Jet>::const_iterator it = jets.begin(); it != jets.end(); it++) {
00060     //std::cout<<"Jet Et: "<<it->et()<<" Eta: "<<fabs(it->eta())<<std::endl;
00061     if(!(pat::Flags::test(*it, pat::Flags::Overlap::Electrons))) continue;
00062     if(it->et()>30. && fabs(it->eta())<2.4) {
00063       seljets.push_back(*it);
00064       nJets++;
00065     }
00066   }
00067   //std::cout<<"selected Jets: "<<nJets<<std::endl;
00068   if(nJets<4) return;
00069 
00070   //sort by Pt
00071   sort(seljets.begin(),seljets.end(),JetwithHigherPt());
00072  
00073   edm::Handle< edm::View<pat::Muon> > muon_handle; 
00074   evt.getByLabel(muons_, muon_handle);
00075   if(!muon_handle.isValid()) return;
00076   const edm::View<pat::Muon> muons = *muon_handle;
00077   int nmuons = 0;
00078   std::vector<pat::Muon> selMuons;
00079   for(edm::View<pat::Muon>::const_iterator it = muons.begin(); it!=muons.end(); it++) {
00080     reco::TrackRef gltr = it->track(); // global track
00081     reco::TrackRef trtr = it->innerTrack(); // tracker track
00082     if(it->pt()>30 && fabs(it->eta())<2.1 && (it->pt()/(it->pt()+it->trackIso()+it->caloIso()))>0.95 && it->isGlobalMuon()){      
00083       if(gltr.isNull()) continue;  //temporary problems with dead trackrefs
00084       if((gltr->chi2()/gltr->ndof())<10 && trtr->numberOfValidHits()>=11) {
00085      
00086         double dRmin = 9999.;
00087         for(std::vector<pat::Jet>::const_iterator ajet = seljets.begin(); ajet != seljets.end(); ajet++) {
00088           math::XYZTLorentzVector jet = ajet->p4();
00089           math::XYZTLorentzVector muon = it->p4();
00090           double tmpdR = DeltaR(muon,jet);
00091           if(tmpdR<dRmin) dRmin = tmpdR;
00092         }
00093         reco::TrackRef trtr = it->track(); // tracker track
00094         if(dRmin>0.3) {   //temporary problems with muon isolation
00095           nmuons++;
00096           selMuons.push_back(*it);
00097         }
00098       }
00099     }
00100   }
00101   //std::cout<<"selected Muons: "<<nleptons<<std::endl;
00102   if(nmuons!=1) return;
00103   
00104   edm::Handle< edm::View<pat::Electron> > electron_handle; 
00105   evt.getByLabel(electrons_, electron_handle);
00106   if(!electron_handle.isValid()) return;
00107   const edm::View<pat::Electron> electrons = *electron_handle;
00108   int nelectrons = 0;
00109   for(edm::View<pat::Electron>::const_iterator it = electrons.begin(); it!=electrons.end(); it++) {
00110     if(it->pt()>30 && fabs(it->eta())<2.4 && (it->pt()/(it->pt()+it->trackIso()+it->caloIso()))>0.95 && it->isElectronIDAvailable("eidTight"))
00111     { 
00112       if(it->electronID("eidTight")==1) nelectrons++;
00113     }
00114   }
00115   if(nelectrons>0) return;
00116   //end of the preselection
00117 
00118   
00119   math::XYZTLorentzVector muon = selMuons.begin()->p4();
00120 
00121   //count the number of selected events
00122   selEv++;
00123   //skip event if enough events are already selected
00124   if(selEv>maxEv_ && maxEv_!=-1) return;
00125 
00126   //calculation of InputVariables
00127   //see TopQuarkAnalysis/TopTools/interface/TtSemiLepSignalSel.h
00128   //                             /src/TtSemiLepSignalSel.cc
00129   //all objects i.e. jets, muons, electrons... which are needed for the calculation
00130   //of the input-variables have to be passed to this class
00131 
00132   TtSemiLepSignalSel selection(seljets,muon,MET);
00133 
00134   //this is only needed for the default example
00135   edm::Handle<TtGenEvent> genEvt;
00136   evt.getByLabel("genEvt", genEvt);
00137 
00138   double weight = 1.0; //standard no weight, i.e. weight=1.0, set this to the corresponding weight if 
00139                        //different weights for different events are available  
00140   if(whatData_==-1) {  //your training-file contains both, signal and background events
00141     bool isSignal;
00142     isSignal = true;//true for signal, false for background this has to be derived in some way
00143     evaluateTtSemiLepSignalSel(mvaComputer, selection, weight, true, isSignal);
00144   }
00145   else {
00146    
00147     if(whatData_==1){ //your tree contains only signal events
00148       //if needed do a special signal selection here
00149       //the following code is for the default example
00150       if(genEvt->isSemiLeptonic() && genEvt->semiLeptonicChannel() == lepChannel_) {
00151         //std::cout<<"a tt_semlep muon event"<<std::endl;
00152         evaluateTtSemiLepSignalSel(mvaComputer, selection, weight, true, true);
00153       }
00154       else selEv--;
00155     }
00156     else if(whatData_==0){
00157       //std::cout<<"a Wjets event"<<std::endl;
00158       evaluateTtSemiLepSignalSel(mvaComputer, selection, weight, true, false);
00159     }
00160     else std::cout<<"Config File Error!! Please check <whatData> in TtSemiLepSignalSelMVATrainer.cfi";
00161   }
00162 }
00163 
00164 void TtSemiLepSignalSelMVATrainer::beginJob(){
00165   selEv = 0;
00166   if(whatData_!=-1 && whatData_!=0 && whatData_!=1){
00167     std::cout<<"Config File Error!! Please check <whatData> in TtSemiLepSignalSelMVATrainer.cfi"<<std::endl;;
00168     return;
00169   }
00170 }
00171 
00172 double TtSemiLepSignalSelMVATrainer::DeltaPhi(math::XYZTLorentzVector v1, math::XYZTLorentzVector v2)
00173 {
00174   double dPhi = fabs(v1.Phi() - v2.Phi());
00175   if (dPhi > TMath::Pi()) dPhi =  2*TMath::Pi() - dPhi;
00176   return dPhi;
00177 }
00178 
00179 double TtSemiLepSignalSelMVATrainer::DeltaR(math::XYZTLorentzVector v1, math::XYZTLorentzVector v2)
00180 {
00181   double dPhi = DeltaPhi(v1,v2);
00182   double dR = TMath::Sqrt((v1.Eta()-v2.Eta())*(v1.Eta()-v2.Eta())+dPhi*dPhi);
00183   return dR;
00184 }
00185 
00186 // implement the plugins for the trainer
00187 // -> defines TtSemiLepSignalSelMVAContainerSaveCondDB
00188 // -> defines TtSemiLepSignalSelMVASaveFile
00189 // -> defines TtSemiLepSignalSelMVATrainerLooper
00190 MVA_TRAINER_IMPLEMENT(TtSemiLepSignalSelMVA);