CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TtFullLepHypGenMatch.cc
Go to the documentation of this file.
4 
6  TtFullLepHypothesis( cfg )
7 {
8 }
9 
11 
12 void
14  const edm::Handle<std::vector<pat::Electron > >& elecs,
15  const edm::Handle<std::vector<pat::Muon> >& mus,
16  const edm::Handle<std::vector<pat::Jet> >& jets,
17  const edm::Handle<std::vector<pat::MET> >& mets,
18  std::vector<int>& match,
19  const unsigned int iComb)
20 {
21  // -----------------------------------------------------
22  // add jets
23  // -----------------------------------------------------
24  for(unsigned idx=0; idx<match.size(); ++idx){
25  if( isValid(match[idx], jets) ){
26  switch(idx){
28  setCandidate(jets, match[idx], b_ , jetCorrectionLevel_); break;
30  setCandidate(jets, match[idx], bBar_, jetCorrectionLevel_); break;
31  }
32  }
33  }
34 
35  // -----------------------------------------------------
36  // add leptons
37  // -----------------------------------------------------
38  // get genEvent
40  evt.getByLabel("genEvt", genEvt);
41 
42  // push back fake indices if no leptons in genevent
43  if( !genEvt->isFullLeptonic() || !genEvt->lepton() || !genEvt->leptonBar() ){
44  match.push_back( -1 );
45  match.push_back( -1 );
46  match.push_back( -1 );
47  match.push_back( -1 );
48  }
49  else if(genEvt->isFullLeptonic(WDecay::kElec, WDecay::kElec) && elecs->size()>=2){
50  //search indices for electrons
51  int iLepBar = findMatchingLepton(genEvt->leptonBar(), elecs);
52  setCandidate(elecs, iLepBar, leptonBar_);
53  match.push_back( iLepBar );
54  int iLep = findMatchingLepton(genEvt->lepton(), elecs);
55  setCandidate(elecs, iLep, lepton_);
56  match.push_back( iLep );
57 
58  // fake indices for muons
59  match.push_back( -1 );
60  match.push_back( -1 );
61  }
62  else if(genEvt->isFullLeptonic(WDecay::kElec, WDecay::kMuon) && !elecs->empty() && !mus->empty() ){
63  if(genEvt->leptonBar()->isElectron()){
64  // push back index for e+
65  int iLepBar = findMatchingLepton(genEvt->leptonBar(), elecs);
66  setCandidate(elecs, iLepBar, leptonBar_);
67  match.push_back( iLepBar );
68  // push back fake indices for e- and mu+
69  match.push_back( -1 );
70  match.push_back( -1 );
71  // push back index for mu-
72  int iLep = findMatchingLepton(genEvt->lepton(), mus);
73  setCandidate(mus, iLep, lepton_);
74  match.push_back( iLep );
75  }
76  else{
77  // push back fake index for e+
78  match.push_back( -1 );
79  // push back index for e-
80  int iLepBar = findMatchingLepton(genEvt->leptonBar(), mus);
81  setCandidate(mus, iLepBar, leptonBar_);
82  match.push_back( iLepBar );
83  // push back index for mu+
84  int iLep = findMatchingLepton(genEvt->lepton(), elecs);
85  setCandidate(elecs, iLep, lepton_);
86  match.push_back( iLep );
87  // push back fake index for mu-
88  match.push_back( -1 );
89  }
90  }
91  else if(genEvt->isFullLeptonic(WDecay::kMuon, WDecay::kMuon) && mus->size()>=2 ){
92  // fake indices for electrons
93  match.push_back( -1 );
94  match.push_back( -1 );
95 
96  //search indices for electrons
97  int iLepBar = findMatchingLepton(genEvt->leptonBar(), mus);
98  setCandidate(mus, iLepBar, leptonBar_);
99  match.push_back( iLepBar );
100  int iLep = findMatchingLepton(genEvt->lepton(), mus);
101  setCandidate(mus, iLep, lepton_);
102  match.push_back( iLep );
103  }
104  else{ //this 'else' should happen if at least one genlepton is a tau
105  match.push_back( -1 );
106  match.push_back( -1 );
107  match.push_back( -1 );
108  match.push_back( -1 );
109  }
110 
111  // -----------------------------------------------------
112  // add met and neutrinos
113  // -----------------------------------------------------
114  if( !mets->empty() ){
115  //setCandidate(mets, 0, met_);
116  buildMatchingNeutrinos(evt, mets);
117  }
118 }
119 
120 template <typename O>
121 int
122 TtFullLepHypGenMatch::findMatchingLepton(const reco::GenParticle* genLep, const edm::Handle<std::vector<O> >& leps)
123 {
124  int idx=-1;
125  double minDR = -1;
126  for(unsigned i=0; i<leps->size(); ++i){
127  double dR = deltaR(genLep->eta(), genLep->phi(), (*leps)[i].eta(), (*leps)[i].phi());
128  if(minDR<0 || dR<minDR){
129  minDR=dR;
130  idx=i;
131  }
132  }
133  return idx;
134 }
135 
136 void
137 TtFullLepHypGenMatch::buildMatchingNeutrinos(edm::Event& evt, const edm::Handle<std::vector<pat::MET> >& mets)
138 {
139  // get genEvent
141  evt.getByLabel("genEvt", genEvt);
142 
143  if( genEvt->isTtBar() && genEvt->isFullLeptonic() && genEvt->neutrino() && genEvt->neutrinoBar() ){
144  double momXNu = genEvt->neutrino() ->px();
145  double momYNu = genEvt->neutrino() ->py();
146  double momXNuBar = genEvt->neutrinoBar()->px();
147  double momYNuBar = genEvt->neutrinoBar()->py();
148 
149  double momXMet = mets->at(0).px();
150  double momYMet = mets->at(0).py();
151 
152  double momXNeutrino = 0.5*(momXNu - momXNuBar + momXMet);
153  double momYNeutrino = 0.5*(momYNu - momYNuBar + momYMet);
154  double momXNeutrinoBar = momXMet - momXNeutrino;
155  double momYNeutrinoBar = momYMet - momYNeutrino;
156 
157  math::XYZTLorentzVector recNuFM(momXNeutrino,momYNeutrino,0,sqrt(momXNeutrino * momXNeutrino + momYNeutrino * momYNeutrino));
158  recNu = new reco::LeafCandidate(0,recNuFM);
159 
160  math::XYZTLorentzVector recNuBarFM(momXNeutrinoBar,momYNeutrinoBar,0,sqrt(momXNeutrinoBar * momXNeutrinoBar + momYNeutrinoBar * momYNeutrinoBar));
161  recNuBar = new reco::LeafCandidate(0,recNuBarFM);
162  }
163 }
int i
Definition: DBlmapReader.cc:9
reco::ShallowClonePtrCandidate * leptonBar_
reco::LeafCandidate * recNuBar
void setCandidate(const edm::Handle< C > &handle, const int &idx, reco::ShallowClonePtrCandidate *&clone)
use one object in a collection to set a ShallowClonePtrCandidate
reco::LeafCandidate * recNu
candidates needed for the genmatch hypothesis
virtual double eta() const
momentum pseudorapidity
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:30
T sqrt(T t)
Definition: SSEVec.h:46
vector< PseudoJet > jets
TtFullLepHypGenMatch(const edm::ParameterSet &)
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
virtual void buildHypo(edm::Event &evt, const edm::Handle< std::vector< pat::Electron > > &elecs, const edm::Handle< std::vector< pat::Muon > > &mus, const edm::Handle< std::vector< pat::Jet > > &jets, const edm::Handle< std::vector< pat::MET > > &mets, std::vector< int > &match, const unsigned int iComb)
build event hypothesis from the reco objects of a semi-leptonic event
void buildMatchingNeutrinos(edm::Event &, const edm::Handle< std::vector< pat::MET > > &)
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
reco::ShallowClonePtrCandidate * lepton_
int findMatchingLepton(const reco::GenParticle *, const edm::Handle< std::vector< O > > &)
reco::ShallowClonePtrCandidate * bBar_
std::string jetCorrectionLevel_
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:6
reco::ShallowClonePtrCandidate * b_
virtual double phi() const
momentum azimuthal angle
bool isValid(const int &idx, const edm::Handle< std::vector< pat::Jet > > &jets)
check if index is in valid range of selected jets