CMS 3D CMS Logo

List of all members | Classes | Public Member Functions | Private Member Functions | Private Attributes
TtSemiLepSignalSelMVATrainer Class Reference

#include <TtSemiLepSignalSelMVATrainer.h>

Inheritance diagram for TtSemiLepSignalSelMVATrainer:
edm::EDAnalyzer edm::EDConsumerBase

Classes

struct  JetwithHigherPt
 

Public Member Functions

 TtSemiLepSignalSelMVATrainer (const edm::ParameterSet &)
 
 ~TtSemiLepSignalSelMVATrainer () override
 
- Public Member Functions inherited from edm::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzer ()
 
SerialTaskQueueglobalLuminosityBlocksQueue ()
 
SerialTaskQueueglobalRunsQueue ()
 
ModuleDescription const & moduleDescription () const
 
std::string workerType () const
 
 ~EDAnalyzer () override
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
virtual ~EDConsumerBase () noexcept(false)
 

Private Member Functions

void analyze (const edm::Event &evt, const edm::EventSetup &setup) override
 
void beginJob () override
 
double DeltaPhi (const math::XYZTLorentzVector &v1, const math::XYZTLorentzVector &v2)
 
double DeltaR (const math::XYZTLorentzVector &v1, const math::XYZTLorentzVector &v2)
 

Private Attributes

edm::EDGetTokenT< edm::View< pat::Electron > > electronsToken_
 
edm::EDGetTokenT< TtGenEventgenEvtToken_
 
edm::EDGetTokenT< std::vector< pat::Jet > > jetsToken_
 
int lepChannel_
 
int maxEv_
 
edm::EDGetTokenT< edm::View< pat::MET > > METsToken_
 
edm::EDGetTokenT< edm::View< pat::Muon > > muonsToken_
 
PhysicsTools::MVAComputerCache mvaComputer
 
int selEv
 
int whatData_
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
static bool wantsGlobalLuminosityBlocks ()
 
static bool wantsGlobalRuns ()
 
static bool wantsStreamLuminosityBlocks ()
 
static bool wantsStreamRuns ()
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Definition at line 28 of file TtSemiLepSignalSelMVATrainer.h.

Constructor & Destructor Documentation

TtSemiLepSignalSelMVATrainer::TtSemiLepSignalSelMVATrainer ( const edm::ParameterSet cfg)
explicit

Definition at line 13 of file TtSemiLepSignalSelMVATrainer.cc.

13  :
16  jetsToken_ (consumes< std::vector<pat::Jet> >(cfg.getParameter<edm::InputTag>("jets"))),
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 }
T getParameter(std::string const &) const
edm::EDGetTokenT< edm::View< pat::MET > > METsToken_
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
edm::EDGetTokenT< std::vector< pat::Jet > > jetsToken_
edm::EDGetTokenT< TtGenEvent > genEvtToken_
edm::EDGetTokenT< edm::View< pat::Muon > > muonsToken_
edm::EDGetTokenT< edm::View< pat::Electron > > electronsToken_
TtSemiLepSignalSelMVATrainer::~TtSemiLepSignalSelMVATrainer ( )
override

Definition at line 25 of file TtSemiLepSignalSelMVATrainer.cc.

26 {
27 }

Member Function Documentation

void TtSemiLepSignalSelMVATrainer::analyze ( const edm::Event evt,
const edm::EventSetup setup 
)
overrideprivate

Definition at line 30 of file TtSemiLepSignalSelMVATrainer.cc.

