CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
TtDilepEvtSolutionMaker Class Reference
Inheritance diagram for TtDilepEvtSolutionMaker:
edm::stream::EDProducer<>

Public Member Functions

void produce (edm::Event &iEvent, const edm::EventSetup &iSetup) override
 
 TtDilepEvtSolutionMaker (const edm::ParameterSet &iConfig)
 
 ~TtDilepEvtSolutionMaker () override
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
 EDProducer (const EDProducer &)=delete
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
const EDProduceroperator= (const EDProducer &)=delete
 

Private Member Functions

bool HasPositiveCharge (const reco::Candidate *) const
 
bool LepDiffCharge (const reco::Candidate *, const reco::Candidate *) const
 
bool PTComp (const reco::Candidate *, const reco::Candidate *) const
 

Private Attributes

bool calcTopMass_
 
bool eeChannel_
 
edm::EDGetTokenT< std::vector< pat::Electron > > electronSourceToken_
 
bool emuChannel_
 
bool etauChannel_
 
edm::EDGetTokenT< TtGenEventevtSourceToken_
 
int jetCorrScheme_
 
edm::EDGetTokenT< std::vector< pat::Jet > > jetSourceToken_
 
bool matchToGenEvt_
 
edm::EDGetTokenT< std::vector< pat::MET > > metSourceToken_
 
bool mumuChannel_
 
edm::EDGetTokenT< std::vector< pat::Muon > > muonSourceToken_
 
bool mutauChannel_
 
TtDilepLRSignalSelObservablesmyLRSignalSelObservables
 
unsigned int nrCombJets_
 
std::vector< double > nupars_
 
TtFullLepKinSolversolver
 
edm::EDGetTokenT< std::vector< pat::Tau > > tauSourceToken_
 
bool tautauChannel_
 
double tmassbegin_
 
double tmassend_
 
double tmassstep_
 
bool useMCforBest_
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer<>
using CacheTypes = CacheContexts< T... >
 
using GlobalCache = typename CacheTypes::GlobalCache
 
using HasAbility = AbilityChecker< T... >
 
using InputProcessBlockCache = typename CacheTypes::InputProcessBlockCache
 
using LuminosityBlockCache = typename CacheTypes::LuminosityBlockCache
 
using LuminosityBlockContext = LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCache >
 
using LuminosityBlockSummaryCache = typename CacheTypes::LuminosityBlockSummaryCache
 
using RunCache = typename CacheTypes::RunCache
 
using RunContext = RunContextT< RunCache, GlobalCache >
 
using RunSummaryCache = typename CacheTypes::RunSummaryCache
 

Detailed Description

Definition at line 18 of file TtDilepEvtSolutionMaker.cc.

Constructor & Destructor Documentation

◆ TtDilepEvtSolutionMaker()

TtDilepEvtSolutionMaker::TtDilepEvtSolutionMaker ( const edm::ParameterSet iConfig)
explicit

Definition at line 59 of file TtDilepEvtSolutionMaker.cc.

References calcTopMass_, eeChannel_, electronSourceToken_, emuChannel_, etauChannel_, evtSourceToken_, edm::ParameterSet::getParameter(), jetCorrScheme_, jetSourceToken_, matchToGenEvt_, metSourceToken_, mumuChannel_, muonSourceToken_, mutauChannel_, myLRSignalSelObservables, nrCombJets_, nupars_, solver, tauSourceToken_, tautauChannel_, tmassbegin_, tmassend_, tmassstep_, and useMCforBest_.

