CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TtDilepEvtSolutionMaker.cc
Go to the documentation of this file.
1 //
2 //
3 
6 
8 
11 
12 #include <memory>
13 #include <vector>
14 
16  // configurables
17  electronSourceToken_ = consumes<std::vector<pat::Electron> >(iConfig.getParameter<edm::InputTag>("electronSource"));
18  muonSourceToken_ = consumes<std::vector<pat::Muon> >(iConfig.getParameter<edm::InputTag>("muonSource"));
19  tauSourceToken_ = consumes<std::vector<pat::Tau> >(iConfig.getParameter<edm::InputTag>("tauSource"));
20  metSourceToken_ = consumes<std::vector<pat::MET> >(iConfig.getParameter<edm::InputTag>("metSource"));
21  jetSourceToken_ = consumes<std::vector<pat::Jet> >(iConfig.getParameter<edm::InputTag>("jetSource"));
22  jetCorrScheme_ = iConfig.getParameter<int>("jetCorrectionScheme");
23  evtSourceToken_ = mayConsume<TtGenEvent>(iConfig.getParameter<edm::InputTag>("evtSource"));
24  nrCombJets_ = iConfig.getParameter<unsigned int>("nrCombJets");
25  matchToGenEvt_ = iConfig.getParameter<bool>("matchToGenEvt");
26  calcTopMass_ = iConfig.getParameter<bool>("calcTopMass");
27  useMCforBest_ = iConfig.getParameter<bool>("bestSolFromMC");
28  eeChannel_ = iConfig.getParameter<bool>("eeChannel");
29  emuChannel_ = iConfig.getParameter<bool>("emuChannel");
30  mumuChannel_ = iConfig.getParameter<bool>("mumuChannel");
31  mutauChannel_ = iConfig.getParameter<bool>("mutauChannel");
32  etauChannel_ = iConfig.getParameter<bool>("etauChannel");
33  tautauChannel_ = iConfig.getParameter<bool>("tautauChannel");
34  tmassbegin_ = iConfig.getParameter<double>("tmassbegin");
35  tmassend_ = iConfig.getParameter<double>("tmassend");
36  tmassstep_ = iConfig.getParameter<double>("tmassstep");
37  nupars_ = iConfig.getParameter<std::vector<double> >("neutrino_parameters");
38 
39  // define what will be produced
40  produces<std::vector<TtDilepEvtSolution> >();
41 
43 }
44 
46 
49 }
50 
53  iEvent.getByToken(tauSourceToken_, taus);
55  iEvent.getByToken(muonSourceToken_, muons);
57  iEvent.getByToken(electronSourceToken_, electrons);
59  iEvent.getByToken(metSourceToken_, mets);
61  iEvent.getByToken(jetSourceToken_, jets);
62 
63  int selMuonp = -1, selMuonm = -1;
64  int selElectronp = -1, selElectronm = -1;
65  int selTaup = -1, selTaum = -1;
66  bool leptonFound = false;
67  bool mumu = false;
68  bool emu = false;
69  bool ee = false;
70  bool etau = false;
71  bool mutau = false;
72  bool tautau = false;
73  bool leptonFoundEE = false;
74  bool leptonFoundMM = false;
75  bool leptonFoundTT = false;
76  bool leptonFoundEpMm = false;
77  bool leptonFoundEmMp = false;
78  bool leptonFoundEpTm = false;
79  bool leptonFoundEmTp = false;
80  bool leptonFoundMpTm = false;
81  bool leptonFoundMmTp = false;
82  bool jetsFound = false;
83  bool METFound = false;
84  std::vector<int> JetVetoByTaus;
85 
86  //select MET (TopMET vector is sorted on ET)
87  if (!mets->empty()) {
88  METFound = true;
89  }
90 
91  // If we have electrons and muons available,
92  // build a solutions with electrons and muons.
93  if (muons->size() + electrons->size() >= 2) {
94  // select leptons
95  if (electrons->empty())
96  mumu = true;
97  else if (muons->empty())
98  ee = true;
99  else if (electrons->size() == 1) {
100  if (muons->size() == 1)
101  emu = true;
102  else if (PTComp(&(*electrons)[0], &(*muons)[1]))
103  emu = true;
104  else
105  mumu = true;
106  } else if (electrons->size() > 1) {
107  if (PTComp(&(*electrons)[1], &(*muons)[0]))
108  ee = true;
109  else if (muons->size() == 1)
110  emu = true;
111  else if (PTComp(&(*electrons)[0], &(*muons)[1]))
112  emu = true;
113  else
114  mumu = true;
115  }
116  if (ee) {
117  if (LepDiffCharge(&(*electrons)[0], &(*electrons)[1])) {
118  leptonFound = true;
119  leptonFoundEE = true;
120  if (HasPositiveCharge(&(*electrons)[0])) {
121  selElectronp = 0;
122  selElectronm = 1;
123  } else {
124  selElectronp = 1;
125  selElectronm = 0;
126  }
127  }
128  } else if (emu) {
129  if (LepDiffCharge(&(*electrons)[0], &(*muons)[0])) {
130  leptonFound = true;
131  if (HasPositiveCharge(&(*electrons)[0])) {
132  leptonFoundEpMm = true;
133  selElectronp = 0;
134  selMuonm = 0;
135  } else {
136  leptonFoundEmMp = true;
137  selMuonp = 0;
138  selElectronm = 0;
139  }
140  }
141  } else if (mumu) {
142  if (LepDiffCharge(&(*muons)[0], &(*muons)[1])) {
143  leptonFound = true;
144  leptonFoundMM = true;
145  if (HasPositiveCharge(&(*muons)[0])) {
146  selMuonp = 0;
147  selMuonm = 1;
148  } else {
149  selMuonp = 1;
150  selMuonm = 0;
151  }
152  }
153  }
154  //select Jets (TopJet vector is sorted on ET)
155  if (jets->size() >= 2) {
156  jetsFound = true;
157  }
158  }
159  // If a tau is needed to have two leptons, then only consider the taus.
160  // This is the minimal modification of the dilept selection that includes taus,
161  // since we are considering taus only when no other solution exist.
162  else if (muons->size() + electrons->size() == 1 && !taus->empty()) {
163  // select leptons
164  if (muons->size() == 1) {
165  mutau = true;
166  // depending on the muon charge, set the right muon index and specify channel
167  int expectedCharge = -muons->begin()->charge();
168  int* tauIdx = nullptr;
169  if (expectedCharge < 0) {
170  selMuonp = 0;
171  tauIdx = &selTaum;
172  leptonFoundMpTm = true;
173  } else {
174  selMuonm = 0;
175  tauIdx = &selTaup;
176  leptonFoundMmTp = true;
177  }
178  // loop over the vector of taus to find the ones
179  // that have the charge opposite to the muon one, and do not match in eta-phi
180  std::vector<std::vector<pat::Tau>::const_iterator> subset1;
181  for (std::vector<pat::Tau>::const_iterator tau = taus->begin(); tau < taus->end(); ++tau) {
182  if (tau->charge() * expectedCharge >= 0 && DeltaR<pat::Particle>()(*tau, *(muons->begin())) > 0.1) {
183  *tauIdx = tau - taus->begin();
184  leptonFound = true;
185  subset1.push_back(tau);
186  }
187  }
188  // if there are more than one tau with ecalIsol==0, take the smallest E/P
189  float iso = 999.;
190  for (std::vector<std::vector<pat::Tau>::const_iterator>::const_iterator tau = subset1.begin();
191  tau < subset1.end();
192  ++tau) {
193  if ((*tau)->isPFTau() && (*tau)->isolationPFChargedHadrCandsPtSum() < iso) {
194  *tauIdx = *tau - taus->begin();
195  iso = (*tau)->isolationPFChargedHadrCandsPtSum();
196  }
197  }
198 
199  // check that one combination has been found
200  if (!leptonFound) {
201  leptonFoundMpTm = false;
202  leptonFoundMmTp = false;
203  }
204  // discard the jet that matches the tau (if one)
205  if (leptonFound) {
206  for (std::vector<pat::Jet>::const_iterator jet = jets->begin(); jet < jets->end(); ++jet) {
207  if (DeltaR<pat::Particle, pat::Jet>()(*(taus->begin() + *tauIdx), *jet) < 0.1) {
208  JetVetoByTaus.push_back(jet - jets->begin());
209  }
210  }
211  }
212  } else {
213  etau = true;
214  // depending on the electron charge, set the right electron index and specify channel
215  int expectedCharge = -electrons->begin()->charge();
216  int* tauIdx = nullptr;
217  if (expectedCharge < 0) {
218  selElectronp = 0;
219  tauIdx = &selTaum;
220  leptonFoundEpTm = true;
221  } else {
222  selElectronm = 0;
223  tauIdx = &selTaup;
224  leptonFoundEmTp = true;
225  }
226  // loop over the vector of taus to find the ones
227  // that have the charge opposite to the muon one, and do not match in eta-phi
228  std::vector<std::vector<pat::Tau>::const_iterator> subset1;
229  for (std::vector<pat::Tau>::const_iterator tau = taus->begin(); tau < taus->end(); ++tau) {
230  if (tau->charge() * expectedCharge >= 0 && DeltaR<pat::Particle>()(*tau, *(electrons->begin())) > 0.1) {
231  *tauIdx = tau - taus->begin();
232  leptonFound = true;
233  subset1.push_back(tau);
234  }
235  }
236  // if there are more than one tau with ecalIsol==0, take the smallest E/P
237  float iso = 999.;
238  for (std::vector<std::vector<pat::Tau>::const_iterator>::const_iterator tau = subset1.begin();
239  tau < subset1.end();
240  ++tau) {
241  if ((*tau)->isPFTau() && (*tau)->isolationPFChargedHadrCandsPtSum() < iso) {
242  *tauIdx = *tau - taus->begin();
243  iso = (*tau)->isolationPFChargedHadrCandsPtSum();
244  }
245  }
246 
247  // check that one combination has been found
248  if (!leptonFound) {
249  leptonFoundEpTm = false;
250  leptonFoundEmTp = false;
251  }
252  // discard the jet that matches the tau (if one)
253  if (leptonFound) {
254  for (std::vector<pat::Jet>::const_iterator jet = jets->begin(); jet < jets->end(); ++jet) {
255  if (DeltaR<pat::Particle, pat::Jet>()(*(taus->begin() + *tauIdx), *jet) < 0.1) {
256  JetVetoByTaus.push_back(jet - jets->begin());
257  }
258  }
259  }
260  }
261  // select Jets (TopJet vector is sorted on ET)
262  jetsFound = ((jets->size() - JetVetoByTaus.size()) >= 2);
263  } else if (taus->size() > 1) {
264  tautau = true;
265  if (LepDiffCharge(&(*taus)[0], &(*taus)[1])) {
266  leptonFound = true;
267  leptonFoundTT = true;
268  if (HasPositiveCharge(&(*taus)[0])) {
269  selTaup = 0;
270  selTaum = 1;
271  } else {
272  selTaup = 1;
273  selTaum = 0;
274  }
275  }
276  for (std::vector<pat::Jet>::const_iterator jet = jets->begin(); jet < jets->end(); ++jet) {
277  if (DeltaR<pat::Particle, pat::Jet>()((*taus)[0], *jet) < 0.1 ||
278  DeltaR<pat::Particle, pat::Jet>()((*taus)[1], *jet) < 0.1) {
279  JetVetoByTaus.push_back(jet - jets->begin());
280  }
281  }
282  // select Jets (TopJet vector is sorted on ET)
283  jetsFound = ((jets->size() - JetVetoByTaus.size()) >= 2);
284  }
285 
286  // Check that the above work makes sense
287  if (int(ee) + int(emu) + int(mumu) + int(etau) + int(mutau) + int(tautau) > 1)
288  std::cout << "[TtDilepEvtSolutionMaker]: "
289  << "Lepton selection criteria uncorrectly defined" << std::endl;
290 
291  bool correctLepton = (leptonFoundEE && eeChannel_) || ((leptonFoundEmMp || leptonFoundEpMm) && emuChannel_) ||
292  (leptonFoundMM && mumuChannel_) || ((leptonFoundMmTp || leptonFoundMpTm) && mutauChannel_) ||
293  ((leptonFoundEmTp || leptonFoundEpTm) && etauChannel_) || (leptonFoundTT && tautauChannel_);
294 
295  std::vector<TtDilepEvtSolution>* evtsols = new std::vector<TtDilepEvtSolution>();
296  if (correctLepton && METFound && jetsFound) {
297  // protect against reading beyond array boundaries while discounting vetoed jets
298  unsigned int nrCombJets = 0;
299  unsigned int numberOfJets = 0;
300  for (; nrCombJets < jets->size() && numberOfJets < nrCombJets_; ++nrCombJets) {
301  if (find(JetVetoByTaus.begin(), JetVetoByTaus.end(), int(nrCombJets)) == JetVetoByTaus.end())
302  ++numberOfJets;
303  }
304  // consider all permutations
305  for (unsigned int ib = 0; ib < nrCombJets; ib++) {
306  // skipped jet vetoed during components-flagging.
307  if (find(JetVetoByTaus.begin(), JetVetoByTaus.end(), int(ib)) != JetVetoByTaus.end())
308  continue;
309  // second loop of the permutations
310  for (unsigned int ibbar = 0; ibbar < nrCombJets; ibbar++) {
311  // avoid the diagonal: b and bbar must be distinct jets
312  if (ib == ibbar)
313  continue;
314  // skipped jet vetoed during components-flagging.
315  if (find(JetVetoByTaus.begin(), JetVetoByTaus.end(), int(ibbar)) != JetVetoByTaus.end())
316  continue;
317  // Build and save a solution
318  TtDilepEvtSolution asol;
320  double xconstraint = 0, yconstraint = 0;
321  // Set e+ in the event
322  if (leptonFoundEE || leptonFoundEpMm || leptonFoundEpTm) {
323  asol.setElectronp(electrons, selElectronp);
324  xconstraint += (*electrons)[selElectronp].px();
325  yconstraint += (*electrons)[selElectronp].py();
326  }
327  // Set e- in the event
328  if (leptonFoundEE || leptonFoundEmMp || leptonFoundEmTp) {
329  asol.setElectronm(electrons, selElectronm);
330  xconstraint += (*electrons)[selElectronm].px();
331  yconstraint += (*electrons)[selElectronm].py();
332  }
333  // Set mu+ in the event
334  if (leptonFoundMM || leptonFoundEmMp || leptonFoundMpTm) {
335  asol.setMuonp(muons, selMuonp);
336  xconstraint += (*muons)[selMuonp].px();
337  yconstraint += (*muons)[selMuonp].py();
338  }
339  // Set mu- in the event
340  if (leptonFoundMM || leptonFoundEpMm || leptonFoundMmTp) {
341  asol.setMuonm(muons, selMuonm);
342  xconstraint += (*muons)[selMuonm].px();
343  yconstraint += (*muons)[selMuonm].py();
344  }
345  // Set tau- in the event
346  if (leptonFoundEpTm || leptonFoundMpTm || leptonFoundTT) {
347  asol.setTaum(taus, selTaum);
348  xconstraint += (*taus)[selTaum].px();
349  yconstraint += (*taus)[selTaum].py();
350  }
351  // Set tau+ in the event
352  if (leptonFoundEmTp || leptonFoundMmTp || leptonFoundTT) {
353  asol.setTaup(taus, selTaup);
354  xconstraint += (*taus)[selTaup].px();
355  yconstraint += (*taus)[selTaup].py();
356  }
357  // Set Jets/MET in the event
358  asol.setB(jets, ib);
359  asol.setBbar(jets, ibbar);
360  asol.setMET(mets, 0);
361  xconstraint += (*jets)[ib].px() + (*jets)[ibbar].px() + (*mets)[0].px();
362  yconstraint += (*jets)[ib].py() + (*jets)[ibbar].py() + (*mets)[0].py();
363  // if asked for, match the event solutions to the gen Event
364  if (matchToGenEvt_) {
366  iEvent.getByToken(evtSourceToken_, genEvt);
367  asol.setGenEvt(genEvt);
368  }
369  // If asked, use the kin fitter to compute the top mass
370  if (calcTopMass_) {
371  solver->SetConstraints(xconstraint, yconstraint);
373  asol = solver->addKinSolInfo(&asol);
374  }
375 
376  // these lines calculate the observables to be used in the TtDilepSignalSelection LR
377  (*myLRSignalSelObservables)(asol, iEvent);
378 
379  evtsols->push_back(asol);
380  }
381  }
382  // flag the best solution (MC matching)
383  if (matchToGenEvt_) {
384  double bestSolDR = 9999.;
385  int bestSol = -1;
386  double dR = 0.;
387  for (size_t s = 0; s < evtsols->size(); s++) {
388  dR = (*evtsols)[s].getJetResidual();
389  if (dR < bestSolDR) {
390  bestSolDR = dR;
391  bestSol = s;
392  }
393  }
394  if (bestSol != -1)
395  (*evtsols)[bestSol].setBestSol(true);
396  }
397  // put the result in the event
398  std::unique_ptr<std::vector<TtDilepEvtSolution> > pOut(evtsols);
399  iEvent.put(std::move(pOut));
400  } else {
401  // no solution: put a dummy solution in the event
402  TtDilepEvtSolution asol;
403  evtsols->push_back(asol);
404  std::unique_ptr<std::vector<TtDilepEvtSolution> > pOut(evtsols);
405  iEvent.put(std::move(pOut));
406  }
407 }
edm::EDGetTokenT< std::vector< pat::Jet > > jetSourceToken_
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
int ib
Definition: cuy.py:661
void setElectronm(const edm::Handle< std::vector< pat::Electron > > &elec, int i)
Definition: deltaR.h:58
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
void setElectronp(const edm::Handle< std::vector< pat::Electron > > &elec, int i)
void setB(const edm::Handle< std::vector< pat::Jet > > &jet, int i)
edm::EDGetTokenT< std::vector< pat::Muon > > muonSourceToken_
bool HasPositiveCharge(const reco::Candidate *) const
void SetConstraints(const double xx=0, const double yy=0)
TtDilepEvtSolutionMaker(const edm::ParameterSet &iConfig)
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
edm::EDGetTokenT< std::vector< pat::Tau > > tauSourceToken_
edm::EDGetTokenT< std::vector< pat::MET > > metSourceToken_
void setMuonm(const edm::Handle< std::vector< pat::Muon > > &muon, int i)
edm::EDGetTokenT< std::vector< pat::Electron > > electronSourceToken_
int iEvent
Definition: GenABIO.cc:224
bool PTComp(const reco::Candidate *, const reco::Candidate *) const
TtDilepEvtSolution addKinSolInfo(TtDilepEvtSolution *asol)
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
vector< PseudoJet > jets
void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override
def move
Definition: eostools.py:511
void setMET(const edm::Handle< std::vector< pat::MET > > &met, int i)
void useWeightFromMC(bool useMC)
void setTaup(const edm::Handle< std::vector< pat::Tau > > &tau, int i)
void setJetCorrectionScheme(int jetCorrScheme)
edm::EDGetTokenT< TtGenEvent > evtSourceToken_
void setGenEvt(const edm::Handle< TtGenEvent > &)
void setMuonp(const edm::Handle< std::vector< pat::Muon > > &muon, int i)
std::vector< double > nupars_
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
bool LepDiffCharge(const reco::Candidate *, const reco::Candidate *) const
double getJetResidual() const
tuple muons
Definition: patZpeak.py:39
void setTaum(const edm::Handle< std::vector< pat::Tau > > &tau, int i)
tuple cout
Definition: gather_cfg.py:144
TtDilepLRSignalSelObservables * myLRSignalSelObservables
void setBbar(const edm::Handle< std::vector< pat::Jet > > &jet, int i)