CMS 3D CMS Logo

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

#include <TtSemiLepHypWMassDeltaTopMass.h>

Inheritance diagram for TtSemiLepHypWMassDeltaTopMass:
TtSemiLepHypothesis edm::EDProducer edm::ProducerBase edm::ProductRegistryHelper

Public Member Functions

 TtSemiLepHypWMassDeltaTopMass (const edm::ParameterSet &)
 
 ~TtSemiLepHypWMassDeltaTopMass ()
 
- 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 ()
 
virtual ~EDProducer ()
 
- Public Member Functions inherited from edm::ProducerBase
 ProducerBase ()
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
boost::function< void(const
BranchDescription &)> 
registrationCallback () const
 used by the fwk to register list of products More...
 
virtual ~ProducerBase ()
 

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...
 

Private Attributes

std::string bTagAlgorithm_
 
double maxBDiscLightJets_
 
int maxNJets_
 
double minBDiscBJets_
 
int neutrinoSolutionType_
 
bool useBTagging_
 
double wMass_
 

Additional Inherited Members

- Public Types inherited from edm::EDProducer
typedef EDProducer ModuleType
 
typedef WorkerT< EDProducerWorkerType
 
- Public Types inherited from edm::ProducerBase
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 
- Static Public Member Functions inherited from edm::EDProducer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from TtSemiLepHypothesis
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
 return key More...
 
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::EDProducer
CurrentProcessingContext const * currentContext () const
 
- Protected Member Functions inherited from edm::ProducerBase
template<class TProducer , class TMethod >
void callWhenNewProductsRegistered (TProducer *iProd, TMethod iMethod)
 
- Protected Attributes inherited from TtSemiLepHypothesis
bool getMatch_
 
reco::ShallowClonePtrCandidatehadronicB_
 
std::string jetCorrectionLevel_
 
edm::InputTag jets_
 input label for all necessary collections More...
 
int key_
 hypothesis key (to be set by the buildKey function) More...
 
edm::InputTag leps_
 
reco::ShallowClonePtrCandidatelepton_
 
reco::ShallowClonePtrCandidateleptonicB_
 
reco::ShallowClonePtrCandidatelightQ_
 
reco::ShallowClonePtrCandidatelightQBar_
 
edm::InputTag match_
 
edm::InputTag mets_
 
reco::ShallowClonePtrCandidateneutrino_
 
int numberOfRealNeutrinoSolutions_
 

Detailed Description

Definition at line 7 of file TtSemiLepHypWMassDeltaTopMass.h.

Constructor & Destructor Documentation

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

Definition at line 4 of file TtSemiLepHypWMassDeltaTopMass.cc.

References edm::hlt::Exception, and maxNJets_.

4  :
5  TtSemiLepHypothesis( cfg ),
6  maxNJets_ (cfg.getParameter<int> ("maxNJets" )),
7  wMass_ (cfg.getParameter<double> ("wMass" )),
8  useBTagging_ (cfg.getParameter<bool> ("useBTagging" )),
9  bTagAlgorithm_ (cfg.getParameter<std::string>("bTagAlgorithm" )),
10  minBDiscBJets_ (cfg.getParameter<double> ("minBDiscBJets" )),
11  maxBDiscLightJets_ (cfg.getParameter<double> ("maxBDiscLightJets" )),
12  neutrinoSolutionType_(cfg.getParameter<int> ("neutrinoSolutionType"))
13 {
14  if(maxNJets_<4 && maxNJets_!=-1)
15  throw cms::Exception("WrongConfig")
16  << "Parameter maxNJets can not be set to " << maxNJets_ << ". \n"
17  << "It has to be larger than 4 or can be set to -1 to take all jets.";
18 }
T getParameter(std::string const &) const
TtSemiLepHypothesis(const edm::ParameterSet &)
default constructor
TtSemiLepHypWMassDeltaTopMass::~TtSemiLepHypWMassDeltaTopMass ( )

Definition at line 20 of file TtSemiLepHypWMassDeltaTopMass.cc.

20 { }

Member Function Documentation

void TtSemiLepHypWMassDeltaTopMass::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 23 of file TtSemiLepHypWMassDeltaTopMass.cc.