59  {
60  // configurables
61  electronSourceToken_ = consumes<std::vector<pat::Electron> >(iConfig.getParameter<edm::InputTag>("electronSource"));
62  muonSourceToken_ = consumes<std::vector<pat::Muon> >(iConfig.getParameter<edm::InputTag>("muonSource"));
63  tauSourceToken_ = consumes<std::vector<pat::Tau> >(iConfig.getParameter<edm::InputTag>("tauSource"));
64  metSourceToken_ = consumes<std::vector<pat::MET> >(iConfig.getParameter<edm::InputTag>("metSource"));
65  jetSourceToken_ = consumes<std::vector<pat::Jet> >(iConfig.getParameter<edm::InputTag>("jetSource"));
66  jetCorrScheme_ = iConfig.getParameter<int>("jetCorrectionScheme");
67  evtSourceToken_ = mayConsume<TtGenEvent>(iConfig.getParameter<edm::InputTag>("evtSource"));
68  nrCombJets_ = iConfig.getParameter<unsigned int>("nrCombJets");
69  matchToGenEvt_ = iConfig.getParameter<bool>("matchToGenEvt");
70  calcTopMass_ = iConfig.getParameter<bool>("calcTopMass");
71  useMCforBest_ = iConfig.getParameter<bool>("bestSolFromMC");
72  eeChannel_ = iConfig.getParameter<bool>("eeChannel");
73  emuChannel_ = iConfig.getParameter<bool>("emuChannel");
74  mumuChannel_ = iConfig.getParameter<bool>("mumuChannel");
75  mutauChannel_ = iConfig.getParameter<bool>("mutauChannel");
76  etauChannel_ = iConfig.getParameter<bool>("etauChannel");
77  tautauChannel_ = iConfig.getParameter<bool>("tautauChannel");
78  tmassbegin_ = iConfig.getParameter<double>("tmassbegin");
79  tmassend_ = iConfig.getParameter<double>("tmassend");
80  tmassstep_ = iConfig.getParameter<double>("tmassstep");
81  nupars_ = iConfig.getParameter<std::vector<double> >("neutrino_parameters");
82 
83  // define what will be produced
84  produces<std::vector<TtDilepEvtSolution> >();
85 
87 
89 }
edm::EDGetTokenT< std::vector< pat::Jet > > jetSourceToken_
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
edm::EDGetTokenT< std::vector< pat::Muon > > muonSourceToken_
edm::EDGetTokenT< std::vector< pat::Tau > > tauSourceToken_
edm::EDGetTokenT< std::vector< pat::MET > > metSourceToken_
edm::EDGetTokenT< std::vector< pat::Electron > > electronSourceToken_
edm::EDGetTokenT< TtGenEvent > evtSourceToken_
std::vector< double > nupars_
TtDilepLRSignalSelObservables * myLRSignalSelObservables

◆ ~TtDilepEvtSolutionMaker()

TtDilepEvtSolutionMaker::~TtDilepEvtSolutionMaker ( )
override

Definition at line 91 of file TtDilepEvtSolutionMaker.cc.

91 {}

Member Function Documentation

◆ HasPositiveCharge()

bool TtDilepEvtSolutionMaker::HasPositiveCharge ( const reco::Candidate l) const
inlineprivate

Definition at line 57 of file TtDilepEvtSolutionMaker.cc.

References MainPageGenerator::l.

Referenced by produce().

57 { return (l->charge() > 0); }

◆ LepDiffCharge()

bool TtDilepEvtSolutionMaker::LepDiffCharge ( const reco::Candidate l1,
const reco::Candidate l2 
) const
inlineprivate

Definition at line 53 of file TtDilepEvtSolutionMaker.cc.

References reco::Candidate::charge().

Referenced by produce().

53  {
54  return (l1->charge() != l2->charge());
55 }
virtual int charge() const =0
electric charge

◆ produce()

void TtDilepEvtSolutionMaker::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
override

Definition at line 93 of file TtDilepEvtSolutionMaker.cc.

