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 
8 { }
9 
11 
12 void
15  const edm::Handle<std::vector<pat::MET> >& mets,
16  const edm::Handle<std::vector<pat::Jet> >& jets,
17  std::vector<int>& match, const unsigned int iComb)
18 {
19  // -----------------------------------------------------
20  // get genEvent (to distinguish between uds and c quarks
21  // and for the lepton matching)
22  // -----------------------------------------------------
24  evt.getByLabel("genEvt", genEvt);
25 
26  // -----------------------------------------------------
27  // add jets
28  // -----------------------------------------------------
29  for(unsigned idx=0; idx<match.size(); ++idx){
30  if( isValid(match[idx], jets) ){
31  switch(idx){
33  if( std::abs(genEvt->hadronicDecayQuark()->pdgId())==4 )
34  setCandidate(jets, match[idx], lightQ_, jetCorrectionLevel("cQuark"));
35  else
36  setCandidate(jets, match[idx], lightQ_, jetCorrectionLevel("udsQuark"));
37  break;
39  if( std::abs(genEvt->hadronicDecayQuarkBar()->pdgId())==4 )
40  setCandidate(jets, match[idx], lightQBar_, jetCorrectionLevel("cQuark"));
41  else
42  setCandidate(jets, match[idx], lightQBar_, jetCorrectionLevel("udsQuark"));
43  break;
45  setCandidate(jets, match[idx], hadronicB_, jetCorrectionLevel("bQuark")); break;
47  setCandidate(jets, match[idx], leptonicB_, jetCorrectionLevel("bQuark")); break;
48  }
49  }
50  }
51 
52  // -----------------------------------------------------
53  // add lepton
54  // -----------------------------------------------------
55  int iLepton = findMatchingLepton(genEvt, leps);
56  if( iLepton<0 )
57  return;
58  setCandidate(leps, iLepton, lepton_);
59  match.push_back( iLepton );
60 
61  // -----------------------------------------------------
62  // add neutrino
63  // -----------------------------------------------------
64  if( mets->empty() )
65  return;
66  if(neutrinoSolutionType_ == -1)
67  setCandidate(mets, 0, neutrino_);
68  else
69  setNeutrino(mets, leps, iLepton, neutrinoSolutionType_);
70 }
71 
73 int
76 {
77  int genIdx=-1;
78 
79  // jump out with -1 when the collection is empty
80  if( leps->empty() ) return genIdx;
81 
82  if( genEvt->isTtBar() && genEvt->isSemiLeptonic( leptonType( &(leps->front()) ) ) && genEvt->singleLepton() ){
83  double minDR=-1;
84  for(unsigned i=0; i<leps->size(); ++i){
85  double dR = deltaR(genEvt->singleLepton()->eta(), genEvt->singleLepton()->phi(), (*leps)[i].eta(), (*leps)[i].phi());
86  if(minDR<0 || dR<minDR){
87  minDR=dR;
88  genIdx=i;
89  }
90  }
91  }
92  return genIdx;
93 }
int i
Definition: DBlmapReader.cc:9
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
reco::ShallowClonePtrCandidate * lepton_
reco::ShallowClonePtrCandidate * lightQBar_
#define abs(x)
Definition: mlp_lapack.h:159
reco::ShallowClonePtrCandidate * neutrino_
vector< PseudoJet > jets
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...
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
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
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:6
int neutrinoSolutionType_
algorithm used to calculate neutrino solutions (see cfi for further details)
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 ...