CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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

Classes

struct  JetwithHigherPt
 

Public Member Functions

 TtSemiLepSignalSelMVATrainer (const edm::ParameterSet &)
 
 ~TtSemiLepSignalSelMVATrainer ()
 
- Public Member Functions inherited from edm::EDAnalyzer
 EDAnalyzer ()
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 

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_
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
typedef WorkerT< EDAnalyzerWorkerType
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::EDAnalyzer
CurrentProcessingContext const * currentContext () const
 

Detailed Description

Definition at line 24 of file TtSemiLepSignalSelMVATrainer.h.

Constructor & Destructor Documentation

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

Definition at line 30 of file TtSemiLepSignalSelMVATrainer.cc.

31 {
32 }

Member Function Documentation

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

Implements edm::EDAnalyzer.

Definition at line 35 of file TtSemiLepSignalSelMVATrainer.cc.

References edm::View< T >::begin(), gather_cfg::cout, DeltaR(), 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, ExpressReco_HICollisions_FallBack::muons, muons_, mvaComputer, elec_selection::selection(), selEv, python.multivaluedict::sort(), pat::Flags::test(), PhysicsTools::MVAComputerCache::update(), CommonMethods::weight(), and whatData_.

36 {
37  //communication with CMSSW CondDB
38  mvaComputer.update<TtSemiLepSignalSelMVARcd>("trainer", setup, "traintreeSaver");
39 
40  // can occur in the last iteration when the
41  // MVATrainer is about to save the result
42  if(!mvaComputer) return;
43 
44 
45  //make your preselection here!!
46  //the following code is for the default example
48  evt.getByLabel(METs_,MET_handle);
49  if(!MET_handle.isValid()) return;
50  const edm::View<pat::MET> MET = *MET_handle;
51 
53  evt.getByLabel(jets_, jet_handle);
54  if(!jet_handle.isValid()) return;
55  const std::vector<pat::Jet> jets = *jet_handle;
56  unsigned int nJets = 0;
57  std::vector<pat::Jet> seljets;
58  //std::cout<<"number of jets: "<<jets.size()<<std::endl;
59  for(std::vector<pat::Jet>::const_iterator it = jets.begin(); it != jets.end(); it++) {
60  //std::cout<<"Jet Et: "<<it->et()<<" Eta: "<<fabs(it->eta())<<std::endl;
62  if(it->et()>30. && fabs(it->eta())<2.4) {
63  seljets.push_back(*it);
64  nJets++;
65  }
66  }
67  //std::cout<<"selected Jets: "<<nJets<<std::endl;
68  if(nJets<4) return;
69 
70  //sort by Pt
71  sort(seljets.begin(),seljets.end(),JetwithHigherPt());
72 
73  edm::Handle< edm::View<pat::Muon> > muon_handle;
74  evt.getByLabel(muons_, muon_handle);
75  if(!muon_handle.isValid()) return;
76  const edm::View<pat::Muon> muons = *muon_handle;
77  int nmuons = 0;
78  std::vector<pat::Muon> selMuons;
79  for(edm::View<pat::Muon>::const_iterator it = muons.begin(); it!=muons.end(); it++) {
80  reco::TrackRef gltr = it->track(); // global track
81  reco::TrackRef trtr = it->innerTrack(); // tracker track
82  if(it->pt()>30 && fabs(it->eta())<2.1 && (it->pt()/(it->pt()+it->trackIso()+it->caloIso()))>0.95 && it->isGlobalMuon()){
83  if(gltr.isNull()) continue; //temporary problems with dead trackrefs
84  if((gltr->chi2()/gltr->ndof())<10 && trtr->numberOfValidHits()>=11) {
85 
86  double dRmin = 9999.;
87  for(std::vector<pat::Jet>::const_iterator ajet = seljets.begin(); ajet != seljets.end(); ajet++) {
88  math::XYZTLorentzVector jet = ajet->p4();
89  math::XYZTLorentzVector muon = it->p4();
90  double tmpdR = DeltaR(muon,jet);
91  if(tmpdR<dRmin) dRmin = tmpdR;
92  }
93  reco::TrackRef trtr = it->track(); // tracker track
94  if(dRmin>0.3) { //temporary problems with muon isolation
95  nmuons++;
96  selMuons.push_back(*it);
97  }
98  }
99  }
100  }
101  //std::cout<<"selected Muons: "<<nleptons<<std::endl;
102  if(nmuons!=1) return;
103 
104  edm::Handle< edm::View<pat::Electron> > electron_handle;
105  evt.getByLabel(electrons_, electron_handle);
106  if(!electron_handle.isValid()) return;
107  const edm::View<pat::Electron> electrons = *electron_handle;
108  int nelectrons = 0;
109  for(edm::View<pat::Electron>::const_iterator it = electrons.begin(); it!=electrons.end(); it++) {
110  if(it->pt()>30 && fabs(it->eta())<2.4 && (it->pt()/(it->pt()+it->trackIso()+it->caloIso()))>0.95 && it->isElectronIDAvailable("eidTight"))
111  {
112  if(it->electronID("eidTight")==1) nelectrons++;
113  }
114  }
115  if(nelectrons>0) return;
116  //end of the preselection
117 
118 
119  math::XYZTLorentzVector muon = selMuons.begin()->p4();
120 
121  //count the number of selected events
122  selEv++;
123  //skip event if enough events are already selected
124  if(selEv>maxEv_ && maxEv_!=-1) return;
125 
126  //calculation of InputVariables
127  //see TopQuarkAnalysis/TopTools/interface/TtSemiLepSignalSel.h
128  // /src/TtSemiLepSignalSel.cc
129  //all objects i.e. jets, muons, electrons... which are needed for the calculation
130  //of the input-variables have to be passed to this class
131 
132  TtSemiLepSignalSel selection(seljets,muon,MET);
133 
134  //this is only needed for the default example
136  evt.getByLabel("genEvt", genEvt);
137 
138  double weight = 1.0; //standard no weight, i.e. weight=1.0, set this to the corresponding weight if
139  //different weights for different events are available
140  if(whatData_==-1) { //your training-file contains both, signal and background events
141  bool isSignal;
142  isSignal = true;//true for signal, false for background this has to be derived in some way
143  evaluateTtSemiLepSignalSel(mvaComputer, selection, weight, true, isSignal);
144  }
145  else {
146 
147  if(whatData_==1){ //your tree contains only signal events
148  //if needed do a special signal selection here
149  //the following code is for the default example
150  if(genEvt->isSemiLeptonic() && genEvt->semiLeptonicChannel() == lepChannel_) {
151  //std::cout<<"a tt_semlep muon event"<<std::endl;
152  evaluateTtSemiLepSignalSel(mvaComputer, selection, weight, true, true);
153  }
154  else selEv--;
155  }
156  else if(whatData_==0){
157  //std::cout<<"a Wjets event"<<std::endl;
158  evaluateTtSemiLepSignalSel(mvaComputer, selection, weight, true, false);
159  }
160  else std::cout<<"Config File Error!! Please check <whatData> in TtSemiLepSignalSelMVATrainer.cfi";
161  }
162 }
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
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:30
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:30
double DeltaR(math::XYZTLorentzVector v1, math::XYZTLorentzVector v2)
bool isNull() const
Checks for null.
Definition: Ref.h:244
bool isValid() const
Definition: HandleBase.h:76
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:359
bool update(const Calibration::MVAComputer *computer)
const_iterator begin() const
tuple cout
Definition: gather_cfg.py:41
const_iterator end() const
void TtSemiLepSignalSelMVATrainer::beginJob ( void  )
privatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 164 of file TtSemiLepSignalSelMVATrainer.cc.