References TtFullLepKinSolver::addKinSolInfo(), calcTopMass_, gather_cfg::cout, HGC3DClusterGenMatchSelector_cfi::dR, eeChannel_, pwdgSkimBPark_cfi::electrons, electronSourceToken_, emuChannel_, etauChannel_, evtSourceToken_, spr::find(), TtGenEvtProducer_cfi::genEvt, HasPositiveCharge(), cuy::ib, iEvent, createfilelist::int, genparticles_cff::iso, metsig::jet, jetCorrScheme_, PDWG_EXODelayedJetMET_cff::jets, jetSourceToken_, LepDiffCharge(), matchToGenEvt_, singleTopDQM_cfi::mets, metSourceToken_, eostools::move(), mumuChannel_, DiMuonV_cfg::muons, muonSourceToken_, mutauChannel_, TtDilepEvtSolProducer_cfi::nrCombJets, nrCombJets_, PTComp(), alignCSCRings::s, TtDilepEvtSolution::setB(), TtDilepEvtSolution::setBbar(), TtFullLepKinSolver::SetConstraints(), TtDilepEvtSolution::setElectronm(), TtDilepEvtSolution::setElectronp(), TtDilepEvtSolution::setGenEvt(), TtDilepEvtSolution::setJetCorrectionScheme(), TtDilepEvtSolution::setMET(), TtDilepEvtSolution::setMuonm(), TtDilepEvtSolution::setMuonp(), TtDilepEvtSolution::setTaum(), TtDilepEvtSolution::setTaup(), solver, metsig::tau, taus_cff::tauIdx, Tau3MuMonitor_cff::taus, tauSourceToken_, tautauChannel_, useMCforBest_, TtFullLepKinSolver::useWeightFromMC(), and trackerHitRTTI::vector.

