CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
TtDilepEvtSolutionMaker Class Reference

#include <TtDilepEvtSolutionMaker.h>

Inheritance diagram for TtDilepEvtSolutionMaker:
edm::EDProducer edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

virtual void beginJob ()
 
virtual void produce (edm::Event &iEvent, const edm::EventSetup &iSetup)
 
 TtDilepEvtSolutionMaker (const edm::ParameterSet &iConfig)
 
 ~TtDilepEvtSolutionMaker ()
 
- Public Member Functions inherited from edm::EDProducer
 EDProducer ()
 
virtual ~EDProducer ()
 
- Public Member Functions inherited from edm::ProducerBase
 ProducerBase ()
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription
const &)> 
registrationCallback () const
 used by the fwk to register list of products More...
 
virtual ~ProducerBase ()
 
- Public Member Functions inherited from edm::EDConsumerBase
 EDConsumerBase ()
 
ProductHolderIndex indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndex > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndex > &) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

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::InputTag electronSource_
 
bool emuChannel_
 
bool etauChannel_
 
edm::InputTag evtSource_
 
int jetCorrScheme_
 
edm::InputTag jetSource_
 
bool matchToGenEvt_
 
edm::InputTag metSource_
 
bool mumuChannel_
 
edm::InputTag muonSource_
 
bool mutauChannel_
 
TtDilepLRSignalSelObservablesmyLRSignalSelObservables
 
unsigned int nrCombJets_
 
std::vector< double > nupars_
 
TtFullLepKinSolversolver
 
edm::InputTag tauSource_
 
bool tautauChannel_
 
double tmassbegin_
 
double tmassend_
 
double tmassstep_
 
bool useMCforBest_
 

Additional Inherited Members

- Public Types inherited from edm::EDProducer
typedef EDProducer ModuleType
 
typedef WorkerT< EDProducerWorkerType
 
- Public Types inherited from edm::ProducerBase
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 
- Static Public Member Functions inherited from edm::EDProducer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::EDProducer
CurrentProcessingContext const * currentContext () const
 
- Protected Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Definition at line 14 of file TtDilepEvtSolutionMaker.h.

Constructor & Destructor Documentation

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

Definition at line 17 of file TtDilepEvtSolutionMaker.cc.

References calcTopMass_, eeChannel_, electronSource_, emuChannel_, etauChannel_, evtSource_, edm::ParameterSet::getParameter(), jetCorrScheme_, jetSource_, matchToGenEvt_, metSource_, mumuChannel_, muonSource_, mutauChannel_, myLRSignalSelObservables, nrCombJets_, nupars_, tauSource_, tautauChannel_, tmassbegin_, tmassend_, tmassstep_, and useMCforBest_.

18 {
19  // configurables
20  electronSource_ = iConfig.getParameter<edm::InputTag>("electronSource");
21  muonSource_ = iConfig.getParameter<edm::InputTag>("muonSource");
22  tauSource_ = iConfig.getParameter<edm::InputTag>("tauSource");
23  metSource_ = iConfig.getParameter<edm::InputTag>("metSource");
24  jetSource_ = iConfig.getParameter<edm::InputTag>("jetSource");
25  jetCorrScheme_ = iConfig.getParameter<int> ("jetCorrectionScheme");
26  evtSource_ = iConfig.getParameter<edm::InputTag>("evtSource");
27  nrCombJets_ = iConfig.getParameter<unsigned int> ("nrCombJets");
28  matchToGenEvt_ = iConfig.getParameter<bool> ("matchToGenEvt");
29  calcTopMass_ = iConfig.getParameter<bool> ("calcTopMass");
30  useMCforBest_ = iConfig.getParameter<bool> ("bestSolFromMC");
31  eeChannel_ = iConfig.getParameter<bool> ("eeChannel");
32  emuChannel_ = iConfig.getParameter<bool> ("emuChannel");
33  mumuChannel_ = iConfig.getParameter<bool> ("mumuChannel");
34  mutauChannel_ = iConfig.getParameter<bool> ("mutauChannel");
35  etauChannel_ = iConfig.getParameter<bool> ("etauChannel");
36  tautauChannel_ = iConfig.getParameter<bool> ("tautauChannel");
37  tmassbegin_ = iConfig.getParameter<double> ("tmassbegin");
38  tmassend_ = iConfig.getParameter<double> ("tmassend");
39  tmassstep_ = iConfig.getParameter<double> ("tmassstep");
40  nupars_ = iConfig.getParameter<std::vector<double> >("neutrino_parameters");
41 
42  // define what will be produced
43  produces<std::vector<TtDilepEvtSolution> >();
44 
46  myLRSignalSelObservables->jetSource(jetSource_);
47 }
T getParameter(std::string const &) const
std::vector< double > nupars_
TtDilepLRSignalSelObservables * myLRSignalSelObservables
TtDilepEvtSolutionMaker::~TtDilepEvtSolutionMaker ( )

