CMS 3D CMS Logo

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

#include <TtSemiLepHypGenMatch.h>

Inheritance diagram for TtSemiLepHypGenMatch:
TtSemiLepHypothesis edm::EDProducer edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

 TtSemiLepHypGenMatch (const edm::ParameterSet &)
 
 ~TtSemiLepHypGenMatch ()
 
- Public Member Functions inherited from TtSemiLepHypothesis
 TtSemiLepHypothesis (const edm::ParameterSet &)
 default constructor More...
 
 ~TtSemiLepHypothesis ()
 default destructor More...
 
- Public Member Functions inherited from edm::EDProducer
 EDProducer ()
 
ModuleDescription const & moduleDescription () const
 
virtual ~EDProducer ()
 
- Public Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 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
std::vector< ConsumesInfoconsumesInfo () const
 
 EDConsumerBase ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesDependentUpon (std::string const &iProcessName, std::string const &iModuleLabel, bool iPrint, std::vector< char const * > &oModuleLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Protected Attributes

edm::EDGetTokenT< TtGenEventgenEvtToken_
 
- Protected Attributes inherited from TtSemiLepHypothesis
bool getMatch_
 
reco::ShallowClonePtrCandidatehadronicB_
 
std::string jetCorrectionLevel_
 
edm::EDGetTokenT< std::vector
< pat::Jet > > 
jetsToken_
 input label for all necessary collections More...
 
int key_
 hypothesis key (to be set by the buildKey function) More...
 
edm::EDGetTokenT< edm::View
< reco::RecoCandidate > > 
lepsToken_
 
reco::ShallowClonePtrCandidatelepton_
 
reco::ShallowClonePtrCandidateleptonicB_
 
reco::ShallowClonePtrCandidatelightQ_
 
reco::ShallowClonePtrCandidatelightQBar_
 
edm::EDGetTokenT< std::vector
< std::vector< int > > > 
matchToken_
 
edm::EDGetTokenT< std::vector
< pat::MET > > 
metsToken_
 
reco::ShallowClonePtrCandidateneutrino_
 
int neutrinoSolutionType_
 algorithm used to calculate neutrino solutions (see cfi for further details) More...
 
edm::EDGetTokenT< int > nJetsConsideredToken_
 
int numberOfRealNeutrinoSolutions_
 

Private Member Functions

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 More...
 
virtual void buildKey ()
 build the event hypothesis key More...
 
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; return -1 if this fails More...
 

Additional Inherited Members

- Public Types inherited from edm::EDProducer
typedef EDProducer ModuleType
 
- Public Types inherited from edm::ProducerBase
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- 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 TtSemiLepHypothesis
void buildHypo (const edm::Handle< edm::View< reco::RecoCandidate > > &leps, const edm::Handle< std::vector< pat::MET > > &mets, const edm::Handle< std::vector< pat::Jet > > &jets, std::vector< int > &jetPartonAssociation)
 minimalistic build function for simple hypotheses More...
 
reco::CompositeCandidate hypo ()
 return event hypothesis More...
 
bool isValid (const int &idx, const edm::Handle< std::vector< pat::Jet > > &jets)
 check if index is in valid range of selected jets More...
 
std::string jetCorrectionLevel (const std::string &quarkType)
 helper function to construct the proper correction level string for corresponding quarkType More...
 
int key () const
 
WDecay::LeptonType leptonType (const reco::RecoCandidate *cand)
 determine lepton type of reco candidate and return a corresponding WDecay::LeptonType; the type is kNone if it is whether a muon nor an electron More...
 
virtual void produce (edm::Event &, const edm::EventSetup &)
 produce the event hypothesis as CompositeCandidate and Key More...
 
void resetCandidates ()
 reset candidate pointers before hypo build process More...
 
template<typename C >
void setCandidate (const edm::Handle< C > &handle, const int &idx, reco::ShallowClonePtrCandidate *&clone)
 use one object in a collection to set a ShallowClonePtrCandidate More...
 
void setCandidate (const edm::Handle< std::vector< pat::Jet > > &handle, const int &idx, reco::ShallowClonePtrCandidate *&clone, const std::string &correctionLevel)
 use one object in a jet collection to set a ShallowClonePtrCandidate with proper jet corrections More...
 
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 More...
 
- 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 7 of file TtSemiLepHypGenMatch.h.

Constructor & Destructor Documentation

TtSemiLepHypGenMatch::TtSemiLepHypGenMatch ( const edm::ParameterSet cfg)
explicit

Definition at line 6 of file TtSemiLepHypGenMatch.cc.

6  :
7  TtSemiLepHypothesis( cfg ),
8  genEvtToken_( consumes<TtGenEvent>( edm::InputTag( "genEvt" ) ) )
9 { }
TtSemiLepHypothesis(const edm::ParameterSet &)
default constructor
edm::EDGetTokenT< TtGenEvent > genEvtToken_
TtSemiLepHypGenMatch::~TtSemiLepHypGenMatch ( )

Definition at line 11 of file TtSemiLepHypGenMatch.cc.

11 { }

Member Function Documentation

void TtSemiLepHypGenMatch::buildHypo ( edm::Event evt,
const edm::Handle< edm::View< reco::RecoCandidate > > &  leps,
const edm::Handle< std::vector< pat::MET > > &  mets,
const edm::Handle< std::vector< pat::Jet > > &  jets,
std::vector< int > &  match,
const unsigned int  iComb 
)
privatevirtual

build event hypothesis from the reco objects of a semi-leptonic event

Implements TtSemiLepHypothesis.

Definition at line 14 of file TtSemiLepHypGenMatch.cc.

References funct::abs(), findMatchingLepton(), TtGenEvtProducer_cfi::genEvt, genEvtToken_, edm::Event::getByToken(), TtSemiLepEvtPartons::HadB, TtSemiLepHypothesis::hadronicB_, customizeTrackingMonitorSeedNumber::idx, TtSemiLepHypothesis::isValid(), TtSemiLepHypothesis::jetCorrectionLevel(), fwrapper::jets, TtSemiLepEvtPartons::LepB, TtSemiLepHypothesis::lepton_, TtSemiLepHypothesis::leptonicB_, TtSemiLepEvtPartons::LightQ, TtSemiLepHypothesis::lightQ_, TtSemiLepEvtPartons::LightQBar, TtSemiLepHypothesis::lightQBar_, TtSemiLepHypothesis::neutrino_, TtSemiLepHypothesis::neutrinoSolutionType_, TtSemiLepHypothesis::setCandidate(), and TtSemiLepHypothesis::setNeutrino().

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 }
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:462
reco::ShallowClonePtrCandidate * lepton_
reco::ShallowClonePtrCandidate * lightQBar_
reco::ShallowClonePtrCandidate * neutrino_
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_
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
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_
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 ...
virtual void TtSemiLepHypGenMatch::buildKey ( )
inlineprivatevirtual