93  {
95  iEvent.getByToken(tauSourceToken_, taus);
97  iEvent.getByToken(muonSourceToken_, muons);
101  iEvent.getByToken(metSourceToken_, mets);
103  iEvent.getByToken(jetSourceToken_, jets);
104 
105  int selMuonp = -1, selMuonm = -1;
106  int selElectronp = -1, selElectronm = -1;
107  int selTaup = -1, selTaum = -1;
108  bool leptonFound = false;
109  bool mumu = false;
110  bool emu = false;
111  bool ee = false;
112  bool etau = false;
113  bool mutau = false;
114  bool tautau = false;
115  bool leptonFoundEE = false;
116  bool leptonFoundMM = false;
117  bool leptonFoundTT = false;
118  bool leptonFoundEpMm = false;
119  bool leptonFoundEmMp = false;
120  bool leptonFoundEpTm = false;
121  bool leptonFoundEmTp = false;
122  bool leptonFoundMpTm = false;
123  bool leptonFoundMmTp = false;
124  bool jetsFound = false;
125  bool METFound = false;
126  std::vector<int> JetVetoByTaus;
127 
128  //select MET (TopMET vector is sorted on ET)
129  if (!mets->empty()) {
130  METFound = true;
131  }
132 
133  // If we have electrons and muons available,
134  // build a solutions with electrons and muons.
135  if (muons->size() + electrons->size() >= 2) {
136  // select leptons
137  if (electrons->empty())
138  mumu = true;
139  else if (muons->empty())
140  ee = true;
141  else if (electrons->size() == 1) {
142  if (muons->size() == 1)
143  emu = true;
144  else if (PTComp(&(*electrons)[0], &(*muons)[1]))
145  emu = true;
146  else
147  mumu = true;
148  } else if (electrons->size() > 1) {
149  if (PTComp(&(*electrons)[1], &(*muons)[0]))
150  ee = true;
151  else if (muons->size() == 1)
152  emu = true;
153  else if (PTComp(&(*electrons)[0], &(*muons)[1]))
154  emu = true;
155  else
156  mumu = true;
157  }
158  if (ee) {
159  if (LepDiffCharge(&(*electrons)[0], &(*electrons)[1])) {
160  leptonFound = true;
161  leptonFoundEE = true;
162  if (HasPositiveCharge(&(*electrons)[0])) {
163  selElectronp = 0;
164  selElectronm = 1;
165  } else {
166  selElectronp = 1;
167  selElectronm = 0;
168  }
169  }
170  } else if (emu) {
171  if (LepDiffCharge(&(*electrons)[0], &(*muons)[0])) {
172  leptonFound = true;
173  if (HasPositiveCharge(&(*electrons)[0])) {
174  leptonFoundEpMm = true;
175  selElectronp = 0;
176  selMuonm = 0;
177  } else {
178  leptonFoundEmMp = true;
179  selMuonp = 0;
180  selElectronm = 0;
181  }
182  }
183  } else if (mumu) {
184  if (LepDiffCharge(&(*muons)[0], &(*muons)[1])) {
185  leptonFound = true;
186  leptonFoundMM = true;
187  if (HasPositiveCharge(&(*muons)[0])) {
188  selMuonp = 0;
189  selMuonm = 1;
190  } else {
191  selMuonp = 1;
192  selMuonm = 0;
193  }
194  }
195  }
196  //select Jets (TopJet vector is sorted on ET)
197  if (jets->size() >= 2) {
198  jetsFound = true;
199  }
200  }
201  // If a tau is needed to have two leptons, then only consider the taus.
202  // This is the minimal modification of the dilept selection that includes taus,
203  // since we are considering taus only when no other solution exist.
204  else if (muons->size() + electrons->size() == 1 && !taus->empty()) {
205  // select leptons
206  if (muons->size() == 1) {
207  mutau = true;
208  // depending on the muon charge, set the right muon index and specify channel
209  int expectedCharge = -muons->begin()->charge();
210  int* tauIdx = nullptr;
211  if (expectedCharge < 0) {
212  selMuonp = 0;
213  tauIdx = &selTaum;
214  leptonFoundMpTm = true;
215  } else {
216  selMuonm = 0;
217  tauIdx = &selTaup;
218  leptonFoundMmTp = true;
219  }
220  // loop over the vector of taus to find the ones
221  // that have the charge opposite to the muon one, and do not match in eta-phi
222  std::vector<std::vector<pat::Tau>::const_iterator> subset1;
223  for (std::vector<pat::Tau>::const_iterator tau = taus->begin(); tau < taus->end(); ++tau) {
224  if (tau->charge() * expectedCharge >= 0 && DeltaR<pat::Particle>()(*tau, *(muons->begin())) > 0.1) {
225  *tauIdx = tau - taus->begin();
226  leptonFound = true;
227  subset1.push_back(tau);
228  }
229  }
230  // if there are more than one tau with ecalIsol==0, take the smallest E/P
231  float iso = 999.;
232  for (std::vector<std::vector<pat::Tau>::const_iterator>::const_iterator tau = subset1.begin();
233  tau < subset1.end();
234  ++tau) {
235  if ((*tau)->isPFTau() && (*tau)->isolationPFChargedHadrCandsPtSum() < iso) {
236  *tauIdx = *tau - taus->begin();
237  iso = (*tau)->isolationPFChargedHadrCandsPtSum();
238  }
239  }
240 
241  // check that one combination has been found
242  if (!leptonFound) {
243  leptonFoundMpTm = false;
244  leptonFoundMmTp = false;
245  }
246  // discard the jet that matches the tau (if one)
247  if (leptonFound) {
248  for (std::vector<pat::Jet>::const_iterator jet = jets->begin(); jet < jets->end(); ++jet) {
249  if (DeltaR<pat::Particle, pat::Jet>()(*(taus->begin() + *tauIdx), *jet) < 0.1) {
250  JetVetoByTaus.push_back(jet - jets->begin());
251  }
252  }
253  }
254  } else {
255  etau = true;
256  // depending on the electron charge, set the right electron index and specify channel
257  int expectedCharge = -electrons->begin()->charge();
258  int* tauIdx = nullptr;
259  if (expectedCharge < 0) {
260  selElectronp = 0;
261  tauIdx = &selTaum;
262  leptonFoundEpTm = true;
263  } else {
264  selElectronm = 0;
265  tauIdx = &selTaup;
266  leptonFoundEmTp = true;
267  }
268  // loop over the vector of taus to find the ones
269  // that have the charge opposite to the muon one, and do not match in eta-phi
270  std::vector<std::vector<pat::Tau>::const_iterator> subset1;
271  for (std::vector<pat::Tau>::const_iterator tau = taus->begin(); tau < taus->end(); ++tau) {
272  if (tau->charge() * expectedCharge >= 0 && DeltaR<pat::Particle>()(*tau, *(electrons->begin())) > 0.1) {
273  *tauIdx = tau - taus->begin();
274  leptonFound = true;
275  subset1.push_back(tau);
276  }
277  }
278  // if there are more than one tau with ecalIsol==0, take the smallest E/P
279  float iso = 999.;
280  for (std::vector<std::vector<pat::Tau>::const_iterator>::const_iterator tau = subset1.begin();
281  tau < subset1.end();
282  ++tau) {
283  if ((*tau)->isPFTau() && (*tau)->isolationPFChargedHadrCandsPtSum() < iso) {
284  *tauIdx = *tau - taus->begin();
285  iso = (*tau)->isolationPFChargedHadrCandsPtSum();
286  }
287  }
288 
289  // check that one combination has been found
290  if (!leptonFound) {
291  leptonFoundEpTm = false;
292  leptonFoundEmTp = false;
293  }
294  // discard the jet that matches the tau (if one)
295  if (leptonFound) {
296  for (std::vector<pat::Jet>::const_iterator jet = jets->begin(); jet < jets->end(); ++jet) {
297  if (DeltaR<pat::Particle, pat::Jet>()(*(taus->begin() + *tauIdx), *jet) < 0.1) {
298  JetVetoByTaus.push_back(jet - jets->begin());
299  }
300  }
301  }
302  }
303  // select Jets (TopJet vector is sorted on ET)
304  jetsFound = ((jets->size() - JetVetoByTaus.size()) >= 2);
305  } else if (taus->size() > 1) {
306  tautau = true;
307  if (LepDiffCharge(&(*taus)[0], &(*taus)[1])) {
308  leptonFound = true;
309  leptonFoundTT = true;
310  if (HasPositiveCharge(&(*taus)[0])) {
311  selTaup = 0;
312  selTaum = 1;
313  } else {
314  selTaup = 1;
315  selTaum = 0;
316  }
317  }
318  for (std::vector<pat::Jet>::const_iterator jet = jets->begin(); jet < jets->end(); ++jet) {
319  if (DeltaR<pat::Particle, pat::Jet>()((*taus)[0], *jet) < 0.1 ||
320  DeltaR<pat::Particle, pat::Jet>()((*taus)[1], *jet) < 0.1) {
321  JetVetoByTaus.push_back(jet - jets->begin());
322  }
323  }
324  // select Jets (TopJet vector is sorted on ET)
325  jetsFound = ((jets->size() - JetVetoByTaus.size()) >= 2);
326  }
327 
328  // Check that the above work makes sense
329  if (int(ee) + int(emu) + int(mumu) + int(etau) + int(mutau) + int(tautau) > 1)
330  std::cout << "[TtDilepEvtSolutionMaker]: "
331  << "Lepton selection criteria uncorrectly defined" << std::endl;
332 
333  bool correctLepton = (leptonFoundEE && eeChannel_) || ((leptonFoundEmMp || leptonFoundEpMm) && emuChannel_) ||
334  (leptonFoundMM && mumuChannel_) || ((leptonFoundMmTp || leptonFoundMpTm) && mutauChannel_) ||
335  ((leptonFoundEmTp || leptonFoundEpTm) && etauChannel_) || (leptonFoundTT && tautauChannel_);
336 
337  std::vector<TtDilepEvtSolution>* evtsols = new std::vector<TtDilepEvtSolution>();
338  if (correctLepton && METFound && jetsFound) {
339  // protect against reading beyond array boundaries while discounting vetoed jets
340  unsigned int nrCombJets = 0;
341  unsigned int numberOfJets = 0;
342  for (; nrCombJets < jets->size() && numberOfJets < nrCombJets_; ++nrCombJets) {
343  if (find(JetVetoByTaus.begin(), JetVetoByTaus.end(), int(nrCombJets)) == JetVetoByTaus.end())
344  ++numberOfJets;
345  }
346  // consider all permutations
347  for (unsigned int ib = 0; ib < nrCombJets; ib++) {
348  // skipped jet vetoed during components-flagging.
349  if (find(JetVetoByTaus.begin(), JetVetoByTaus.end(), int(ib)) != JetVetoByTaus.end())
350  continue;
351  // second loop of the permutations
352  for (unsigned int ibbar = 0; ibbar < nrCombJets; ibbar++) {
353  // avoid the diagonal: b and bbar must be distinct jets
354  if (ib == ibbar)
355  continue;
356  // skipped jet vetoed during components-flagging.
357  if (find(JetVetoByTaus.begin(), JetVetoByTaus.end(), int(ibbar)) != JetVetoByTaus.end())
358  continue;
359  // Build and save a solution
360  TtDilepEvtSolution asol;
362  double xconstraint = 0, yconstraint = 0;
363  // Set e+ in the event
364  if (leptonFoundEE || leptonFoundEpMm || leptonFoundEpTm) {
365  asol.setElectronp(electrons, selElectronp);
366  xconstraint += (*electrons)[selElectronp].px();
367  yconstraint += (*electrons)[selElectronp].py();
368  }
369  // Set e- in the event
370  if (leptonFoundEE || leptonFoundEmMp || leptonFoundEmTp) {
371  asol.setElectronm(electrons, selElectronm);
372  xconstraint += (*electrons)[selElectronm].px();
373  yconstraint += (*electrons)[selElectronm].py();
374  }
375  // Set mu+ in the event
376  if (leptonFoundMM || leptonFoundEmMp || leptonFoundMpTm) {
377  asol.setMuonp(muons, selMuonp);
378  xconstraint += (*muons)[selMuonp].px();
379  yconstraint += (*muons)[selMuonp].py();
380  }
381  // Set mu- in the event
382  if (leptonFoundMM || leptonFoundEpMm || leptonFoundMmTp) {
383  asol.setMuonm(muons, selMuonm);
384  xconstraint += (*muons)[selMuonm].px();
385  yconstraint += (*muons)[selMuonm].py();
386  }
387  // Set tau- in the event
388  if (leptonFoundEpTm || leptonFoundMpTm || leptonFoundTT) {
389  asol.setTaum(taus, selTaum);
390  xconstraint += (*taus)[selTaum].px();
391  yconstraint += (*taus)[selTaum].py();
392  }
393  // Set tau+ in the event
394  if (leptonFoundEmTp || leptonFoundMmTp || leptonFoundTT) {
395  asol.setTaup(taus, selTaup);
396  xconstraint += (*taus)[selTaup].px();
397  yconstraint += (*taus)[selTaup].py();
398  }
399  // Set Jets/MET in the event
400  asol.setB(jets, ib);
401  asol.setBbar(jets, ibbar);
402  asol.setMET(mets, 0);
403  xconstraint += (*jets)[ib].px() + (*jets)[ibbar].px() + (*mets)[0].px();
404  yconstraint += (*jets)[ib].py() + (*jets)[ibbar].py() + (*mets)[0].py();
405  // if asked for, match the event solutions to the gen Event
406  if (matchToGenEvt_) {
408  iEvent.getByToken(evtSourceToken_, genEvt);
409  asol.setGenEvt(genEvt);
410  }
411  // If asked, use the kin fitter to compute the top mass
412  if (calcTopMass_) {
413  solver->SetConstraints(xconstraint, yconstraint);
415  asol = solver->addKinSolInfo(&asol);
416  }
417 
418  // these lines calculate the observables to be used in the TtDilepSignalSelection LR
419  (*myLRSignalSelObservables)(asol, iEvent);
420 
421  evtsols->push_back(asol);
422  }
423  }
424  // flag the best solution (MC matching)
425  if (matchToGenEvt_) {
426  double bestSolDR = 9999.;
427  int bestSol = -1;
428  double dR = 0.;
429  for (size_t s = 0; s < evtsols->size(); s++) {
430  dR = (*evtsols)[s].getJetResidual();
431  if (dR < bestSolDR) {
432  bestSolDR = dR;
433  bestSol = s;
434  }
435  }
436  if (bestSol != -1)
437  (*evtsols)[bestSol].setBestSol(true);
438  }
439  // put the result in the event
440  std::unique_ptr<std::vector<TtDilepEvtSolution> > pOut(evtsols);
441  iEvent.put(std::move(pOut));
442  } else {
443  // no solution: put a dummy solution in the event
444  TtDilepEvtSolution asol;
445  evtsols->push_back(asol);
446  std::unique_ptr<std::vector<TtDilepEvtSolution> > pOut(evtsols);
447  iEvent.put(std::move(pOut));
448  }
449 }
edm::EDGetTokenT< std::vector< pat::Jet > > jetSourceToken_
bool HasPositiveCharge(const reco::Candidate *) const
void setElectronm(const edm::Handle< std::vector< pat::Electron > > &elec, int i)
Definition: deltaR.h:58
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_
void SetConstraints(const double xx=0, const double yy=0)
bool LepDiffCharge(const reco::Candidate *, const reco::Candidate *) const
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
muons
the two sets of parameters below are mutually exclusive, depending if RECO or ALCARECO is used the us...
Definition: DiMuonV_cfg.py:214
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
TtDilepEvtSolution addKinSolInfo(TtDilepEvtSolution *asol)
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)
bool PTComp(const reco::Candidate *, const reco::Candidate *) const
void setTaum(const edm::Handle< std::vector< pat::Tau > > &tau, int i)
def move(src, dest)
Definition: eostools.py:511
ib
Definition: cuy.py:661
void setBbar(const edm::Handle< std::vector< pat::Jet > > &jet, int i)

