CMS 3D CMS Logo

All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TtSemiLepHypGeom.cc
Go to the documentation of this file.
3 
5  TtSemiLepHypothesis( cfg ),
6  neutrinoSolutionType_(cfg.getParameter<int>("neutrinoSolutionType"))
7 { }
8 
10 
11 void
14  const edm::Handle<std::vector<pat::MET> >& mets,
15  const edm::Handle<std::vector<pat::Jet> >& jets,
16  std::vector<int>& match, const unsigned int iComb)
17 {
18  // -----------------------------------------------------
19  // add jets
20  // -----------------------------------------------------
21  for(unsigned idx=0; idx<match.size(); ++idx){
22  if( isValid(match[idx], jets) ){
23  switch(idx){
25  setCandidate(jets, match[idx], lightQ_, jetCorrectionLevel("wQuarkMix")); break;
27  setCandidate(jets, match[idx], lightQBar_, jetCorrectionLevel("wQuarkMix")); break;
29  setCandidate(jets, match[idx], hadronicB_, jetCorrectionLevel("bQuark")); break;
31  setCandidate(jets, match[idx], leptonicB_, jetCorrectionLevel("bQuark")); break;
32  }
33  }
34  }
35 
36  // -----------------------------------------------------
37  // add lepton
38  // -----------------------------------------------------
39  if( leps->empty() )
40  return;
41  setCandidate(leps, 0, lepton_);
42  match.push_back( 0 );
43 
44  // -----------------------------------------------------
45  // add neutrino
46  // -----------------------------------------------------
47  if( mets->empty() )
48  return;
49  if(neutrinoSolutionType_ == -1)
50  setCandidate(mets, 0, neutrino_);
51  else
52  setNeutrino(mets, leps, 0, neutrinoSolutionType_);
53 }
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_
TtSemiLepHypGeom(const edm::ParameterSet &)
reco::ShallowClonePtrCandidate * lightQBar_
reco::ShallowClonePtrCandidate * neutrino_
reco::ShallowClonePtrCandidate * hadronicB_
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
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
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 ...