References gather_cfg::cout, selEv, and whatData_.

164  {
165  selEv = 0;
166  if(whatData_!=-1 && whatData_!=0 && whatData_!=1){
167  std::cout<<"Config File Error!! Please check <whatData> in TtSemiLepSignalSelMVATrainer.cfi"<<std::endl;;
168  return;
169  }
170 }
tuple cout
Definition: gather_cfg.py:41
double TtSemiLepSignalSelMVATrainer::DeltaPhi ( math::XYZTLorentzVector  v1,
math::XYZTLorentzVector  v2 
)
private

Definition at line 172 of file TtSemiLepSignalSelMVATrainer.cc.

References dPhi(), and Pi.

Referenced by DeltaR().

173 {
174  double dPhi = fabs(v1.Phi() - v2.Phi());
175  if (dPhi > TMath::Pi()) dPhi = 2*TMath::Pi() - dPhi;
176  return dPhi;
177 }
const double Pi
double dPhi(double phi1, double phi2)
Definition: JetUtil.h:30
double TtSemiLepSignalSelMVATrainer::DeltaR ( math::XYZTLorentzVector  v1,
math::XYZTLorentzVector  v2 
)
private

Definition at line 179 of file TtSemiLepSignalSelMVATrainer.cc.

References DeltaPhi(), and dPhi().

Referenced by analyze().

180 {
181  double dPhi = DeltaPhi(v1,v2);
182  double dR = TMath::Sqrt((v1.Eta()-v2.Eta())*(v1.Eta()-v2.Eta())+dPhi*dPhi);
183  return dR;
184 }
double dPhi(double phi1, double phi2)
Definition: JetUtil.h:30
double DeltaPhi(math::XYZTLorentzVector v1, math::XYZTLorentzVector v2)

Member Data Documentation

edm::InputTag TtSemiLepSignalSelMVATrainer::electrons_
private

Definition at line 47 of file TtSemiLepSignalSelMVATrainer.h.

Referenced by analyze().

edm::InputTag TtSemiLepSignalSelMVATrainer::jets_
private

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().

edm::InputTag TtSemiLepSignalSelMVATrainer::METs_
private

Definition at line 49 of file TtSemiLepSignalSelMVATrainer.h.

Referenced by analyze().

edm::InputTag TtSemiLepSignalSelMVATrainer::muons_
private

Definition at line 46 of file TtSemiLepSignalSelMVATrainer.h.

Referenced by analyze().

PhysicsTools::MVAComputerCache TtSemiLepSignalSelMVATrainer::mvaComputer
private

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().