◆ PTComp()

bool TtDilepEvtSolutionMaker::PTComp ( const reco::Candidate l1,
const reco::Candidate l2 
) const
inlineprivate

Definition at line 49 of file TtDilepEvtSolutionMaker.cc.

References reco::Candidate::pt().

Referenced by produce().

49  {
50  return (l1->pt() > l2->pt());
51 }
virtual double pt() const =0
transverse momentum

Member Data Documentation

◆ calcTopMass_

bool TtDilepEvtSolutionMaker::calcTopMass_
private

Definition at line 40 of file TtDilepEvtSolutionMaker.cc.

Referenced by produce(), and TtDilepEvtSolutionMaker().

◆ eeChannel_

bool TtDilepEvtSolutionMaker::eeChannel_
private

Definition at line 41 of file TtDilepEvtSolutionMaker.cc.

Referenced by produce(), and TtDilepEvtSolutionMaker().

◆ electronSourceToken_

edm::EDGetTokenT<std::vector<pat::Electron> > TtDilepEvtSolutionMaker::electronSourceToken_
private

Definition at line 32 of file TtDilepEvtSolutionMaker.cc.

Referenced by produce(), and TtDilepEvtSolutionMaker().

◆ emuChannel_

bool TtDilepEvtSolutionMaker::emuChannel_
private

