CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TtSemiLepSignalSelMVAComputer.cc
Go to the documentation of this file.
2 
5 
7 
12 
14  : mvaToken_(esConsumes()),
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  putToken_(produces("DiscSel")) {}
20 
22 
24  mvaComputer.update(&setup.getData(mvaToken_), "ttSemiLepSignalSelMVA");
25 
26  //make your preselection! This must!! be the same one as in TraintreeSaver.cc
28  evt.getByToken(METsToken_, MET_handle);
29  if (!MET_handle.isValid())
30  return;
31  const edm::View<pat::MET> MET = *MET_handle;
32 
34  evt.getByToken(jetsToken_, jet_handle);
35  if (!jet_handle.isValid())
36  return;
37  const std::vector<pat::Jet> jets = *jet_handle;
38  unsigned int nJets = 0;
39  std::vector<pat::Jet> seljets;
40  for (std::vector<pat::Jet>::const_iterator it = jets.begin(); it != jets.end(); it++) {
42  continue;
43  if (it->pt() > 30. && fabs(it->eta()) < 2.4) {
44  seljets.push_back(*it);
45  nJets++;
46  }
47  }
48 
50  evt.getByToken(muonsToken_, muon_handle);
51  if (!muon_handle.isValid())
52  return;
53  const edm::View<pat::Muon> muons = *muon_handle;
54  int nmuons = 0;
55  std::vector<pat::Muon> selMuons;
56  for (edm::View<pat::Muon>::const_iterator it = muons.begin(); it != muons.end(); it++) {
57  reco::TrackRef gltr = it->track(); // global track
58  reco::TrackRef trtr = it->innerTrack(); // tracker track
59  if (it->pt() > 30 && fabs(it->eta()) < 2.1 && (it->pt() / (it->pt() + it->trackIso() + it->caloIso())) > 0.95 &&
60  it->isGlobalMuon()) {
61  if (gltr.isNull())
62  continue; //temporary problems with dead trackrefs
63  if ((gltr->chi2() / gltr->ndof()) < 10 && trtr->numberOfValidHits() > 11) {
64  double dRmin = 9999.;
65  for (std::vector<pat::Jet>::const_iterator ajet = seljets.begin(); ajet != seljets.end(); ajet++) {
66  math::XYZTLorentzVector jet = ajet->p4();
67  math::XYZTLorentzVector muon = it->p4();
68  double tmpdR = DeltaR(muon, jet);
69  if (tmpdR < dRmin)
70  dRmin = tmpdR;
71  }
72  if (dRmin > 0.3) { //temporary problems with muon isolation
73  nmuons++;
74  selMuons.push_back(*it);
75  }
76  }
77  }
78  }
79 
80  edm::Handle<edm::View<pat::Electron> > electron_handle;
81  evt.getByToken(electronsToken_, electron_handle);
82  if (!electron_handle.isValid())
83  return;
84  const edm::View<pat::Electron> electrons = *electron_handle;
85  int nelectrons = 0;
86  for (edm::View<pat::Electron>::const_iterator it = electrons.begin(); it != electrons.end(); it++) {
87  if (it->pt() > 30 && fabs(it->eta()) < 2.4 && (it->pt() / (it->pt() + it->trackIso() + it->caloIso())) > 0.95 &&
88  it->isElectronIDAvailable("eidTight")) {
89  if (it->electronID("eidTight") == 1)
90  nelectrons++;
91  }
92  }
93 
94  double discrim;
95  // discriminator output for events which do not pass the preselection is set to -1
96  if (nmuons != 1 || nJets < 4 || nelectrons > 0)
97  discrim = -1.; //std::cout<<"nJets: "<<seljets.size()<<" numLeptons: "<<nleptons<<std::endl;}
98  else {
99  //check wheter a event was already selected (problem with duplicated events)
100  math::XYZTLorentzVector muon = selMuons.begin()->p4();
101 
102  TtSemiLepSignalSel selection(seljets, muon, MET);
103 
104  discrim = evaluateTtSemiLepSignalSel(mvaComputer, selection);
105  }
106 
107  evt.emplace(putToken_, discrim);
108 }
109 
111  double dPhi = fabs(v1.Phi() - v2.Phi());
112  if (dPhi > TMath::Pi())
113  dPhi = 2 * TMath::Pi() - dPhi;
114  return dPhi;
115 }
116 
118  double dPhi = DeltaPhi(v1, v2);
119  double dR = TMath::Sqrt((v1.Eta() - v2.Eta()) * (v1.Eta() - v2.Eta()) + dPhi * dPhi);
120  return dR;
121 }
122 
123 // implement the plugins for the computer container
124 // -> register TtSemiLepSignalSelMVARcd
125 // -> define TtSemiLepSignalSelMVAFileSource
126 MVA_COMPUTER_CONTAINER_IMPLEMENT(TtSemiLepSignalSelMVA);
const double Pi
tuple cfg
Definition: looper.py:296
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
double DeltaPhi(const math::XYZTLorentzVector &v1, const math::XYZTLorentzVector &v2)
static bool test(uint32_t val, uint32_t mask)
Definition: Flags.h:28
edm::EDGetTokenT< edm::View< pat::Electron > > electronsToken_
selection
main part
Definition: corrVsCorr.py:100
PhysicsTools::MVAComputerCache mvaComputer
bool getData(T &iHolder) const
Definition: EventSetup.h:128
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
const_iterator begin() const
edm::EDGetTokenT< edm::View< pat::Muon > > muonsToken_
double evaluateTtSemiLepSignalSel(PhysicsTools::MVAComputerCache &mvaComputer, const TtSemiLepSignalSel &sigsel, float weight=1., const bool isSignal=false)
vector< PseudoJet > jets
double DeltaR(const math::XYZTLorentzVector &v1, const math::XYZTLorentzVector &v2)
bool isValid() const
Definition: HandleBase.h:70
bool isNull() const
Checks for null.
Definition: Ref.h:235
edm::EDGetTokenT< edm::View< pat::MET > > METsToken_
bool update(const Calibration::MVAComputer *computer)
OrphanHandle< PROD > emplace(EDPutTokenT< PROD > token, Args &&...args)
puts a new product
Definition: Event.h:433
constexpr char Jet[]
Definition: modules.cc:9
tuple muons
Definition: patZpeak.py:39
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
edm::EDGetTokenT< std::vector< pat::Jet > > jetsToken_
TtSemiLepSignalSelMVAComputer(const edm::ParameterSet &)
const_iterator end() const
void produce(edm::Event &evt, const edm::EventSetup &setup) override
constexpr char Electron[]
Definition: modules.cc:12
ESGetTokenH3DDVariant esConsumes(std::string const &Reccord, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
edm::ESGetToken< PhysicsTools::Calibration::MVAComputerContainer, TtSemiLepSignalSelMVARcd > mvaToken_
#define MVA_COMPUTER_CONTAINER_IMPLEMENT(N)
Definition: HelperMacros.h:53