Definition at line 49 of file TtDilepEvtSolutionMaker.cc.

50 {
51 }

Member Function Documentation

void TtDilepEvtSolutionMaker::beginJob ( void  )
virtual
bool TtDilepEvtSolutionMaker::HasPositiveCharge ( const reco::Candidate l) const
inlineprivate

Definition at line 60 of file TtDilepEvtSolutionMaker.h.

References reco::Candidate::charge().

Referenced by produce().

61 {
62  return (l->charge() > 0);
63 }
virtual int charge() const =0
electric charge
bool TtDilepEvtSolutionMaker::LepDiffCharge ( const reco::Candidate l1,
const reco::Candidate l2 
) const
inlineprivate

Definition at line 55 of file TtDilepEvtSolutionMaker.h.

References reco::Candidate::charge().

Referenced by produce().

56 {
57  return (l1->charge() != l2->charge());
58 }
virtual int charge() const =0
electric charge
void TtDilepEvtSolutionMaker::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
virtual

Implements edm::EDProducer.

Definition at line 58 of file TtDilepEvtSolutionMaker.cc.

References TtFullLepKinSolver::addKinSolInfo(), calcTopMass_, gather_cfg::cout, PFRecoTauDiscriminationAgainstElectronDeadECAL_cfi::dR, eeChannel_, HI_PhotonSkim_cff::electrons, electronSource_, emuChannel_, etauChannel_, evtSource_, spr::find(), TtGenEvtProducer_cfi::genEvt, edm::Event::getByLabel(), HasPositiveCharge(), cuy::ib, iEvent, metsig::jet, jetCorrScheme_, fwrapper::jets, jetSource_, LepDiffCharge(), matchToGenEvt_, metSource_, mumuChannel_, patZpeak::muons, muonSource_, mutauChannel_, nrCombJets_, NULL, PTComp(), edm::Event::put(), 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, tauSource_, tautauChannel_, useMCforBest_, and TtFullLepKinSolver::useWeightFromMC().