Definition at line 41 of file TtDilepEvtSolutionMaker.cc.

Referenced by produce(), and TtDilepEvtSolutionMaker().

◆ etauChannel_

bool TtDilepEvtSolutionMaker::etauChannel_
private

Definition at line 41 of file TtDilepEvtSolutionMaker.cc.

Referenced by produce(), and TtDilepEvtSolutionMaker().

◆ evtSourceToken_

edm::EDGetTokenT<TtGenEvent> TtDilepEvtSolutionMaker::evtSourceToken_
private

Definition at line 37 of file TtDilepEvtSolutionMaker.cc.

Referenced by produce(), and TtDilepEvtSolutionMaker().

◆ jetCorrScheme_

int TtDilepEvtSolutionMaker::jetCorrScheme_
private

Definition at line 38 of file TtDilepEvtSolutionMaker.cc.

Referenced by produce(), and TtDilepEvtSolutionMaker().

◆ jetSourceToken_

edm::EDGetTokenT<std::vector<pat::Jet> > TtDilepEvtSolutionMaker::jetSourceToken_
private

Definition at line 36 of file TtDilepEvtSolutionMaker.cc.

Referenced by produce(), and TtDilepEvtSolutionMaker().

◆ matchToGenEvt_

bool TtDilepEvtSolutionMaker::matchToGenEvt_
private

