CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
JetFlavourIdentifier.cc
Go to the documentation of this file.
1 //
2 // Translation of BTag MCJetFlavour tool to identify real flavour of a jet
3 // work with CaloJet objects
4 // Store Infos by Values in JetFlavour.h
5 // Author: Attilio
6 // Date: 05.10.2007
7 //
8 //
9 // \class JetFlavourIdentifier
10 //
11 // \brief Interface to pull out the proper flavour identifier from a
12 // jet->parton matching collection.
13 //
14 // In detail, the definitions are as follows:
15 //
16 // Definitions:
17 // The default behavior is that the "definition" is NULL_DEF,
18 // so the software will fall back to either "physics" or "algorithmic" definition
19 // as per the "physDefinition" switch.
20 //
21 // If the user specifies "definition", then that definition is taken. However,
22 // if the requested definition is not defined, the software reverts back to
23 // either "physics" or "algorithmic" based on the "physDefinition" switch.
24 // For example, if the user specifies "heaviest" as a flavor ID, and there
25 // are no bottom, charm, or top quarks in the event that match to the jet,
26 // then the software will fall back to the "physics" or "algorithmic" definition.
27 //
28 // Modifications:
29 //
30 // 09.03.2008: Sal Rappoccio.
31 // Added capability for all methods of identification, not just
32 // "physics" or "algorithmic". If the requested method does not exist
33 // (i.e. is unphysical), the "physics" or "algorithmic" definition
34 // is defaulted to.
35 
36 //=======================================================================
37 
38 // user include files
43 
49 
50 //#include "SimGeneral/HepPDTRecord/interface/ParticleDataTable.h"
51 
52 
60 
64 
67 
68 #include <memory>
69 #include <string>
70 #include <iostream>
71 #include <vector>
72 #include <Math/VectorUtil.h>
73 #include <TMath.h>
74 
75 using namespace std;
76 using namespace reco;
77 using namespace edm;
78 using namespace ROOT::Math::VectorUtil;
79 
80 namespace reco { namespace modules {
81 
82 //--------------------------------------------------------------------------
83 //
84 //--------------------------------------------------------------------------
86 {
87  public:
88  enum DEFINITION_T { PHYSICS=0, ALGO, NEAREST_STATUS2, NEAREST_STATUS3, HEAVIEST,
90  NULL_DEF};
91 
94 
95  private:
96  virtual void produce(edm::Event&, const edm::EventSetup& ) override;
97 
98  JetFlavour::Leptons findLeptons(const GenParticleRef &);
99  std::vector<const reco::Candidate*> findCandidates(const reco::Candidate*, int);
100  void fillLeptons(const std::vector<const reco::Candidate*> &, JetFlavour::Leptons &, int, int);
101  static int heaviestFlavour(int);
102 
109 
110 };
111 } }
113 
114 //=========================================================================
115 
116 JetFlavourIdentifier::JetFlavourIdentifier( const edm::ParameterSet& iConfig )
117 {
118  produces<JetFlavourMatchingCollection>();
119  sourceByReferToken_ = consumes<JetMatchedPartonsCollection>(iConfig.getParameter<InputTag>("srcByReference"));
120  physDefinition = iConfig.getParameter<bool>("physicsDefinition");
121  leptonInfo_ = iConfig.exists("leptonInfo") ? iConfig.getParameter<bool>("leptonInfo") : true;
122  // If we have a definition of which parton to identify, use it,
123  // otherwise we default to the "old" behavior of either "physics" or "algorithmic".
124  // Furthermore, if the specified definition is not sensible for the given jet,
125  // then the "physDefinition" switch is used to identify the flavour of the jet.
126  if ( iConfig.exists("definition") ) {
127  definition = static_cast<DEFINITION_T>( iConfig.getParameter<int>("definition") );
128  } else {
129  definition = NULL_DEF;
130  }
131 }
132 
133 //=========================================================================
134 
135 JetFlavourIdentifier::~JetFlavourIdentifier()
136 {
137 }
138 
139 // ------------ method called to produce the data ------------
140 
141 void JetFlavourIdentifier::produce( Event& iEvent, const EventSetup& iEs )
142 {
143  // Get the JetMatchedPartons
144  iEvent.getByToken (sourceByReferToken_, theTagByRef);
145 
146  // Create a JetFlavourMatchingCollection
148  if (!theTagByRef->empty()) {
149  RefToBase<Jet> jj = theTagByRef->begin()->first;
151  } else {
152  jfmc = new JetFlavourMatchingCollection();
153  }
154  auto_ptr<reco::JetFlavourMatchingCollection> jetFlavMatching(jfmc);
155 
156  // Loop over the matched partons and see which match.
157  for ( JetMatchedPartonsCollection::const_iterator j = theTagByRef->begin();
158  j != theTagByRef->end();
159  j ++ ) {
160 
161 
162  // Consider this match.
163  const MatchedPartons aMatch = (*j).second;
164 
165  // This will hold the 4-vector, vertex, flavour and the leptonian decays (0: no lepton, xyz: x leptons in layer 2, y in layer 1 and z in layer 0) of the requested object.
166  math::XYZTLorentzVector thePartonLorentzVector(0,0,0,0);
167  math::XYZPoint thePartonVertex(0,0,0);
168  int thePartonFlavour = 0;
169  JetFlavour::Leptons theLeptons;
170 
171  // get the partons based on which definition to use.
172  switch (definition) {
173  case PHYSICS: {
174  const GenParticleRef aPartPhy = aMatch.physicsDefinitionParton();
175  if (aPartPhy.isNonnull()) {
176  thePartonLorentzVector = aPartPhy.get()->p4();
177  thePartonVertex = aPartPhy.get()->vertex();
178  thePartonFlavour = aPartPhy.get()->pdgId();
179  if (leptonInfo_) theLeptons = findLeptons(aPartPhy);
180  }
181  break;
182  }
183  case ALGO: {
184  const GenParticleRef aPartAlg = aMatch.algoDefinitionParton();
185  if (aPartAlg.isNonnull()) {
186  thePartonLorentzVector = aPartAlg.get()->p4();
187  thePartonVertex = aPartAlg.get()->vertex();
188  thePartonFlavour = aPartAlg.get()->pdgId();
189  if (leptonInfo_) theLeptons = findLeptons(aPartAlg);
190  }
191  break;
192  }
193  case NEAREST_STATUS2 : {
194  const GenParticleRef aPartN2 = aMatch.nearest_status2();
195  if (aPartN2.isNonnull()) {
196  thePartonLorentzVector = aPartN2.get()->p4();
197  thePartonVertex = aPartN2.get()->vertex();
198  thePartonFlavour = aPartN2.get()->pdgId();
199  if (leptonInfo_) theLeptons = findLeptons(aPartN2);
200  }
201  break;
202  }
203  case NEAREST_STATUS3: {
204  const GenParticleRef aPartN3 = aMatch.nearest_status3();
205  if (aPartN3.isNonnull()) {
206  thePartonLorentzVector = aPartN3.get()->p4();
207  thePartonVertex = aPartN3.get()->vertex();
208  thePartonFlavour = aPartN3.get()->pdgId();
209  if (leptonInfo_) theLeptons = findLeptons(aPartN3);
210  }
211  break;
212  }
213  case HEAVIEST: {
214  const GenParticleRef aPartHeaviest = aMatch.heaviest();
215  if (aPartHeaviest.isNonnull()) {
216  thePartonLorentzVector = aPartHeaviest.get()->p4();
217  thePartonVertex = aPartHeaviest.get()->vertex();
218  thePartonFlavour = aPartHeaviest.get()->pdgId();
219  if (leptonInfo_) theLeptons = findLeptons(aPartHeaviest);
220  }
221  break;
222  }
223  // Default case is backwards-compatible
224  default:{
225  if (physDefinition) {
226  const GenParticleRef aPartPhy = aMatch.physicsDefinitionParton();
227  if (aPartPhy.isNonnull()) {
228  thePartonLorentzVector = aPartPhy.get()->p4();
229  thePartonVertex = aPartPhy.get()->vertex();
230  thePartonFlavour = aPartPhy.get()->pdgId();
231  if (leptonInfo_) theLeptons = findLeptons(aPartPhy);
232  }
233  } else {
234  const GenParticleRef aPartAlg = aMatch.algoDefinitionParton();
235  if (aPartAlg.isNonnull()) {
236  thePartonLorentzVector = aPartAlg.get()->p4();
237  thePartonVertex = aPartAlg.get()->vertex();
238  thePartonFlavour = aPartAlg.get()->pdgId();
239  if (leptonInfo_) theLeptons = findLeptons(aPartAlg);
240  }
241  }
242  } break;
243  } // end switch on definition
244 
245  // Now make sure we have a match. If the user specified "heaviest", for instance,
246  // and there is no b- or c-quarks, then fall back to the "physDefinition" switch.
247 
248  if (thePartonFlavour == 0) {
249  if (physDefinition) {
250  const GenParticleRef aPartPhy = aMatch.physicsDefinitionParton();
251  if (aPartPhy.isNonnull()) {
252  thePartonLorentzVector = aPartPhy.get()->p4();
253  thePartonVertex = aPartPhy.get()->vertex();
254  thePartonFlavour = aPartPhy.get()->pdgId();
255  if (leptonInfo_) theLeptons = findLeptons(aPartPhy);
256  }
257  } else {
258  const GenParticleRef aPartAlg = aMatch.algoDefinitionParton();
259  if (aPartAlg.isNonnull()) {
260  thePartonLorentzVector = aPartAlg.get()->p4();
261  thePartonVertex = aPartAlg.get()->vertex();
262  thePartonFlavour = aPartAlg.get()->pdgId();
263  if (leptonInfo_) theLeptons = findLeptons(aPartAlg);
264  }
265  }
266  }
267 
268 /*
269  std::cout << "Leptons of " <<thePartonFlavour << " Jet: " << std::endl;
270  std::cout << " electrons: " <<theLeptons.electron << std::endl;
271  std::cout << " muons : " <<theLeptons.muon << std::endl;
272  std::cout << " tau : " <<theLeptons.tau << std::endl;
273 */
274  // Add the jet->flavour match to the map.
275  (*jetFlavMatching)[(*j).first] = JetFlavour(thePartonLorentzVector, thePartonVertex, thePartonFlavour, theLeptons);
276  }// end loop over jets
277 
278 
279  // Put the object into the event.
280  iEvent.put( jetFlavMatching );
281 
282 }
283 
284 JetFlavour::Leptons JetFlavourIdentifier::findLeptons(const GenParticleRef &parton)
285 {
286  JetFlavour::Leptons theLeptons;
287 
288  thePartonLV = parton->p4();
289 
291  const reco::Candidate *mcstring = parton->daughter(0);
292  int partonFlavour = std::abs(parton->pdgId());
293 // std::cout << "parton DeltaR: " << DeltaR(thePartonLV, parton->p4()) << std::endl;
294 
296  std::vector<const reco::Candidate*> candidates = findCandidates(mcstring, partonFlavour);
297 // std::cout << "Candidates are:" << std::endl;
298 // for(unsigned int j = 0; j < candidates.size(); j++) std::cout << " --> " << candidates[j]->pdgId() << std::endl;
299 
301  fillLeptons(candidates, theLeptons, 1, partonFlavour);
302 
303  return theLeptons;
304 }
305 
306 std::vector<const reco::Candidate*> JetFlavourIdentifier::findCandidates(const reco::Candidate *cand, int partonFlavour)
307 {
308  std::vector<const reco::Candidate*> cands;
309 
310  for(unsigned int i = 0; i < cand->numberOfDaughters(); i++) {
311 /*
312  std::cout << "DeltaR - " << partonFlavour << " ";
313  if (DeltaR(thePartonLV, cand->daughter(i)->p4()) > 0.7) std::cout << "(";
314  std::cout << cand->daughter(i)->pdgId() << ": " << DeltaR(thePartonLV, cand->daughter(i)->p4());
315  if (DeltaR(thePartonLV, cand->daughter(i)->p4()) > 0.7) std::cout << ")";
316  std::cout << std::endl;
317 */
318 
319  if (DeltaR(thePartonLV, cand->daughter(i)->p4()) < 0.7) {
320  int pdgId = std::abs(cand->daughter(i)->pdgId());
321  int flavour = heaviestFlavour(pdgId);
322  if (flavour == partonFlavour ||
323  (flavour >= 10 && partonFlavour >= 10)) {
324 // std::cout << "<------- " << std::endl;
325  std::vector<const reco::Candidate*> newcands = findCandidates(cand->daughter(i), partonFlavour);
326 // std::cout << " ------->" << std::endl;
327  std::copy(newcands.begin(), newcands.end(), std::back_inserter(cands));
328  }
329  if (partonFlavour >= 10)
330  cands.push_back(cand->daughter(i));
331  }
332  }
333 
334  if (cands.empty() && std::abs(cand->pdgId()) > 110 &&
335  !(partonFlavour >= 4 && partonFlavour < 10 &&
336  heaviestFlavour(cand->pdgId()) < 4))
337  cands.push_back(cand);
338 
339  return cands;
340 }
341 
342 void JetFlavourIdentifier::fillLeptons(const std::vector<const reco::Candidate*> &cands, JetFlavour::Leptons &leptons, int rank, int flavour)
343 {
344  for(unsigned int j = 0; j < cands.size(); j++) {
345  for(unsigned int i = 0; i < cands[j]->numberOfDaughters(); i++) {
346  int pdgId = std::abs(cands[j]->daughter(i)->pdgId());
347 
348 // for(int z = 1; z <= rank; z *= 10) std::cout << " ------ ";
349 // std::cout << pdgId << std::endl;
350 
352  if (pdgId == 12)
353  leptons.electron += rank;
354  else if (pdgId == 14)
355  leptons.muon += rank;
356  else if (pdgId == 16)
357  leptons.tau += rank;
358  else {
359  int heaviest = heaviestFlavour(pdgId);
360  int heaviest_ = heaviest < 10 ? heaviest : 0;
361  if (!heaviest || (flavour < 4 ? (heaviest_ < 4) : (heaviest >= 4))) {
362  std::vector<const reco::Candidate*> newcands = findCandidates(cands[j]->daughter(i), heaviest);
363  if (pdgId <= 110) newcands.push_back(cands[j]->daughter(i));
364  fillLeptons(newcands, leptons, rank * 10, std::max(heaviest_, flavour));
365  }
366  }
367  }
368  }
369 }
370 
371 int JetFlavourIdentifier::heaviestFlavour(int pdgId)
372 {
373  int flavour = 0;
374 
375  pdgId = std::abs(pdgId) % 100000;
376  if (pdgId > 110) {
377  while(pdgId % 10 > 0 && pdgId % 10 < 6) {
378  pdgId /= 10;
379  if (pdgId % 10 > flavour)
380  flavour = pdgId % 10;
381  }
382  } else
383  flavour = pdgId;
384 
385  return flavour;
386 }
387 
388 
389 //define this as a plug-in
T getParameter(std::string const &) const
Handle< JetMatchedPartonsCollection > theTagByRef
int i
Definition: DBlmapReader.cc:9
virtual const Candidate * daughter(size_type i) const =0
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
Definition: deltaR.h:79
const GenParticleRef & nearest_status3() const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
bool exists(std::string const &parameterName) const
checks if a parameter exists
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:250
virtual size_type numberOfDaughters() const =0
number of daughters
int iEvent
Definition: GenABIO.cc:243
const T & max(const T &a, const T &b)
const GenParticleRef & algoDefinitionParton() const
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int j
Definition: DBlmapReader.cc:9
const GenParticleRef & nearest_status2() const
virtual int pdgId() const =0
PDG identifier.
lepton info struct
Definition: JetFlavour.h:25
const GenParticleRef heaviest() const
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
const GenParticleRef & physicsDefinitionParton() const
T const * get() const
Returns C++ pointer to the item.
Definition: Ref.h:242
EDGetTokenT< JetMatchedPartonsCollection > sourceByReferToken_
int flavour(const Candidate &part)
Definition: pdgIdUtils.h:31
virtual const LorentzVector & p4() const =0
four-momentum Lorentz vector