build the event hypothesis key

Implements TtSemiLepHypothesis.

Definition at line 17 of file TtSemiLepHypGenMatch.h.

References TtSemiLepHypothesis::key_, and TtEvent::kGenMatch.

int key_
hypothesis key (to be set by the buildKey function)
int TtSemiLepHypGenMatch::findMatchingLepton ( const edm::Handle< TtGenEvent > &  genEvt,
const edm::Handle< edm::View< reco::RecoCandidate > > &  leps 
)
private

find index of the candidate nearest to the singleLepton of the generator event in the collection; return -1 if this fails

Definition at line 75 of file TtSemiLepHypGenMatch.cc.

References deltaR(), PFRecoTauDiscriminationAgainstElectronDeadECAL_cfi::dR, i, TtSemiLepHypothesis::leptonType(), and HLT_FULL_cff::minDR.

Referenced by buildHypo().

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
WDecay::LeptonType leptonType(const reco::RecoCandidate *cand)
determine lepton type of reco candidate and return a corresponding WDecay::LeptonType; the type is kN...
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17

Member Data Documentation

edm::EDGetTokenT<TtGenEvent> TtSemiLepHypGenMatch::genEvtToken_
protected

Definition at line 30 of file TtSemiLepHypGenMatch.h.

Referenced by buildHypo().