Definition at line 40 of file TtDilepEvtSolutionMaker.cc.

Referenced by produce(), and TtDilepEvtSolutionMaker().

◆ metSourceToken_

edm::EDGetTokenT<std::vector<pat::MET> > TtDilepEvtSolutionMaker::metSourceToken_
private

Definition at line 35 of file TtDilepEvtSolutionMaker.cc.

Referenced by produce(), and TtDilepEvtSolutionMaker().

◆ mumuChannel_

bool TtDilepEvtSolutionMaker::mumuChannel_
private

Definition at line 41 of file TtDilepEvtSolutionMaker.cc.

Referenced by produce(), and TtDilepEvtSolutionMaker().

◆ muonSourceToken_

edm::EDGetTokenT<std::vector<pat::Muon> > TtDilepEvtSolutionMaker::muonSourceToken_
private

Definition at line 33 of file TtDilepEvtSolutionMaker.cc.

Referenced by produce(), and TtDilepEvtSolutionMaker().

◆ mutauChannel_

bool TtDilepEvtSolutionMaker::mutauChannel_
private

Definition at line 41 of file TtDilepEvtSolutionMaker.cc.

Referenced by produce(), and TtDilepEvtSolutionMaker().

◆ myLRSignalSelObservables

TtDilepLRSignalSelObservables* TtDilepEvtSolutionMaker::myLRSignalSelObservables
private