59 {
61  iEvent.getByLabel(tauSource_, taus);
63  iEvent.getByLabel(muonSource_, muons);
65  iEvent.getByLabel(electronSource_, electrons);
67  iEvent.getByLabel(metSource_, mets);
69  iEvent.getByLabel(jetSource_, jets);
70 
71  int selMuonp = -1, selMuonm = -1;
72  int selElectronp = -1, selElectronm = -1;
73  int selTaup = -1, selTaum = -1;
74  bool leptonFound = false;
75  bool mumu = false;
76  bool emu = false;
77  bool ee = false;
78  bool etau = false;
79  bool mutau = false;
80  bool tautau = false;
81  bool leptonFoundEE = false;
82  bool leptonFoundMM = false;
83  bool leptonFoundTT = false;
84  bool leptonFoundEpMm = false;
85  bool leptonFoundEmMp = false;
86  bool leptonFoundEpTm = false;
87  bool leptonFoundEmTp = false;
88  bool leptonFoundMpTm = false;
89  bool leptonFoundMmTp = false;
90  bool jetsFound = false;
91  bool METFound = false;
92  std::vector<int> JetVetoByTaus;
93 
94  //select MET (TopMET vector is sorted on ET)
95  if(mets->size()>=1) { METFound = true; }
96 
97  // If we have electrons and muons available,
98  // build a solutions with electrons and muons.
99  if (muons->size() + electrons->size() >=2) {
100  // select leptons
101  if (electrons->size() == 0) mumu = true;
102  else if (muons->size() == 0) ee = true;
103  else if (electrons->size() == 1) {
104  if (muons->size() == 1) emu = true;
105  else if (PTComp(&(*electrons)[0], &(*muons)[1])) emu = true;
106  else mumu = true;
107  }
108  else if (electrons->size() > 1) {
109  if (PTComp(&(*electrons)[1], &(*muons)[0])) ee = true;
110  else if (muons->size() == 1) emu = true;
111  else if (PTComp(&(*electrons)[0], &(*muons)[1])) emu = true;
112  else mumu = true;
113  }
114  if (ee) {
115  if (LepDiffCharge(&(*electrons)[0], &(*electrons)[1])) {
116  leptonFound = true;
117  leptonFoundEE = true;
118  if (HasPositiveCharge(&(*electrons)[0])) {
119  selElectronp = 0;
120  selElectronm = 1;
121  } else {
122  selElectronp = 1;
123  selElectronm = 0;
124  }
125  }
126  }
127  else if (emu) {
128  if (LepDiffCharge(&(*electrons)[0], &(*muons)[0])) {
129  leptonFound = true;
130  if (HasPositiveCharge(&(*electrons)[0])) {
131  leptonFoundEpMm = true;
132  selElectronp = 0;
133  selMuonm = 0;
134  } else {
135  leptonFoundEmMp = true;
136  selMuonp = 0;
137  selElectronm = 0;
138  }
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) { jetsFound = true; }
156  }
157  // If a tau is needed to have two leptons, then only consider the taus.
158  // This is the minimal modification of the dilept selection that includes taus,
159  // since we are considering taus only when no other solution exist.
160  else if(muons->size() + electrons->size()==1 && taus->size()>0) {
161  // select leptons
162  if(muons->size()==1) {
163  mutau = true;
164  // depending on the muon charge, set the right muon index and specify channel
165  int expectedCharge = - muons->begin()->charge();
166  int* tauIdx = NULL;
167  if (expectedCharge<0) {
168  selMuonp = 0;
169  tauIdx = &selTaum;
170  leptonFoundMpTm = true;
171  } else {
172  selMuonm = 0;
173  tauIdx = &selTaup;
174  leptonFoundMmTp = true;
175  }
176  // loop over the vector of taus to find the ones
177  // that have the charge opposite to the muon one, and do not match in eta-phi
178  std::vector<std::vector<pat::Tau>::const_iterator> subset1;
179  for(std::vector<pat::Tau>::const_iterator tau = taus->begin(); tau < taus->end(); ++tau ) {
180  if(tau->charge()*expectedCharge>=0 && DeltaR<pat::Particle>()(*tau,*(muons->begin()))>0.1) {
181  *tauIdx = tau-taus->begin();
182  leptonFound = true;
183  subset1.push_back(tau);
184  }
185  }
186  // if there are more than one tau with ecalIsol==0, take the smallest E/P
187  float iso = 999.;
188  for(std::vector<std::vector<pat::Tau>::const_iterator>::const_iterator tau = subset1.begin(); tau < subset1.end(); ++tau) {
189  if((*tau)->isCaloTau() && (*tau)->isolationTracksPtSum()<iso) {
190  *tauIdx = *tau - taus->begin();
191  iso = (*tau)->isolationTracksPtSum();
192  }
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) { leptonFoundMpTm = false; leptonFoundMmTp = false; }
201  // discard the jet that matches the tau (if one)
202  if(leptonFound) {
203  for(std::vector<pat::Jet>::const_iterator jet = jets->begin(); jet<jets->end(); ++jet) {
204  if(DeltaR<pat::Particle, pat::Jet>()(*(taus->begin()+*tauIdx),*jet)<0.1) {
205  JetVetoByTaus.push_back(jet-jets->begin());
206  }
207  }
208  }
209  }
210  else {
211  etau = true;
212  // depending on the electron charge, set the right electron index and specify channel
213  int expectedCharge = - electrons->begin()->charge();
214  int* tauIdx = NULL;
215  if (expectedCharge<0) {
216  selElectronp = 0;
217  tauIdx = &selTaum;
218  leptonFoundEpTm = true;
219  } else {
220  selElectronm = 0;
221  tauIdx = &selTaup;
222  leptonFoundEmTp = true;
223  }
224  // loop over the vector of taus to find the ones
225  // that have the charge opposite to the muon one, and do not match in eta-phi
226  std::vector<std::vector<pat::Tau>::const_iterator> subset1;
227  for(std::vector<pat::Tau>::const_iterator tau = taus->begin(); tau < taus->end(); ++tau ) {
228  if(tau->charge()*expectedCharge>=0 && DeltaR<pat::Particle>()(*tau,*(electrons->begin()))>0.1) {
229  *tauIdx = tau-taus->begin();
230  leptonFound = true;
231  subset1.push_back(tau);
232  }
233  }
234  // if there are more than one tau with ecalIsol==0, take the smallest E/P
235  float iso = 999.;
236  for(std::vector<std::vector<pat::Tau>::const_iterator>::const_iterator tau = subset1.begin(); tau < subset1.end(); ++tau) {
237  if((*tau)->isCaloTau() && (*tau)->isolationTracksPtSum()<iso) {
238  *tauIdx = *tau - taus->begin();
239  iso = (*tau)->isolationTracksPtSum();
240  }
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) { leptonFoundEpTm = false; leptonFoundEmTp = false; }
249  // discard the jet that matches the tau (if one)
250  if(leptonFound) {
251  for(std::vector<pat::Jet>::const_iterator jet = jets->begin(); jet<jets->end(); ++jet) {
252  if(DeltaR<pat::Particle, pat::Jet>()(*(taus->begin()+*tauIdx),*jet)<0.1) {
253  JetVetoByTaus.push_back(jet-jets->begin());
254  }
255  }
256  }
257  }
258  // select Jets (TopJet vector is sorted on ET)
259  jetsFound = ((jets->size()-JetVetoByTaus.size())>=2);
260  } else if(taus->size()>1) {
261  tautau = true;
262  if(LepDiffCharge(&(*taus)[0],&(*taus)[1])) {
263  leptonFound = true;
264  leptonFoundTT = true;
265  if(HasPositiveCharge(&(*taus)[0])) {
266  selTaup = 0;
267  selTaum = 1;
268  }
269  else {
270  selTaup = 1;
271  selTaum = 0;
272  }
273  }
274  for(std::vector<pat::Jet>::const_iterator jet = jets->begin(); jet<jets->end(); ++jet) {
275  if(DeltaR<pat::Particle, pat::Jet>()((*taus)[0],*jet)<0.1 || DeltaR<pat::Particle, pat::Jet>()((*taus)[1],*jet)<0.1) {
276  JetVetoByTaus.push_back(jet-jets->begin());
277  }
278  }
279  // select Jets (TopJet vector is sorted on ET)
280  jetsFound = ((jets->size()-JetVetoByTaus.size())>=2);
281  }
282 
283  // Check that the above work makes sense
284  if(int(ee)+int(emu)+int(mumu)+int(etau)+int(mutau)+int(tautau)>1)
285  std::cout << "[TtDilepEvtSolutionMaker]: "
286  << "Lepton selection criteria uncorrectly defined" << std::endl;
287 
288  bool correctLepton = (leptonFoundEE && eeChannel_) ||
289  ((leptonFoundEmMp || leptonFoundEpMm) && emuChannel_) ||
290  (leptonFoundMM && mumuChannel_) ||
291  ((leptonFoundMmTp || leptonFoundMpTm) && mutauChannel_)||
292  ((leptonFoundEmTp || leptonFoundEpTm) && etauChannel_) ||
293  (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()) ++numberOfJets;
302  }
303  // consider all permutations
304  for (unsigned int ib = 0; ib < nrCombJets; ib++) {
305  // skipped jet vetoed during components-flagging.
306  if(find(JetVetoByTaus.begin(),JetVetoByTaus.end(),int(ib))!=JetVetoByTaus.end())continue;
307  // second loop of the permutations
308  for (unsigned int ibbar = 0; ibbar < nrCombJets; ibbar++) {
309  // avoid the diagonal: b and bbar must be distinct jets
310  if(ib==ibbar) continue;
311  // skipped jet vetoed during components-flagging.
312  if(find(JetVetoByTaus.begin(),JetVetoByTaus.end(),int(ibbar))!=JetVetoByTaus.end())continue;
313  // Build and save a solution
314  TtDilepEvtSolution asol;
316  double xconstraint = 0, yconstraint = 0;
317  // Set e+ in the event
318  if (leptonFoundEE || leptonFoundEpMm || leptonFoundEpTm) {
319  asol.setElectronp(electrons, selElectronp);
320  xconstraint += (*electrons)[selElectronp].px();
321  yconstraint += (*electrons)[selElectronp].py();
322  }
323  // Set e- in the event
324  if (leptonFoundEE || leptonFoundEmMp || leptonFoundEmTp) {
325  asol.setElectronm(electrons, selElectronm);
326  xconstraint += (*electrons)[selElectronm].px();
327  yconstraint += (*electrons)[selElectronm].py();
328  }
329  // Set mu+ in the event
330  if (leptonFoundMM || leptonFoundEmMp || leptonFoundMpTm) {
331  asol.setMuonp(muons, selMuonp);
332  xconstraint += (*muons)[selMuonp].px();
333  yconstraint += (*muons)[selMuonp].py();
334  }
335  // Set mu- in the event
336  if (leptonFoundMM || leptonFoundEpMm || leptonFoundMmTp) {
337  asol.setMuonm(muons, selMuonm);
338  xconstraint += (*muons)[selMuonm].px();
339  yconstraint += (*muons)[selMuonm].py();
340  }
341  // Set tau- in the event
342  if (leptonFoundEpTm || leptonFoundMpTm || leptonFoundTT) {
343  asol.setTaum(taus, selTaum);
344  xconstraint += (*taus)[selTaum].px();
345  yconstraint += (*taus)[selTaum].py();
346  }
347  // Set tau+ in the event
348  if (leptonFoundEmTp || leptonFoundMmTp || leptonFoundTT) {
349  asol.setTaup(taus, selTaup);
350  xconstraint += (*taus)[selTaup].px();
351  yconstraint += (*taus)[selTaup].py();
352  }
353  // Set Jets/MET in the event
354  asol.setB(jets, ib);
355  asol.setBbar(jets, ibbar);
356  asol.setMET(mets, 0);
357  xconstraint += (*jets)[ib].px() + (*jets)[ibbar].px() + (*mets)[0].px();
358  yconstraint += (*jets)[ib].py() + (*jets)[ibbar].py() + (*mets)[0].py();
359  // if asked for, match the event solutions to the gen Event
360  if(matchToGenEvt_){
362  iEvent.getByLabel (evtSource_,genEvt);
363  asol.setGenEvt(genEvt);
364  }
365  // If asked, use the kin fitter to compute the top mass
366  if (calcTopMass_) {
367  solver->SetConstraints(xconstraint, yconstraint);
369  asol = solver->addKinSolInfo(&asol);
370  }
371 
372  // these lines calculate the observables to be used in the TtDilepSignalSelection LR
373  (*myLRSignalSelObservables)(asol, iEvent);
374 
375  evtsols->push_back(asol);
376  }
377  }
378  // flag the best solution (MC matching)
379  if(matchToGenEvt_){
380  double bestSolDR = 9999.;
381  int bestSol = -1;
382  double dR = 0.;
383  for(size_t s=0; s<evtsols->size(); s++) {
384  dR = (*evtsols)[s].getJetResidual();
385  if(dR<bestSolDR) { bestSolDR = dR; bestSol = s; }
386  }
387  if(bestSol!=-1) (*evtsols)[bestSol].setBestSol(true);
388  }
389  // put the result in the event
390  std::auto_ptr<std::vector<TtDilepEvtSolution> > pOut(evtsols);
391  iEvent.put(pOut);
392  } else {
393  // no solution: put a dummy solution in the event
394  TtDilepEvtSolution asol;
395  evtsols->push_back(asol);
396  std::auto_ptr<std::vector<TtDilepEvtSolution> > pOut(evtsols);
397  iEvent.put(pOut);
398  }
399 }
int ib
Definition: cuy.py:660
void setElectronm(const edm::Handle< std::vector< pat::Electron > > &elec, int i)
Definition: deltaR.h:79
void setElectronp(const edm::Handle< std::vector< pat::Electron > > &elec, int i)
void setB(const edm::Handle< std::vector< pat::Jet > > &jet, int i)
bool HasPositiveCharge(const reco::Candidate *) const
void SetConstraints(const double xx=0, const double yy=0)
#define NULL
Definition: scimark2.h:8
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
void setMuonm(const edm::Handle< std::vector< pat::Muon > > &muon, int i)
int iEvent
Definition: GenABIO.cc:243
bool PTComp(const reco::Candidate *, const reco::Candidate *) const
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:94
TtDilepEvtSolution addKinSolInfo(TtDilepEvtSolution *asol)
vector< PseudoJet > jets
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)
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
void setJetCorrectionScheme(int jetCorrScheme)
void setGenEvt(const edm::Handle< TtGenEvent > &)
void setMuonp(const edm::Handle< std::vector< pat::Muon > > &muon, int i)
bool LepDiffCharge(const reco::Candidate *, const reco::Candidate *) const
tuple muons
Definition: patZpeak.py:38
void setTaum(const edm::Handle< std::vector< pat::Tau > > &tau, int i)
tuple cout
Definition: gather_cfg.py:121
void setBbar(const edm::Handle< std::vector< pat::Jet > > &jet, int i)
bool TtDilepEvtSolutionMaker::PTComp ( const reco::Candidate l1,
const reco::Candidate l2 
) const
inlineprivate

