CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TtSemiLepHypothesis.cc
Go to the documentation of this file.
5 
7 
10  jets_(cfg.getParameter<edm::InputTag>("jets")),
11  leps_(cfg.getParameter<edm::InputTag>("leps")),
12  mets_(cfg.getParameter<edm::InputTag>("mets")),
13  numberOfRealNeutrinoSolutions_(-1),
14  lightQ_(0), lightQBar_(0), hadronicB_(0),
15  leptonicB_(0), neutrino_(0), lepton_(0)
16 {
17  getMatch_ = false;
18  if( cfg.exists("match") ) {
19  getMatch_ = true;
20  match_ = cfg.getParameter<edm::InputTag>("match");
21  }
22  if( cfg.exists("jetCorrectionLevel") ) {
23  jetCorrectionLevel_ = cfg.getParameter<std::string>("jetCorrectionLevel");
24  }
25  produces<std::vector<std::pair<reco::CompositeCandidate, std::vector<int> > > >();
26  produces<int>("Key");
27  produces<int>("NumberOfRealNeutrinoSolutions");
28 }
29 
32 {
33  if( lightQ_ ) delete lightQ_;
34  if( lightQBar_ ) delete lightQBar_;
35  if( hadronicB_ ) delete hadronicB_;
36  if( leptonicB_ ) delete leptonicB_;
37  if( neutrino_ ) delete neutrino_;
38  if( lepton_ ) delete lepton_;
39 }
40 
42 void
44 {
46  evt.getByLabel(jets_, jets);
47 
49  evt.getByLabel(leps_, leps);
50 
52  evt.getByLabel(mets_, mets);
53 
54  std::vector<std::vector<int> > matchVec;
55  if( getMatch_ ) {
57  evt.getByLabel(match_, matchHandle);
58  matchVec = *matchHandle;
59  }
60  else {
61  std::vector<int> dummyMatch;
62  for(unsigned int i = 0; i < 4; ++i)
63  dummyMatch.push_back( -1 );
64  matchVec.push_back( dummyMatch );
65  }
66 
67  // declare auto_ptr for products
68  std::auto_ptr<std::vector<std::pair<reco::CompositeCandidate, std::vector<int> > > >
69  pOut( new std::vector<std::pair<reco::CompositeCandidate, std::vector<int> > > );
70  std::auto_ptr<int> pKey(new int);
71  std::auto_ptr<int> pNeutrinoSolutions(new int);
72 
73  // go through given vector of jet combinations
74  unsigned int idMatch = 0;
75  typedef std::vector<std::vector<int> >::iterator MatchVecIterator;
76  for(MatchVecIterator match = matchVec.begin(); match != matchVec.end(); ++match) {
77  // reset pointers
79  // build hypothesis
80  buildHypo(evt, leps, mets, jets, *match, idMatch++);
81  pOut->push_back( std::make_pair(hypo(), *match) );
82  }
83  // feed out hyps and matches
84  evt.put(pOut);
85 
86  // build and feed out key
87  buildKey();
88  *pKey=key();
89  evt.put(pKey, "Key");
90 
91  // feed out number of real neutrino solutions
92  *pNeutrinoSolutions=numberOfRealNeutrinoSolutions_;
93  evt.put(pNeutrinoSolutions, "NumberOfRealNeutrinoSolutions");
94 }
95 
97 void
99 {
101  lightQ_ = 0;
102  lightQBar_ = 0;
103  hadronicB_ = 0;
104  leptonicB_ = 0;
105  neutrino_ = 0;
106  lepton_ = 0;
107 }
108 
112 {
113  // check for sanity of the hypothesis
114  if( !lightQ_ || !lightQBar_ || !hadronicB_ ||
115  !leptonicB_ || !neutrino_ || !lepton_ )
116  return reco::CompositeCandidate();
117 
118  // setup transient references
119  reco::CompositeCandidate hyp, hadTop, hadW, lepTop, lepW;
120 
121  AddFourMomenta addFourMomenta;
122  // build up the top branch that decays leptonically
125  addFourMomenta.set( lepW );
126  lepTop.addDaughter( lepW, TtSemiLepDaughter::LepW );
128  addFourMomenta.set( lepTop );
129 
130  // build up the top branch that decays hadronically
133  addFourMomenta.set( hadW );
134  hadTop.addDaughter( hadW, TtSemiLepDaughter::HadW );
136  addFourMomenta.set( hadTop );
137 
138  // build ttbar hypotheses
139  hyp.addDaughter( lepTop, TtSemiLepDaughter::LepTop );
140  hyp.addDaughter( hadTop, TtSemiLepDaughter::HadTop );
141  addFourMomenta.set( hyp );
142 
143  return hyp;
144 }
145 
149 {
150  // check whetherwe are dealing with a reco muon or a reco electron
152  if( dynamic_cast<const reco::Muon*>(cand) ){
153  type = WDecay::kMuon;
154  }
155  else if( dynamic_cast<const reco::GsfElectron*>(cand) ){
156  type = WDecay::kElec;
157  }
158  return type;
159 }
160 
162 std::string
163 TtSemiLepHypothesis::jetCorrectionLevel(const std::string& quarkType)
164 {
165  // jetCorrectionLevel was not configured
166  if(jetCorrectionLevel_.empty())
167  throw cms::Exception("Configuration")
168  << "Unconfigured jetCorrectionLevel. Please use an appropriate, non-empty string.\n";
169 
170  // quarkType is unknown
171  if( !(quarkType=="wQuarkMix" ||
172  quarkType=="udsQuark" ||
173  quarkType=="cQuark" ||
174  quarkType=="bQuark") )
175  throw cms::Exception("Configuration")
176  << quarkType << " is unknown as a quarkType for the jetCorrectionLevel.\n";
177 
178  // combine correction level; start with a ':' even if
179  // there is no flavor tag to be added, as it is needed
180  // by setCandidate to disentangle the correction tag
181  // from a potential flavor tag, which can be empty
182  std::string level=jetCorrectionLevel_+":";
183  if( level=="L5Flavor:" || level=="L6UE:" || level=="L7Parton:" ){
184  if(quarkType=="wQuarkMix"){level+="wMix";}
185  if(quarkType=="udsQuark" ){level+="uds";}
186  if(quarkType=="cQuark" ){level+="charm";}
187  if(quarkType=="bQuark" ){level+="bottom";}
188  }
189  else{
190  level+="none";
191  }
192  return level;
193 }
194 
196 void
197 TtSemiLepHypothesis::setCandidate(const edm::Handle<std::vector<pat::Jet> >& handle, const int& idx,
198  reco::ShallowClonePtrCandidate*& clone, const std::string& correctionLevel)
199 {
201  // disentangle the correction from the potential flavor tag
202  // by the separating ':'; the flavor tag can be empty though
203  std::string step = correctionLevel.substr(0,correctionLevel.find(":"));
204  std::string flavor = correctionLevel.substr(1+correctionLevel.find(":"));
205  float corrFactor = 1.;
206  if(flavor=="wMix")
207  corrFactor = 0.75*ptr->jecFactor(step, "uds") + 0.25*ptr->jecFactor(step, "charm");
208  else
209  corrFactor = ptr->jecFactor(step, flavor);
210  clone = new reco::ShallowClonePtrCandidate( ptr, ptr->charge(), ptr->p4()*corrFactor, ptr->vertex() );
211 }
212 
214 void TtSemiLepHypothesis::setNeutrino(const edm::Handle<std::vector<pat::MET> >& met,
215  const edm::Handle<edm::View<reco::RecoCandidate> >& leps, const int& idx, const int& type)
216 {
218  MEzCalculator mez;
219  mez.SetMET( *(met->begin()) );
220  if( leptonType( &(leps->front()) ) == WDecay::kMuon )
221  mez.SetLepton( (*leps)[idx], true );
222  else if( leptonType( &(leps->front()) ) == WDecay::kElec )
223  mez.SetLepton( (*leps)[idx], false );
224  else
225  throw cms::Exception("UnimplementedFeature") << "Type of lepton given together with MET for solving neutrino kinematics is neither muon nor electron.\n";
226  double pz = mez.Calculate(type);
227  numberOfRealNeutrinoSolutions_ = mez.IsComplex() ? 0 : 2;
228  const math::XYZTLorentzVector p4( ptr->px(), ptr->py(), pz, sqrt(ptr->px()*ptr->px() + ptr->py()*ptr->py() + pz*pz) );
229  neutrino_ = new reco::ShallowClonePtrCandidate( ptr, ptr->charge(), p4, ptr->vertex() );
230 }
std::string jetCorrectionLevel_
type
Definition: HCALResponse.h:22
T getParameter(std::string const &) const
static const std::string HadB
int i
Definition: DBlmapReader.cc:9
static const std::string LepTop
static const std::string Lep
static const std::string LepW
WDecay::LeptonType leptonType(const reco::RecoCandidate *cand)
determine lepton type of reco candidate and return a corresponding WDecay::LeptonType; the type is kN...
reco::ShallowClonePtrCandidate * lepton_
virtual void buildHypo(edm::Event &event, const edm::Handle< edm::View< reco::RecoCandidate > > &lepton, const edm::Handle< std::vector< pat::MET > > &neutrino, const edm::Handle< std::vector< pat::Jet > > &jets, std::vector< int > &jetPartonAssociation, const unsigned int iComb)=0
build event hypothesis from the reco objects of a semi-leptonic event
reco::ShallowClonePtrCandidate * lightQBar_
bool exists(std::string const &parameterName) const
checks if a parameter exists
static const std::string HadP
~TtSemiLepHypothesis()
default destructor
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:30
TtSemiLepHypothesis(const edm::ParameterSet &)
default constructor
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:84
reco::ShallowClonePtrCandidate * neutrino_
T sqrt(T t)
Definition: SSEVec.h:28
double p4[4]
Definition: TauolaWrapper.h:92
void resetCandidates()
reset candidate pointers before hypo build process
tuple handle
Definition: patZpeak.py:22
edm::InputTag jets_
input label for all necessary collections
virtual void buildKey()=0
build the event hypothesis key
static const std::string HadTop
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:359
void addDaughter(const Candidate &, const std::string &s="")
add a clone of the passed candidate as daughter
reco::ShallowClonePtrCandidate * hadronicB_
virtual void produce(edm::Event &, const edm::EventSetup &)
produce the event hypothesis as CompositeCandidate and Key
static const std::string HadQ
T * clone(const T *tp)
Definition: Ptr.h:42
static const std::string HadW
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
reco::CompositeCandidate hypo()
return event hypothesis
static const std::string Nu
int key() const
return key
void setCandidate(const edm::Handle< C > &handle, const int &idx, reco::ShallowClonePtrCandidate *&clone)
use one object in a collection to set a ShallowClonePtrCandidate
static const std::string LepB
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:6
void set(reco::Candidate &c) const
set up a candidate
tuple level
Definition: testEve_cfg.py:34
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 ...