Definition at line 45 of file TtDilepEvtSolutionMaker.cc.

Referenced by TtDilepEvtSolutionMaker().

◆ nrCombJets_

unsigned int TtDilepEvtSolutionMaker::nrCombJets_
private

Definition at line 39 of file TtDilepEvtSolutionMaker.cc.

Referenced by produce(), and TtDilepEvtSolutionMaker().

◆ nupars_

std::vector<double> TtDilepEvtSolutionMaker::nupars_
private

Definition at line 43 of file TtDilepEvtSolutionMaker.cc.

Referenced by TtDilepEvtSolutionMaker().

◆ solver

TtFullLepKinSolver* TtDilepEvtSolutionMaker::solver
private

Definition at line 46 of file TtDilepEvtSolutionMaker.cc.

Referenced by produce(), and TtDilepEvtSolutionMaker().

◆ tauSourceToken_

edm::EDGetTokenT<std::vector<pat::Tau> > TtDilepEvtSolutionMaker::tauSourceToken_
private

Definition at line 34 of file TtDilepEvtSolutionMaker.cc.

Referenced by produce(), and TtDilepEvtSolutionMaker().

◆ tautauChannel_

bool TtDilepEvtSolutionMaker::tautauChannel_
private

Definition at line 41 of file TtDilepEvtSolutionMaker.cc.

Referenced by produce(), and TtDilepEvtSolutionMaker().

◆ tmassbegin_

double TtDilepEvtSolutionMaker::tmassbegin_
private

Definition at line 42 of file TtDilepEvtSolutionMaker.cc.

Referenced by TtDilepEvtSolutionMaker().

◆ tmassend_

double TtDilepEvtSolutionMaker::tmassend_
private

Definition at line 42 of file TtDilepEvtSolutionMaker.cc.

Referenced by TtDilepEvtSolutionMaker().

◆ tmassstep_

double TtDilepEvtSolutionMaker::tmassstep_
private

Definition at line 42 of file TtDilepEvtSolutionMaker.cc.

Referenced by TtDilepEvtSolutionMaker().

◆ useMCforBest_

bool TtDilepEvtSolutionMaker::useMCforBest_
private

Definition at line 40 of file TtDilepEvtSolutionMaker.cc.

Referenced by produce(), and TtDilepEvtSolutionMaker().