Definition at line 50 of file TtDilepEvtSolutionMaker.h.

References reco::Candidate::pt().

Referenced by produce().

51 {
52  return (l1->pt() > l2->pt());
53 }
virtual float pt() const =0
transverse momentum

Member Data Documentation

bool TtDilepEvtSolutionMaker::calcTopMass_
private

Definition at line 41 of file TtDilepEvtSolutionMaker.h.

Referenced by produce(), and TtDilepEvtSolutionMaker().

bool TtDilepEvtSolutionMaker::eeChannel_
private

Definition at line 42 of file TtDilepEvtSolutionMaker.h.

Referenced by produce(), and TtDilepEvtSolutionMaker().

edm::InputTag TtDilepEvtSolutionMaker::electronSource_
private

Definition at line 33 of file TtDilepEvtSolutionMaker.h.

Referenced by produce(), and TtDilepEvtSolutionMaker().

bool TtDilepEvtSolutionMaker::emuChannel_
private

Definition at line 42 of file TtDilepEvtSolutionMaker.h.

Referenced by produce(), and TtDilepEvtSolutionMaker().

bool TtDilepEvtSolutionMaker::etauChannel_
private

Definition at line 42 of file TtDilepEvtSolutionMaker.h.

Referenced by produce(), and TtDilepEvtSolutionMaker().

