CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
JetPartonMatcher.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 JetPartonMatcher
10 //
11 // \brief Interface to create a specified "matched" parton to jet match based
12 // on the user interpretation.
13 //
14 // Algorithm for definitions:
15 //
16 // If the particle is matched to any item in the priority list (given by user),
17 // then count that prioritized particle as the match.
18 //
19 // If not resort to default behavior:
20 //
21 //
22 // 1) Algorithmic Definition:
23 //
24 // A particle is a parton if its daughter is a string(code=92) or a cluster(code=93)
25 // If (parton is within the cone defined by coneSizeToAssociate) then:
26 // if (parton is a b) then associatedParton is the b
27 // else if (associatedParton =! b and parton is a c) then associatedParton is the c
28 // else if (associatedParton =! b and associatedParton =! c) then associatedParton is the one with the highest pT
29 // associatedParton can be -1 --> no partons in the cone
30 // True Flavour of the jet --> flavour of the associatedParton
31 //
32 // ToDo: if more than one b(c) in the cone --> the selected parton is not always the one with highest pT
33 //
34 //
35 // 2) Physics Definition:
36 //
37 // A particle is a parton if its daughter is a string(code=92) or a cluster(code=93)
38 // if( only one initialParticle within the cone defined by coneSizeToAssociate) associatedInitialParticle is the initialParticle
39 // TheBiggerConeSize = 0.7 --> it's hard coded!
40 // if( a parton not coming from associatedInitialParticle in TheBiggerConeSize is a b or a c) reject the association
41 //
42 // associatedInitialParticle can be -1 --> no initialParticle in the cone or rejected association
43 // True Flavour of the jet --> flavour of the associatedInitialParticle
44 //
45 //
46 // Modifications:
47 //
48 // 09.03.2008: Sal Rappoccio.
49 // Added capability for "priority" list.
50 //
51 
52 //=======================================================================
53 
54 // user include files
59 
65 
73 
77 //#include "DataFormats/Candidate/interface/Candidate.h"
78 //#include "DataFormats/Candidate/interface/CandidateFwd.h"
81 
82 #include <memory>
83 #include <string>
84 #include <iostream>
85 #include <vector>
86 #include <Math/VectorUtil.h>
87 #include <TMath.h>
88 
89 using namespace std;
90 using namespace reco;
91 using namespace edm;
92 using namespace ROOT::Math::VectorUtil;
93 
95 {
96  public:
99 
100  private:
101  virtual void produce(edm::Event&, const edm::EventSetup& ) override;
102 
103  int fillAlgoritDefinition( const Jet& );
104  int fillPhysicsDefinition( const Jet& );
109 
111 
117  vector<int> priorityList;
118 
119 };
120 
121 //=========================================================================
122 
124  theHeaviest(0),
125  theNearest2(0),
126  theNearest3(0),
127  theHardest(0)
128 {
129  produces<JetMatchedPartonsCollection>();
130  m_jetsSrcToken = consumes<edm::View <reco::Jet> >(iConfig.getParameter<edm::InputTag>("jets"));
131  m_ParticleSrcToken = consumes<GenParticleRefVector>(iConfig.getParameter<edm::InputTag>("partons"));
132  coneSizeToAssociate = iConfig.getParameter<double>("coneSizeToAssociate");
133  if ( iConfig.exists("doPriority") ) {
134  doPriority = iConfig.getParameter<bool>("doPriority");
135  priorityList = iConfig.getParameter<vector<int> >("priorityList");
136  } else {
137  doPriority = false;
138  priorityList.clear();
139  }
140 }
141 
142 //=========================================================================
143 
145 {
146 }
147 
148 // ------------ method called to produce the data ------------
149 
151 {
153  iEvent.getByToken(m_jetsSrcToken, jets_h );
155 
156  edm::LogVerbatim("JetPartonMatcher") << "=== Partons size:" << particles->size();
157 
158  for( size_t m = 0; m != particles->size(); ++ m ) {
159  const GenParticle & aParton = *(particles->at(m).get());
160  edm::LogVerbatim("JetPartonMatcher") << aParton.status() << " " <<
161  aParton.pdgId() << " " <<
162  aParton.pt() << " " <<
163  aParton.eta() << " " <<
164  aParton.phi() << endl;
165  }
166 
167  auto_ptr<reco::JetMatchedPartonsCollection> jetMatchedPartons( new JetMatchedPartonsCollection(reco::JetRefBaseProd(jets_h)));
168 
169  for (size_t j = 0; j < jets_h->size(); j++) {
170 
171  const int theMappedPartonAlg = fillAlgoritDefinition( (*jets_h)[j] );
172  const int theMappedPartonPhy = fillPhysicsDefinition( (*jets_h)[j] );
173 
174  GenParticleRef pHV;
175  GenParticleRef pN2;
176  GenParticleRef pN3;
177  GenParticleRef pPH;
178  GenParticleRef pAL;
179 
180  if(theHeaviest>=0) pHV = particles->at( theHeaviest );
181  if(theNearest2>=0) pN2 = particles->at( theNearest2 );
182  if(theNearest3>=0) pN3 = particles->at( theNearest3 );
183  if(theMappedPartonPhy>=0) pPH = particles->at( theMappedPartonPhy );
184  if(theMappedPartonAlg>=0) pAL = particles->at( theMappedPartonAlg );
185 
186  (*jetMatchedPartons)[jets_h->refAt(j)]=MatchedPartons(pHV,pN2,pN3,pPH,pAL);
187  }
188 
189  iEvent.put( jetMatchedPartons );
190 
191 }
192 
193 //
194 // Algorithmic Definition:
195 // Output: define one associatedParton
196 // Loop on all particle.
197 //
198 // (Note: This part added by Salvatore Rappoccio)
199 // [begin]
200 // If the particle is matched to any item in the priority list (given by user),
201 // then count that prioritized particle as the match.
202 //
203 // If not resort to default behavior:
204 // [end]
205 //
206 // A particle is a parton if its daughter is a string(code=92) or a cluster(code=93)
207 // If (parton is within the cone defined by coneSizeToAssociate) then:
208 // if (parton is a b) then associatedParton is the b
209 // else if (associatedParton =! b and parton is a c) then associatedParton is the c
210 // else if (associatedParton =! b and associatedParton =! c) then associatedParton is the one with the highest pT
211 // associatedParton can be -1 --> no partons in the cone
212 // True Flavour of the jet --> flavour of the associatedParton
213 //
214 // ToDo: if more than one b(c) in the cone --> the selected parton is not always the one with highest pT
215 //
217 
218  int tempParticle = -1;
219  int tempPartonHighestPt = -1;
220  int tempNearest = -1;
221  float maxPt = 0;
222  float minDr = 1000;
223  bool foundPriority = false;
224 
225  // Loop over the particles in question until we find a priority
226  // "hit", or if we find none in the priority list (or we don't want
227  // to consider priority), then loop through all particles and fill
228  // standard definition.
229 
230  //
231  // Matching:
232  //
233  // 1) First try to match by hand. The "priority list" is given
234  // by the user. The algorithm finds any particles in that
235  // "priority list" that are within the specified cone size.
236  // If it finds one, it counts the object as associating to that
237  // particle.
238  // NOTE! The objects in the priority list are given in order
239  // of priority. So if the user specifies:
240  // 6, 24, 21,
241  // then it will first look for top quarks, then W bosons, then gluons.
242  // 2) If no priority items are found, do the default "standard"
243  // matching.
244  for( size_t m = 0; m != particles->size() && !foundPriority; ++ m ) {
245  const Candidate & aParton = *(particles->at(m).get());
246 
247  // "Priority" behavoir:
248  // Associate to the first particle found in the priority list, regardless
249  // of delta R.
250  if ( doPriority ) {
251  int ipdgId = aParton.pdgId();
252  vector<int>::const_iterator ipriority = find( priorityList.begin(), priorityList.end(), abs(ipdgId) );
253  // Check to see if this particle is in our priority list
254  if ( ipriority != priorityList.end() ) {
255  // This particle is on our priority list. If it matches,
256  // we break, since we take in order of priority, not deltaR
257  double dist = DeltaR( theJet.p4(), aParton.p4() );
258  if( dist <= coneSizeToAssociate ) {
259  tempParticle = m;
260  foundPriority = true;
261  }
262  }
263  }
264  // Here we do not want to do priority matching. Ensure the "foundPriority" swtich
265  // is turned off.
266  else {
267  foundPriority = false;
268  }
269 
270  // Default behavior:
271  // Look for partons before the color string to associate.
272  // Do this if we don't want to do priority matching, or if
273  // we didn't find a priority object.
274 
275  if( !foundPriority && aParton.status() != 3 ) {
276  double dist = DeltaR( theJet.p4(), aParton.p4() );
277  if( dist <= coneSizeToAssociate ) {
278  if( dist < minDr ) {
279  minDr = dist;
280  tempNearest = m;
281  }
282  if( tempParticle == -1 && ( abs( aParton.pdgId() ) == 4 ) ) tempParticle = m;
283  if( abs( aParton.pdgId() ) == 5 ) tempParticle = m;
284  if( aParton.pt() > maxPt ) {
285  maxPt = aParton.pt();
286  tempPartonHighestPt = m;
287  }
288  }
289  }
290  }
291 
292  if ( foundPriority ) {
293  theHeaviest = tempParticle; // The priority-matched particle
294  theHardest = -1; // set the rest to -1
295  theNearest2 = -1; // " "
296  } else {
297  theHeaviest = tempParticle;
298  theHardest = tempPartonHighestPt;
299  theNearest2 = tempNearest;
300  if ( tempParticle == -1 ) tempParticle = tempPartonHighestPt;
301  }
302  return tempParticle;
303 }
304 
305 //
306 // Physics Definition:
307 // A initialParticle is a particle with status=3
308 // Output: define one associatedInitialParticle
309 // Loop on all particles
310 //
311 // (Note: This part added by Salvatore Rappoccio)
312 // [begin]
313 // If the particle is matched to any item in the priority list (given by user),
314 // then count that prioritized particle as the match.
315 //
316 // If not resort to default behavior:
317 // [end]
318 //
319 // A particle is a parton if its daughter is a string(code=92) or a cluster(code=93)
320 // if( only one initialParticle within the cone defined by coneSizeToAssociate) associatedInitialParticle is the initialParticle
321 // TheBiggerConeSize = 0.7 --> it's hard coded!
322 // if( a parton not coming from associatedInitialParticle in TheBiggerConeSize is a b or a c) reject the association
323 //
324 // associatedInitialParticle can be -1 --> no initialParticle in the cone or rejected association
325 // True Flavour of the jet --> flavour of the associatedInitialParticle
326 //
328 
329  float TheBiggerConeSize = 0.7; // In HepMC it's 0.3 --> it's a mistake: value has to be 0.7
330  int tempParticle = -1;
331  int nInTheCone = 0;
332  int tempNearest = -1;
333  float minDr = 1000;
334  bool foundPriority = false;
335 
336  vector<const reco::Candidate *> theContaminations;
337  theContaminations.clear();
338 
339  // Loop over the particles in question until we find a priority
340  // "hit", or if we find none in the priority list (or we don't want
341  // to consider priority), then loop through all particles and fill
342  // standard definition.
343 
344  //
345  // Matching:
346  //
347  // 1) First try to match by hand. The "priority list" is given
348  // by the user. The algorithm finds any particles in that
349  // "priority list" that are within the specified cone size.
350  // If it finds one, it counts the object as associating to that
351  // particle.
352  // NOTE! The objects in the priority list are given in order
353  // of priority. So if the user specifies:
354  // 6, 24, 21,
355  // then it will first look for top quarks, then W bosons, then gluons.
356  // 2) If no priority items are found, do the default "standard"
357  // matching.
358  for( size_t m = 0; m != particles->size() && !foundPriority; ++ m ) {
359 
360  const Candidate & aParticle = *(particles->at(m).get());
361 
362  // "Priority" behavoir:
363  // Associate to the first particle found in the priority list, regardless
364  // of delta R.
365  if ( doPriority ) {
366  int ipdgId = aParticle.pdgId();
367  vector<int>::const_iterator ipriority = find( priorityList.begin(), priorityList.end(), abs(ipdgId) );
368  // Check to see if this particle is in our priority list
369  if ( ipriority != priorityList.end() ) {
370  // This particle is on our priority list. If it matches,
371  // we break, since we take in order of priority, not deltaR
372  double dist = DeltaR( theJet.p4(), aParticle.p4() );
373  if( dist <= coneSizeToAssociate ) {
374  tempParticle = m;
375  foundPriority = true;
376  }
377  }
378  }
379  // Here we do not want to do priority matching. Ensure the "foundPriority" swtich
380  // is turned off.
381  else {
382  foundPriority = false;
383  }
384 
385  // Default behavior:
386  if ( !foundPriority ) {
387 
388  // skipping all particle but udscbg (is this correct/enough?!?!)
389  bool isAParton = false;
390  int flavour = abs(aParticle.pdgId());
391  if(flavour == 1 ||
392  flavour == 2 ||
393  flavour == 3 ||
394  flavour == 4 ||
395  flavour == 5 ||
396  flavour == 21 ) isAParton = true;
397  if(!isAParton) continue;
398  double dist = DeltaR( theJet.p4(), aParticle.p4() );
399  if( aParticle.status() == 3 && dist < minDr ) {
400  minDr = dist;
401  tempNearest = m;
402  }
403  if( aParticle.status() == 3 && dist <= coneSizeToAssociate ) {
404  //cout << "particle in small cone=" << aParticle.pdgId() << endl;
405  tempParticle = m;
406  nInTheCone++;
407  }
408  // Look for heavy partons in TheBiggerConeSize now
409  if( aParticle.numberOfDaughters() > 0 && aParticle.status() != 3 ) {
410  if( flavour == 1 ||
411  flavour == 2 ||
412  flavour == 3 ||
413  flavour == 21 ) continue;
414  if( dist < TheBiggerConeSize ) theContaminations.push_back( &aParticle );
415  }
416  }
417  }
418 
419 
420  // Here's the default behavior for assignment if there is no priority.
421  if ( !foundPriority ) {
422  theNearest3 = tempNearest;
423 
424  if(nInTheCone != 1) return -1; // rejected --> only one initialParton requested
425  if(theContaminations.size() == 0 ) return tempParticle; //no contamination
426  int initialPartonFlavour = abs( (particles->at(tempParticle).get()) ->pdgId() );
427 
428  vector<const Candidate *>::const_iterator itCont = theContaminations.begin();
429  for( ; itCont != theContaminations.end(); itCont++ ) {
430  int contaminatingFlavour = abs( (*itCont)->pdgId() );
431  if( (*itCont)->numberOfMothers()>0 && (*itCont)->mother(0) == particles->at(tempParticle).get() ) continue; // mother is the initialParton --> OK
432  if( initialPartonFlavour == 4 ) {
433  if( contaminatingFlavour == 4 ) continue; // keep association --> the initialParton is a c --> the contaminated parton is a c
434  tempParticle = -1; // all the other cases reject!
435  return tempParticle;
436  }
437  }
438  }
439  // If there is priority, then just set the heaviest to priority, the rest are -1.
440  else {
441  theHeaviest = tempParticle; // Set the heaviest to tempParticle
442  theNearest2 = -1; // Set the rest to -1
443  theNearest3 = -1; // " "
444  theHardest = -1; // " "
445  }
446 
447  return tempParticle;
448 }
449 
450 //define this as a plug-in
452 
T getParameter(std::string const &) const
virtual int pdgId() const
PDG identifier.
virtual float pt() const
transverse momentum
int fillPhysicsDefinition(const Jet &)
Definition: deltaR.h:79
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
virtual int status() const
status word
virtual float phi() const
momentum azimuthal angle
Base class for all types of Jets.
Definition: Jet.h:20
virtual int status() const =0
status word
bool exists(std::string const &parameterName) const
checks if a parameter exists
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
int fillAlgoritDefinition(const Jet &)
virtual size_type numberOfDaughters() const =0
number of daughters
edm::EDGetTokenT< GenParticleRefVector > m_ParticleSrcToken
virtual float pt() const =0
transverse momentum
int iEvent
Definition: GenABIO.cc:230
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
virtual float eta() const
momentum pseudorapidity
Handle< GenParticleRefVector > particles
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int j
Definition: DBlmapReader.cc:9
JetPartonMatcher(const edm::ParameterSet &)
virtual int pdgId() const =0
PDG identifier.
vector< int > priorityList
virtual void produce(edm::Event &, const edm::EventSetup &) override
edm::EDGetTokenT< edm::View< reco::Jet > > m_jetsSrcToken
int flavour(const Candidate &part)
Definition: pdgIdUtils.h:31
virtual const LorentzVector & p4() const
four-momentum Lorentz vector
virtual const LorentzVector & p4() const =0
four-momentum Lorentz vector