CMS 3D CMS Logo

TtSemiLepSignalSelMVAComputer.cc
Go to the documentation of this file.
2 
5 
7 
12 
13 
15  muonsToken_ (consumes< edm::View<pat::Muon> >(cfg.getParameter<edm::InputTag>("muons"))),
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  electronsToken_ (consumes< edm::View<pat::Electron> >(cfg.getParameter<edm::InputTag>("elecs")))
19 {
20  produces< double >("DiscSel");
21 }
22 
23 
24 
26 {
27 }
28 
29 void
31 {
32  std::unique_ptr< double > pOutDisc (new double);
33 
34  mvaComputer.update<TtSemiLepSignalSelMVARcd>(setup, "ttSemiLepSignalSelMVA");
35 
36  // read name of the last processor in the MVA calibration
37  // (to be used as meta information)
39  setup.get<TtSemiLepSignalSelMVARcd>().get( calibContainer );
40  std::vector<PhysicsTools::Calibration::VarProcessor*> processors
41  = (calibContainer->find("ttSemiLepSignalSelMVA")).getProcessors();
42 
43  //make your preselection! This must!! be the same one as in TraintreeSaver.cc
45  evt.getByToken(METsToken_,MET_handle);
46  if(!MET_handle.isValid()) return;
47  const edm::View<pat::MET> MET = *MET_handle;
48 
50  evt.getByToken(jetsToken_, jet_handle);
51  if(!jet_handle.isValid()) return;
52  const std::vector<pat::Jet> jets = *jet_handle;
53  unsigned int nJets = 0;
54  std::vector<pat::Jet> seljets;
55  for(std::vector<pat::Jet>::const_iterator it = jets.begin(); it != jets.end(); it++) {
57  if(it->pt()>30. && fabs(it->eta())<2.4) {
58  seljets.push_back(*it);
59  nJets++;
60  }
61  }
62 
64  evt.getByToken(muonsToken_, muon_handle);
65  if(!muon_handle.isValid()) return;
66  const edm::View<pat::Muon> muons = *muon_handle;
67  int nmuons = 0;
68  std::vector<pat::Muon> selMuons;
69  for(edm::View<pat::Muon>::const_iterator it = muons.begin(); it!=muons.end(); it++) {
70  reco::TrackRef gltr = it->track(); // global track
71  reco::TrackRef trtr = it->innerTrack(); // tracker track
72  if(it->pt()>30 && fabs(it->eta())<2.1 && (it->pt()/(it->pt()+it->trackIso()+it->caloIso()))>0.95 && it->isGlobalMuon()){
73  if(gltr.isNull()) continue; //temporary problems with dead trackrefs
74  if((gltr->chi2()/gltr->ndof())<10 && trtr->numberOfValidHits()>11) {
75  double dRmin = 9999.;
76  for(std::vector<pat::Jet>::const_iterator ajet = seljets.begin(); ajet != seljets.end(); ajet++) {
77  math::XYZTLorentzVector jet = ajet->p4();
78  math::XYZTLorentzVector muon = it->p4();
79  double tmpdR = DeltaR(muon,jet);
80  if(tmpdR<dRmin) dRmin = tmpdR;
81  }
82  if(dRmin>0.3) { //temporary problems with muon isolation
83  nmuons++;
84  selMuons.push_back(*it);
85  }
86  }
87  }
88  }
89 
90  edm::Handle< edm::View<pat::Electron> > electron_handle;
91  evt.getByToken(electronsToken_, electron_handle);
92  if(!electron_handle.isValid()) return;
93  const edm::View<pat::Electron> electrons = *electron_handle;
94  int nelectrons = 0;
95  for(edm::View<pat::Electron>::const_iterator it = electrons.begin(); it!=electrons.end(); it++) {
96  if(it->pt()>30 && fabs(it->eta())<2.4 && (it->pt()/(it->pt()+it->trackIso()+it->caloIso()))>0.95 && it->isElectronIDAvailable("eidTight"))
97  {
98  if(it->electronID("eidTight")==1) nelectrons++;
99  }
100  }
101 
102 
103  double discrim;
104  // discriminator output for events which do not pass the preselection is set to -1
105  if( nmuons!=1 ||
106  nJets < 4 ||
107  nelectrons > 0 ) discrim = -1.; //std::cout<<"nJets: "<<seljets.size()<<" numLeptons: "<<nleptons<<std::endl;}
108  else {
109  //check wheter a event was already selected (problem with duplicated events)
110  math::XYZTLorentzVector muon = selMuons.begin()->p4();
111 
112  TtSemiLepSignalSel selection(seljets,muon,MET);
113 
114  discrim = evaluateTtSemiLepSignalSel(mvaComputer, selection);
115  }
116 
117  *pOutDisc = discrim;
118 
119  evt.put(std::move(pOutDisc), "DiscSel");
120 
121  DiscSel = discrim;
122 
123 }
124 
125 void
127 {
128 }
129 
130 void
132 {
133 }
134 
136 {
137  double dPhi = fabs(v1.Phi() - v2.Phi());
138  if (dPhi > TMath::Pi()) dPhi = 2*TMath::Pi() - dPhi;
139  return dPhi;
140 }
141 
143 {
144  double dPhi = DeltaPhi(v1,v2);
145  double dR = TMath::Sqrt((v1.Eta()-v2.Eta())*(v1.Eta()-v2.Eta())+dPhi*dPhi);
146  return dR;
147 }
148 
149 // implement the plugins for the computer container
150 // -> register TtSemiLepSignalSelMVARcd
151 // -> define TtSemiLepSignalSelMVAFileSource
152 MVA_COMPUTER_CONTAINER_IMPLEMENT(TtSemiLepSignalSelMVA);
const double Pi
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:122
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:460
double evaluateTtSemiLepSignalSel(PhysicsTools::MVAComputerCache &mvaComputer, const TtSemiLepSignalSel &sigsel, float weight=1., const bool training=false, const bool isSignal=false)
double DeltaPhi(const math::XYZTLorentzVector &v1, const math::XYZTLorentzVector &v2)
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
PhysicsTools::MVAComputerCache mvaComputer
Definition: HeavyIon.h:7
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
const_iterator begin() const
Definition: Muon.py:1
Definition: Jet.py:1
vector< PseudoJet > jets
void get(HolderT &iHolder) const
double DeltaR(const math::XYZTLorentzVector &v1, const math::XYZTLorentzVector &v2)
edm::EDGetTokenT< edm::View< pat::Muon > > muonsToken_
virtual void produce(edm::Event &evt, const edm::EventSetup &setup)
bool isValid() const
Definition: HandleBase.h:74
bool isNull() const
Checks for null.
Definition: Ref.h:249
edm::EDGetTokenT< edm::View< pat::MET > > METsToken_
bool update(const Calibration::MVAComputer *computer)
edm::EDGetTokenT< std::vector< pat::Jet > > jetsToken_
#define MVA_COMPUTER_CONTAINER_IMPLEMENT(N)
Definition: HelperMacros.h:46
HLT enums.
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
TtSemiLepSignalSelMVAComputer(const edm::ParameterSet &)
edm::EDGetTokenT< edm::View< pat::Electron > > electronsToken_
const_iterator end() const
def move(src, dest)
Definition: eostools.py:510