References bTagAlgorithm_, TtSemiLepEvtPartons::HadB, TtSemiLepHypothesis::hadronicB_, i, TtSemiLepHypothesis::isValid(), TtSemiLepHypothesis::jetCorrectionLevel(), analyzePatCleaning_cfg::jets, TtSemiLepEvtPartons::LepB, TtSemiLepEvtPartons::Lepton, TtSemiLepHypothesis::lepton_, TtSemiLepHypothesis::leptonicB_, TtSemiLepEvtPartons::LightQ, TtSemiLepHypothesis::lightQ_, TtSemiLepEvtPartons::LightQBar, TtSemiLepHypothesis::lightQBar_, maxBDiscLightJets_, maxNJets_, minBDiscBJets_, TtSemiLepHypothesis::neutrino_, neutrinoSolutionType_, reco::LeafCandidate::p4(), TtSemiLepHypothesis::setCandidate(), TtSemiLepHypothesis::setNeutrino(), useBTagging_, and wMass_.

28 {
29  if(leps->empty() || mets->empty() || jets->size()<4){
30  // create empty hypothesis
31  return;
32  }
33 
34  int maxNJets = maxNJets_;
35  if(maxNJets_ == -1 || (int)jets->size() < maxNJets_) maxNJets = jets->size();
36 
37  std::vector<bool> isBJet;
38  std::vector<bool> isLJet;
39  int cntBJets = 0;
40  if(useBTagging_) {
41  for(int idx=0; idx<maxNJets; ++idx) {
42  isBJet.push_back( ((*jets)[idx].bDiscriminator(bTagAlgorithm_) > minBDiscBJets_ ) );
43  isLJet.push_back( ((*jets)[idx].bDiscriminator(bTagAlgorithm_) < maxBDiscLightJets_) );
44  if((*jets)[idx].bDiscriminator(bTagAlgorithm_) > minBDiscBJets_ )cntBJets++;
45  }
46  }
47 
48  match.clear();
49  for(unsigned int i=0; i<5; ++i)
50  match.push_back(-1);
51 
52  // -----------------------------------------------------
53  // add lepton
54  // -----------------------------------------------------
55  setCandidate(leps, 0, lepton_);
57 
58  // -----------------------------------------------------
59  // add neutrino
60  // -----------------------------------------------------
61  if(neutrinoSolutionType_ == -1)
62  setCandidate(mets, 0, neutrino_);
63  else
64  setNeutrino(mets, leps, 0, neutrinoSolutionType_);
65 
66  // -----------------------------------------------------
67  // associate those jets that get closest to the W mass
68  // with their invariant mass to the hadronic W boson
69  // -----------------------------------------------------
70  double wDist =-1.;
71  std::vector<int> closestToWMassIndices;
72  closestToWMassIndices.push_back(-1);
73  closestToWMassIndices.push_back(-1);
74  for(int idx=0; idx<maxNJets; ++idx){
75  if(useBTagging_ && (!isLJet[idx] || (cntBJets<=2 && isBJet[idx]))) continue;
76  for(int jdx=(idx+1); jdx<maxNJets; ++jdx){
77  if(useBTagging_ && (!isLJet[jdx] || (cntBJets<=2 && isBJet[jdx]) || (cntBJets==3 && isBJet[idx] && isBJet[jdx]))) continue;
79  (*jets)[idx].p4()+
80  (*jets)[jdx].p4();
81  if( wDist<0. || wDist>fabs(sum.mass()-wMass_) ){
82  wDist=fabs(sum.mass()-wMass_);
83  closestToWMassIndices.clear();
84  closestToWMassIndices.push_back(idx);
85  closestToWMassIndices.push_back(jdx);
86  }
87  }
88  }
89 
90  // -----------------------------------------------------
91  // associate those jets to the hadronic and the leptonic
92  // b quark that minimize the difference between
93  // hadronic and leptonic top mass
94  // -----------------------------------------------------
95  double deltaTop=-1.;
96  int hadB=-1;
97  int lepB=-1;
98  if( isValid(closestToWMassIndices[0], jets) && isValid(closestToWMassIndices[1], jets)) {
101  (*jets)[closestToWMassIndices[0]].p4()+
102  (*jets)[closestToWMassIndices[1]].p4();
103  // find hadronic b candidate
104  for(int idx=0; idx<maxNJets; ++idx){
105  if(useBTagging_ && !isBJet[idx]) continue;
106  // make sure it's not used up already from the hadronic W
107  if( idx!=closestToWMassIndices[0] && idx!=closestToWMassIndices[1] ){
108  reco::Particle::LorentzVector hadTop = hadW + (*jets)[idx].p4();
109  // find leptonic b candidate
110  for(int jdx=0; jdx<maxNJets; ++jdx){
111  if(useBTagging_ && !isBJet[jdx]) continue;
112  // make sure it's not used up already from the hadronic branch
113  if( jdx!=closestToWMassIndices[0] && jdx!=closestToWMassIndices[1] && jdx!=idx ){
114  reco::Particle::LorentzVector lepTop = lepW + (*jets)[jdx].p4();
115  if( deltaTop<0. || deltaTop>fabs(hadTop.mass()-lepTop.mass()) ){
116  deltaTop=fabs(hadTop.mass()-lepTop.mass());
117  hadB=idx;
118  lepB=jdx;
119  }
120  }
121  }
122  }
123  }
124  }
125 
126  // -----------------------------------------------------
127  // add jets
128  // -----------------------------------------------------
129  if( isValid(closestToWMassIndices[0], jets) ){
130  setCandidate(jets, closestToWMassIndices[0], lightQ_, jetCorrectionLevel("wQuarkMix"));
131  match[TtSemiLepEvtPartons::LightQ] = closestToWMassIndices[0];
132  }
133 
134  if( isValid(closestToWMassIndices[1], jets) ){
135  setCandidate(jets, closestToWMassIndices[1], lightQBar_, jetCorrectionLevel("wQuarkMix"));
136  match[TtSemiLepEvtPartons::LightQBar] = closestToWMassIndices[1];
137  }
138 
139  if( isValid(hadB, jets) ){
140  setCandidate(jets, hadB, hadronicB_, jetCorrectionLevel("bQuark"));
142  }
143 
144  if( isValid(lepB, jets) ){
145  setCandidate(jets, lepB, leptonicB_, jetCorrectionLevel("bQuark"));
147  }
148 
149 }
int i
Definition: DBlmapReader.cc:9
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_
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 const LorentzVector & p4() const
four-momentum Lorentz vector
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:26
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 TtSemiLepHypWMassDeltaTopMass::buildKey ( )
inlineprivatevirtual