edm::InputTag TtDilepEvtSolutionMaker::evtSource_
private

Definition at line 38 of file TtDilepEvtSolutionMaker.h.

Referenced by produce(), and TtDilepEvtSolutionMaker().

int TtDilepEvtSolutionMaker::jetCorrScheme_
private

Definition at line 39 of file TtDilepEvtSolutionMaker.h.

Referenced by produce(), and TtDilepEvtSolutionMaker().

edm::InputTag TtDilepEvtSolutionMaker::jetSource_
private

Definition at line 37 of file TtDilepEvtSolutionMaker.h.

Referenced by produce(), and TtDilepEvtSolutionMaker().

bool TtDilepEvtSolutionMaker::matchToGenEvt_
private

Definition at line 41 of file TtDilepEvtSolutionMaker.h.

Referenced by produce(), and TtDilepEvtSolutionMaker().

edm::InputTag TtDilepEvtSolutionMaker::metSource_
private

Definition at line 36 of file TtDilepEvtSolutionMaker.h.

Referenced by produce(), and TtDilepEvtSolutionMaker().

bool TtDilepEvtSolutionMaker::mumuChannel_
private

Definition at line 42 of file TtDilepEvtSolutionMaker.h.

Referenced by produce(), and TtDilepEvtSolutionMaker().

