CMS 3D CMS Logo

TtSemiLepSignalSelMVATrainer.cc
Go to the documentation of this file.
1 #include "TMath.h"
2 #include <algorithm>
3 
6 
9 
11 
12 
14  muonsToken_ (consumes< edm::View<pat::Muon> >(cfg.getParameter<edm::InputTag>("muons"))),
15  electronsToken_ (consumes< edm::View<pat::Electron> >(cfg.getParameter<edm::InputTag>("elecs"))),
16  jetsToken_ (consumes< std::vector<pat::Jet> >(cfg.getParameter<edm::InputTag>("jets"))),
17  METsToken_ (consumes<edm::View<pat::MET> >(cfg.getParameter<edm::InputTag>("mets"))),
18  genEvtToken_ (consumes<TtGenEvent>(edm::InputTag("genEvt"))),
19  lepChannel_(cfg.getParameter<int>("lepChannel")),
20  whatData_ (cfg.getParameter<int>("whatData")),
21  maxEv_ (cfg.getParameter<int>("maxEv"))
22 {
23 }
24 
26 {
27 }
28 
29 void
31 {
32  //communication with CMSSW CondDB
33  mvaComputer.update<TtSemiLepSignalSelMVARcd>("trainer", setup, "traintreeSaver");
34 
35  // can occur in the last iteration when the
36  // MVATrainer is about to save the result
37  if(!mvaComputer) return;
38 
39 
40  //make your preselection here!!
41  //the following code is for the default example
43  evt.getByToken(METsToken_,MET_handle);
44  if(!MET_handle.isValid()) return;
45  const edm::View<pat::MET> MET = *MET_handle;
46 
48  evt.getByToken(jetsToken_, jet_handle);
49  if(!jet_handle.isValid()) return;
50  const std::vector<pat::Jet> jets = *jet_handle;
51  unsigned int nJets = 0;
52  std::vector<pat::Jet> seljets;
53  //std::cout<<"number of jets: "<<jets.size()<<std::endl;
54  for(std::vector<pat::Jet>::const_iterator it = jets.begin(); it != jets.end(); it++) {
55  //std::cout<<"Jet Et: "<<it->et()<<" Eta: "<<fabs(it->eta())<<std::endl;
57  if(it->et()>30. && fabs(it->eta())<2.4) {
58  seljets.push_back(*it);
59  nJets++;
60  }
61  }
62  //std::cout<<"selected Jets: "<<nJets<<std::endl;
63  if(nJets<4) return;
64 
65  //sort by Pt
66  sort(seljets.begin(),seljets.end(),JetwithHigherPt());
67 
69  evt.getByToken(muonsToken_, muon_handle);
70  if(!muon_handle.isValid()) return;
71  const edm::View<pat::Muon> muons = *muon_handle;
72  int nmuons = 0;
73  std::vector<pat::Muon> selMuons;
74  for(edm::View<pat::Muon>::const_iterator it = muons.begin(); it!=muons.end(); it++) {
75  reco::TrackRef gltr = it->track(); // global track
76  reco::TrackRef trtr = it->innerTrack(); // tracker track
77  if(it->pt()>30 && fabs(it->eta())<2.1 && (it->pt()/(it->pt()+it->trackIso()+it->caloIso()))>0.95 && it->isGlobalMuon()){
78  if(gltr.isNull()) continue; //temporary problems with dead trackrefs
79  if((gltr->chi2()/gltr->ndof())<10 && trtr->numberOfValidHits()>=11) {
80 
81  double dRmin = 9999.;
82  for(std::vector<pat::Jet>::const_iterator ajet = seljets.begin(); ajet != seljets.end(); ajet++) {
83  math::XYZTLorentzVector jet = ajet->p4();
84  math::XYZTLorentzVector muon = it->p4();
85  double tmpdR = DeltaR(muon,jet);
86  if(tmpdR<dRmin) dRmin = tmpdR;
87  }
88  if(dRmin>0.3) { //temporary problems with muon isolation
89  nmuons++;
90  selMuons.push_back(*it);
91  }
92  }
93  }
94  }
95  //std::cout<<"selected Muons: "<<nleptons<<std::endl;
96  if(nmuons!=1) return;
97 
98  edm::Handle< edm::View<pat::Electron> > electron_handle;
99  evt.getByToken(electronsToken_, electron_handle);
100  if(!electron_handle.isValid()) return;
101  const edm::View<pat::Electron> electrons = *electron_handle;
102  int nelectrons = 0;
103  for(edm::View<pat::Electron>::const_iterator it = electrons.begin(); it!=electrons.end(); it++) {
104  if(it->pt()>30 && fabs(it->eta())<2.4 && (it->pt()/(it->pt()+it->trackIso()+it->caloIso()))>0.95 && it->isElectronIDAvailable("eidTight"))
105  {
106  if(it->electronID("eidTight")==1) nelectrons++;
107  }
108  }
109  if(nelectrons>0) return;
110  //end of the preselection
111 
112 
113  math::XYZTLorentzVector muon = selMuons.begin()->p4();
114 
115  //count the number of selected events
116  selEv++;
117  //skip event if enough events are already selected
118  if(selEv>maxEv_ && maxEv_!=-1) return;
119 
120  //calculation of InputVariables
121  //see TopQuarkAnalysis/TopTools/interface/TtSemiLepSignalSel.h
122  // /src/TtSemiLepSignalSel.cc
123  //all objects i.e. jets, muons, electrons... which are needed for the calculation
124  //of the input-variables have to be passed to this class
125 
126  TtSemiLepSignalSel selection(seljets,muon,MET);
127 
128  //this is only needed for the default example
130  evt.getByToken(genEvtToken_, genEvt);
131 
132  double weight = 1.0; //standard no weight, i.e. weight=1.0, set this to the corresponding weight if
133  //different weights for different events are available
134  if(whatData_==-1) { //your training-file contains both, signal and background events
135  bool isSignal;
136  isSignal = true;//true for signal, false for background this has to be derived in some way
137  evaluateTtSemiLepSignalSel(mvaComputer, selection, weight, true, isSignal);
138  }
139  else {
140 
141  if(whatData_==1){ //your tree contains only signal events
142  //if needed do a special signal selection here
143  //the following code is for the default example
144  if(genEvt->isSemiLeptonic() && genEvt->semiLeptonicChannel() == lepChannel_) {
145  //std::cout<<"a tt_semlep muon event"<<std::endl;
146  evaluateTtSemiLepSignalSel(mvaComputer, selection, weight, true, true);
147  }
148  else selEv--;
149  }
150  else if(whatData_==0){
151  //std::cout<<"a Wjets event"<<std::endl;
152  evaluateTtSemiLepSignalSel(mvaComputer, selection, weight, true, false);
153  }
154  else std::cout<<"Config File Error!! Please check <whatData> in TtSemiLepSignalSelMVATrainer.cfi";
155  }
156 }
157 
159  selEv = 0;
160  if(whatData_!=-1 && whatData_!=0 && whatData_!=1){
161  std::cout<<"Config File Error!! Please check <whatData> in TtSemiLepSignalSelMVATrainer.cfi"<<std::endl;;
162  return;
163  }
164 }
165 
167 {
168  double dPhi = fabs(v1.Phi() - v2.Phi());
169  if (dPhi > TMath::Pi()) dPhi = 2*TMath::Pi() - dPhi;
170  return dPhi;
171 }
172 
174 {
175  double dPhi = DeltaPhi(v1,v2);
176  double dR = TMath::Sqrt((v1.Eta()-v2.Eta())*(v1.Eta()-v2.Eta())+dPhi*dPhi);
177  return dR;
178 }
179 
180 // implement the plugins for the trainer
181 // -> defines TtSemiLepSignalSelMVAContainerSaveCondDB
182 // -> defines TtSemiLepSignalSelMVASaveFile
183 // -> defines TtSemiLepSignalSelMVATrainerLooper
184 MVA_TRAINER_IMPLEMENT(TtSemiLepSignalSelMVA);
const double Pi
edm::EDGetTokenT< edm::View< pat::MET > > METsToken_
bool isSemiLeptonic(bool excludeTauLeptons=false) const
check if the event can be classified as semi-laptonic
Definition: TtGenEvent.h:38
void analyze(const edm::Event &evt, const edm::EventSetup &setup) override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
double evaluateTtSemiLepSignalSel(PhysicsTools::MVAComputerCache &mvaComputer, const TtSemiLepSignalSel &sigsel, float weight=1., const bool training=false, const bool isSignal=false)
PhysicsTools::MVAComputerCache mvaComputer
static bool test(uint32_t val, uint32_t mask)
Definition: Flags.h:28
#define MVA_TRAINER_IMPLEMENT(N)
Definition: HelperMacros.h:40
selection
main part
Definition: corrVsCorr.py:98
def setup(process, global_tag, zero_tesla=False)
Definition: GeneralSetup.py:1
Definition: weight.py:1
TtSemiLepSignalSelMVATrainer(const edm::ParameterSet &)
Definition: HeavyIon.h:7
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
double DeltaPhi(const math::XYZTLorentzVector &v1, const math::XYZTLorentzVector &v2)
edm::EDGetTokenT< std::vector< pat::Jet > > jetsToken_
const_iterator begin() const
edm::EDGetTokenT< TtGenEvent > genEvtToken_
Class derived from the TopGenEvent for ttbar events.
Definition: TtGenEvent.h:18
Definition: Muon.py:1
WDecay::LeptonType semiLeptonicChannel() const
return decay channel; all leptons including taus are allowed
Definition: TtGenEvent.cc:39
Definition: Jet.py:1
vector< PseudoJet > jets
edm::EDGetTokenT< edm::View< pat::Muon > > muonsToken_
bool isValid() const
Definition: HandleBase.h:74
bool isNull() const
Checks for null.
Definition: Ref.h:250
bool update(const Calibration::MVAComputer *computer)
HLT enums.
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
const_iterator end() const
double DeltaR(const math::XYZTLorentzVector &v1, const math::XYZTLorentzVector &v2)
edm::EDGetTokenT< edm::View< pat::Electron > > electronsToken_