build the event hypothesis key

Implements TtSemiLepHypothesis.

Definition at line 17 of file TtSemiLepHypWMassDeltaTopMass.h.

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

int key_
hypothesis key (to be set by the buildKey function)

Member Data Documentation

std::string TtSemiLepHypWMassDeltaTopMass::bTagAlgorithm_
private

Definition at line 30 of file TtSemiLepHypWMassDeltaTopMass.h.

Referenced by buildHypo().

double TtSemiLepHypWMassDeltaTopMass::maxBDiscLightJets_
private

Definition at line 32 of file TtSemiLepHypWMassDeltaTopMass.h.

Referenced by buildHypo().

int TtSemiLepHypWMassDeltaTopMass::maxNJets_
private

Definition at line 27 of file TtSemiLepHypWMassDeltaTopMass.h.

Referenced by buildHypo(), and TtSemiLepHypWMassDeltaTopMass().

double TtSemiLepHypWMassDeltaTopMass::minBDiscBJets_
private

Definition at line 31 of file TtSemiLepHypWMassDeltaTopMass.h.

Referenced by buildHypo().

int TtSemiLepHypWMassDeltaTopMass::neutrinoSolutionType_
private

Definition at line 33 of file TtSemiLepHypWMassDeltaTopMass.h.

Referenced by buildHypo().

bool TtSemiLepHypWMassDeltaTopMass::useBTagging_
private

Definition at line 29 of file TtSemiLepHypWMassDeltaTopMass.h.

Referenced by buildHypo().

double TtSemiLepHypWMassDeltaTopMass::wMass_
private

Definition at line 28 of file TtSemiLepHypWMassDeltaTopMass.h.

Referenced by buildHypo().