edm::InputTag TtDilepEvtSolutionMaker::muonSource_
private

Definition at line 34 of file TtDilepEvtSolutionMaker.h.

Referenced by produce(), and TtDilepEvtSolutionMaker().

bool TtDilepEvtSolutionMaker::mutauChannel_
private

Definition at line 42 of file TtDilepEvtSolutionMaker.h.

Referenced by produce(), and TtDilepEvtSolutionMaker().

TtDilepLRSignalSelObservables* TtDilepEvtSolutionMaker::myLRSignalSelObservables
private

Definition at line 46 of file TtDilepEvtSolutionMaker.h.

Referenced by TtDilepEvtSolutionMaker().

unsigned int TtDilepEvtSolutionMaker::nrCombJets_
private

Definition at line 40 of file TtDilepEvtSolutionMaker.h.

Referenced by produce(), and TtDilepEvtSolutionMaker().

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

Definition at line 44 of file TtDilepEvtSolutionMaker.h.

Referenced by beginJob(), and TtDilepEvtSolutionMaker().

TtFullLepKinSolver* TtDilepEvtSolutionMaker::solver
private

Definition at line 47 of file TtDilepEvtSolutionMaker.h.

Referenced by beginJob(), and produce().

edm::InputTag TtDilepEvtSolutionMaker::tauSource_
private

Definition at line 35 of file TtDilepEvtSolutionMaker.h.

Referenced by produce(), and TtDilepEvtSolutionMaker().

bool TtDilepEvtSolutionMaker::tautauChannel_
private

Definition at line 42 of file TtDilepEvtSolutionMaker.h.

Referenced by produce(), and TtDilepEvtSolutionMaker().

double TtDilepEvtSolutionMaker::tmassbegin_
private

Definition at line 43 of file TtDilepEvtSolutionMaker.h.

Referenced by beginJob(), and TtDilepEvtSolutionMaker().

double TtDilepEvtSolutionMaker::tmassend_
private

Definition at line 43 of file TtDilepEvtSolutionMaker.h.

Referenced by beginJob(), and TtDilepEvtSolutionMaker().

double TtDilepEvtSolutionMaker::tmassstep_
private

Definition at line 43 of file TtDilepEvtSolutionMaker.h.

Referenced by beginJob(), and TtDilepEvtSolutionMaker().

bool TtDilepEvtSolutionMaker::useMCforBest_
private

Definition at line 41 of file TtDilepEvtSolutionMaker.h.

Referenced by produce(), and TtDilepEvtSolutionMaker().