CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TtSemiLepHypGenMatch.cc
Go to the documentation of this file.
5 
7  TtSemiLepHypothesis( cfg ),
8  genEvtToken_( consumes<TtGenEvent>( edm::InputTag( "genEvt" ) ) )
9 { }
10 
12 
13 void
16  const edm::Handle<std::vector<pat::MET> >& mets,
17  const edm::Handle<std::vector<pat::Jet> >& jets,
18  std::vector<int>& match, const unsigned int iComb)
19 {
20  // -----------------------------------------------------
21  // get genEvent (to distinguish between uds and c quarks
22  // and for the lepton matching)
23  // -----------------------------------------------------
25  evt.getByToken(genEvtToken_, genEvt);
26 
27  // -----------------------------------------------------
28  // add jets
29  // -----------------------------------------------------
30  for(unsigned idx=0; idx<match.size(); ++idx){
31  if( isValid(match[idx], jets) ){
32  switch(idx){
34  if( std::abs(genEvt->hadronicDecayQuark()->pdgId())==4 )
35  setCandidate(jets, match[idx], lightQ_, jetCorrectionLevel("cQuark"));
36  else
37  setCandidate(jets, match[idx], lightQ_, jetCorrectionLevel("udsQuark"));
38  break;
40  if( std::abs(genEvt->hadronicDecayQuarkBar()->pdgId())==4 )
41  setCandidate(jets, match[idx], lightQBar_, jetCorrectionLevel("cQuark"));
42  else
43  setCandidate(jets, match[idx], lightQBar_, jetCorrectionLevel("udsQuark"));
44  break;
46  setCandidate(jets, match[idx], hadronicB_, jetCorrectionLevel("bQuark")); break;
48  setCandidate(jets, match[idx], leptonicB_, jetCorrectionLevel("bQuark")); break;
49  }
50  }
51  }
52 
53  // -----------------------------------------------------
54  // add lepton
55  // -----------------------------------------------------
56  int iLepton = findMatchingLepton(genEvt, leps);
57  if( iLepton<0 )
58  return;
59  setCandidate(leps, iLepton, lepton_);
60  match.push_back( iLepton );
61 
62  // -----------------------------------------------------
63  // add neutrino
64  // -----------------------------------------------------
65  if( mets->empty() )
66  return;
67  if(neutrinoSolutionType_ == -1)
68  setCandidate(mets, 0, neutrino_);
69  else
70  setNeutrino(mets, leps, iLepton, neutrinoSolutionType_);
71 }
72 
74 int
77 {
78  int genIdx=-1;
79 
80  // jump out with -1 when the collection is empty
81  if( leps->empty() ) return genIdx;
82 
83  if( genEvt->isTtBar() && genEvt->isSemiLeptonic( leptonType( &(leps->front()) ) ) && genEvt->singleLepton() ){
84  double minDR=-1;
85  for(unsigned i=0; i<leps->size(); ++i){
86  double dR = deltaR(genEvt->singleLepton()->eta(), genEvt->singleLepton()->phi(), (*leps)[i].eta(), (*leps)[i].phi());
87  if(minDR<0 || dR<minDR){
88  minDR=dR;
89  genIdx=i;
90  }
91  }
92  }
93  return genIdx;
94 }
int i
Definition: DBlmapReader.cc:9
tuple cfg
Definition: looper.py:259
WDecay::LeptonType leptonType(const reco::RecoCandidate *cand)
determine lepton type of reco candidate and return a corresponding WDecay::LeptonType; the type is kN...
bool isValid(const int &idx, const edm::Handle< std::vector< pat::Jet > > &jets)
check if index is in valid range of selected jets
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
reco::ShallowClonePtrCandidate * lepton_
reco::ShallowClonePtrCandidate * lightQBar_
Class derived from the TopGenEvent for ttbar events.
Definition: TtGenEvent.h:18
reco::ShallowClonePtrCandidate * neutrino_
vector< PseudoJet > jets
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int findMatchingLepton(const edm::Handle< TtGenEvent > &genEvt, const edm::Handle< edm::View< reco::RecoCandidate > > &)
find index of the candidate nearest to the singleLepton of the generator event in the collection; ret...
reco::ShallowClonePtrCandidate * hadronicB_
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
virtual void buildHypo(edm::Event &, const edm::Handle< edm::View< reco::RecoCandidate > > &, const edm::Handle< std::vector< pat::MET > > &, const edm::Handle< std::vector< pat::Jet > > &, std::vector< int > &, const unsigned int iComb)
build event hypothesis from the reco objects of a semi-leptonic event
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
TtSemiLepHypGenMatch(const edm::ParameterSet &)
void setNeutrino(const edm::Handle< std::vector< pat::MET > > &met, const edm::Handle< edm::View< reco::RecoCandidate > > &leps, const int &idx, const int &type)
set neutrino, using mW = 80.4 to calculate the neutrino pz
void setCandidate(const edm::Handle< C > &handle, const int &idx, reco::ShallowClonePtrCandidate *&clone)
use one object in a collection to set a ShallowClonePtrCandidate
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
int neutrinoSolutionType_
algorithm used to calculate neutrino solutions (see cfi for further details)
edm::EDGetTokenT< TtGenEvent > genEvtToken_
genEvtToken_(mayConsume< TtGenEvent >(genEvt_))
reco::ShallowClonePtrCandidate * lightQ_
reco::ShallowClonePtrCandidate * leptonicB_
std::string jetCorrectionLevel(const std::string &quarkType)
helper function to construct the proper correction level string for corresponding quarkType ...