References edm::View< T >::begin(), gather_cfg::cout, DeltaR(), nano_cff::electrons, pat::Flags::Overlap::Electrons, electronsToken_, edm::View< T >::end(), evaluateTtSemiLepSignalSel(), TtGenEvtProducer_cfi::genEvt, genEvtToken_, edm::Event::getByToken(), edm::Ref< C, T, F >::isNull(), TtGenEvent::isSemiLeptonic(), edm::HandleBase::isValid(), metsig::jet, fwrapper::jets, jetsToken_, lepChannel_, maxEv_, METsToken_, extraflags_cff::muons, muonsToken_, mvaComputer, corrVsCorr::selection, selEv, TtGenEvent::semiLeptonicChannel(), GeneralSetup::setup(), pat::Flags::test(), PhysicsTools::MVAComputerCache::update(), and whatData_.

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 }
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
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
selection
main part
Definition: corrVsCorr.py:98
def setup(process, global_tag, zero_tesla=False)
Definition: GeneralSetup.py:1
Definition: weight.py:1
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
edm::EDGetTokenT< std::vector< pat::Jet > > jetsToken_
const_iterator begin() const
edm::EDGetTokenT< TtGenEvent > genEvtToken_
WDecay::LeptonType semiLeptonicChannel() const
return decay channel; all leptons including taus are allowed
Definition: TtGenEvent.cc:39
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)
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_
void TtSemiLepSignalSelMVATrainer::beginJob ( void  )
overrideprivatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 158 of file TtSemiLepSignalSelMVATrainer.cc.

References gather_cfg::cout, selEv, and whatData_.

158  {
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 }
double TtSemiLepSignalSelMVATrainer::DeltaPhi ( const math::XYZTLorentzVector v1,
const math::XYZTLorentzVector v2 
)
private

Definition at line 166 of file TtSemiLepSignalSelMVATrainer.cc.

References particleFlow_cfi::dPhi, and Pi.

Referenced by DeltaR().

167 {
168  double dPhi = fabs(v1.Phi() - v2.Phi());
169  if (dPhi > TMath::Pi()) dPhi = 2*TMath::Pi() - dPhi;
170  return dPhi;
171 }
const double Pi
double TtSemiLepSignalSelMVATrainer::DeltaR ( const math::XYZTLorentzVector v1,
const math::XYZTLorentzVector v2 
)
private

Definition at line 173 of file TtSemiLepSignalSelMVATrainer.cc.

References DeltaPhi(), particleFlow_cfi::dPhi, and PFRecoTauDiscriminationAgainstElectronDeadECAL_cfi::dR.

Referenced by analyze().

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 }
double DeltaPhi(const math::XYZTLorentzVector &v1, const math::XYZTLorentzVector &v2)

Member Data Documentation

edm::EDGetTokenT< edm::View<pat::Electron> > TtSemiLepSignalSelMVATrainer::electronsToken_
private

Definition at line 51 of file TtSemiLepSignalSelMVATrainer.h.

Referenced by analyze().

edm::EDGetTokenT<TtGenEvent> TtSemiLepSignalSelMVATrainer::genEvtToken_
private

Definition at line 54 of file TtSemiLepSignalSelMVATrainer.h.

Referenced by analyze().

edm::EDGetTokenT< std::vector<pat::Jet> > TtSemiLepSignalSelMVATrainer::jetsToken_
private

Definition at line 52 of file TtSemiLepSignalSelMVATrainer.h.

Referenced by analyze().

int TtSemiLepSignalSelMVATrainer::lepChannel_
private

Definition at line 56 of file TtSemiLepSignalSelMVATrainer.h.

Referenced by analyze().

int TtSemiLepSignalSelMVATrainer::maxEv_
private

Definition at line 58 of file TtSemiLepSignalSelMVATrainer.h.

Referenced by analyze().

edm::EDGetTokenT<edm::View<pat::MET> > TtSemiLepSignalSelMVATrainer::METsToken_
private

Definition at line 53 of file TtSemiLepSignalSelMVATrainer.h.

Referenced by analyze().

edm::EDGetTokenT< edm::View<pat::Muon> > TtSemiLepSignalSelMVATrainer::muonsToken_
private

Definition at line 50 of file TtSemiLepSignalSelMVATrainer.h.

Referenced by analyze().

PhysicsTools::MVAComputerCache TtSemiLepSignalSelMVATrainer::mvaComputer
private

Definition at line 61 of file TtSemiLepSignalSelMVATrainer.h.

Referenced by analyze().

int TtSemiLepSignalSelMVATrainer::selEv
private

Definition at line 59 of file TtSemiLepSignalSelMVATrainer.h.

Referenced by analyze(), and beginJob().

int TtSemiLepSignalSelMVATrainer::whatData_
private

Definition at line 57 of file TtSemiLepSignalSelMVATrainer.h.

Referenced by analyze(), and beginJob().