CMS 3D CMS Logo

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:
98  ~JetPartonMatcher() override;
99 
100  private:
101  void produce(edm::StreamID, edm::Event&, const edm::EventSetup& ) const override;
102 
105  int theHeaviest = 0;
106  int theNearest2 =0;
107  int theNearest3 = 0;
108  int theHardest = 0;
109  };
110 
111  int fillAlgoritDefinition( const Jet&, WorkingVariables& ) const;
112  int fillPhysicsDefinition( const Jet&, WorkingVariables& ) const;
113 
119  vector<int> priorityList;
120 
121 };
122 
123 //=========================================================================
124 
126 {
127  produces<JetMatchedPartonsCollection>();
128  m_jetsSrcToken = consumes<edm::View <reco::Jet> >(iConfig.getParameter<edm::InputTag>("jets"));
129  m_ParticleSrcToken = consumes<GenParticleRefVector>(iConfig.getParameter<edm::InputTag>("partons"));
130  coneSizeToAssociate = iConfig.getParameter<double>("coneSizeToAssociate");
131  if ( iConfig.exists("doPriority") ) {
132  doPriority = iConfig.getParameter<bool>("doPriority");
133  priorityList = iConfig.getParameter<vector<int> >("priorityList");
134  } else {
135  doPriority = false;
136  priorityList.clear();
137  }
138 }
139 
140 //=========================================================================
141 
143 {
144 }
145 
146 // ------------ method called to produce the data ------------
147 
149 {
150  WorkingVariables wv;
152  iEvent.getByToken(m_jetsSrcToken, jets_h );
153  iEvent.getByToken(m_ParticleSrcToken, wv.particles );
154 
155  edm::LogVerbatim("JetPartonMatcher") << "=== Partons size:" << wv.particles->size();
156 
157  for( size_t m = 0; m < wv.particles->size(); ++ m ) {
158  const GenParticle & aParton = *(wv.particles->at(m).get());
159  edm::LogVerbatim("JetPartonMatcher") << aParton.status() << " " <<
160  aParton.pdgId() << " " <<
161  aParton.pt() << " " <<
162  aParton.eta() << " " <<
163  aParton.phi() << endl;
164  }
165 
166  auto jetMatchedPartons = std::make_unique<JetMatchedPartonsCollection>(reco::JetRefBaseProd(jets_h));
167 
168  for (size_t j = 0; j < jets_h->size(); j++) {
169 
170  const int theMappedPartonAlg = fillAlgoritDefinition( (*jets_h)[j], wv );
171  const int theMappedPartonPhy = fillPhysicsDefinition( (*jets_h)[j], wv );
172 
173  GenParticleRef pHV;
174  GenParticleRef pN2;
175  GenParticleRef pN3;
176  GenParticleRef pPH;
177  GenParticleRef pAL;
178 
179  if(wv.theHeaviest>=0) pHV = wv.particles->at( wv.theHeaviest );
180  if(wv.theNearest2>=0) pN2 = wv.particles->at( wv.theNearest2 );
181  if(wv.theNearest3>=0) pN3 = wv.particles->at( wv.theNearest3 );
182  if(theMappedPartonPhy>=0) pPH = wv.particles->at( theMappedPartonPhy );
183  if(theMappedPartonAlg>=0) pAL = wv.particles->at( theMappedPartonAlg );
184 
185  (*jetMatchedPartons)[jets_h->refAt(j)]=MatchedPartons(pHV,pN2,pN3,pPH,pAL);
186  }
187 
188  iEvent.put(std::move(jetMatchedPartons));
189 
190 }
191 
192 //
193 // Algorithmic Definition:
194 // Output: define one associatedParton
195 // Loop on all particle.
196 //
197 // (Note: This part added by Salvatore Rappoccio)
198 // [begin]
199 // If the particle is matched to any item in the priority list (given by user),
200 // then count that prioritized particle as the match.
201 //
202 // If not resort to default behavior:
203 // [end]
204 //
205 // A particle is a parton if its daughter is a string(code=92) or a cluster(code=93)
206 // If (parton is within the cone defined by coneSizeToAssociate) then:
207 // if (parton is a b) then associatedParton is the b
208 // else if (associatedParton =! b and parton is a c) then associatedParton is the c
209 // else if (associatedParton =! b and associatedParton =! c) then associatedParton is the one with the highest pT
210 // associatedParton can be -1 --> no partons in the cone
211 // True Flavour of the jet --> flavour of the associatedParton
212 //
213 // ToDo: if more than one b(c) in the cone --> the selected parton is not always the one with highest pT
214 //
216 
217  int tempParticle = -1;
218  int tempPartonHighestPt = -1;
219  int tempNearest = -1;
220  float maxPt = 0;
221  float minDr = 1000;
222  bool foundPriority = false;
223 
224  // Loop over the particles in question until we find a priority
225  // "hit", or if we find none in the priority list (or we don't want
226  // to consider priority), then loop through all particles and fill
227  // standard definition.
228 
229  //
230  // Matching:
231  //
232  // 1) First try to match by hand. The "priority list" is given
233  // by the user. The algorithm finds any particles in that
234  // "priority list" that are within the specified cone size.
235  // If it finds one, it counts the object as associating to that
236  // particle.
237  // NOTE! The objects in the priority list are given in order
238  // of priority. So if the user specifies:
239  // 6, 24, 21,
240  // then it will first look for top quarks, then W bosons, then gluons.
241  // 2) If no priority items are found, do the default "standard"
242  // matching.
243  for( size_t m = 0; m < wv.particles->size() && !foundPriority; ++ m ) {
244  // const Candidate & aParton = *(wv.particles->at(m).get());
245  const GenParticle & aParton = *(wv.particles->at(m).get());
246  int flavour = abs(aParton.pdgId());
247 
248  // "Priority" behavoir:
249  // Associate to the first particle found in the priority list, regardless
250  // of delta R.
251  if ( doPriority ) {
252  vector<int>::const_iterator ipriority = find( priorityList.begin(), priorityList.end(), flavour );
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 && ( flavour == 4 ) ) tempParticle = m;
283  if( flavour == 5 ) tempParticle = m;
284  if( aParton.pt() > maxPt ) {
285  maxPt = aParton.pt();
286  tempPartonHighestPt = m;
287  }
288  }
289  }
290  }
291 
292  if ( foundPriority ) {
293  wv.theHeaviest = tempParticle; // The priority-matched particle
294  wv.theHardest = -1; // set the rest to -1
295  wv.theNearest2 = -1; // " "
296  } else {
297  wv.theHeaviest = tempParticle;
298  wv.theHardest = tempPartonHighestPt;
299  wv.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 < wv.particles->size() && !foundPriority; ++ m ) {
359 
360  // const Candidate & aParticle = *(wv.particles->at(m).get());
361  const GenParticle & aParticle = *(wv.particles->at(m).get());
362  int flavour = abs(aParticle.pdgId());
363 
364  // "Priority" behavoir:
365  // Associate to the first particle found in the priority list, regardless
366  // of delta R.
367  if ( doPriority ) {
368  vector<int>::const_iterator ipriority = find( priorityList.begin(), priorityList.end(), flavour );
369  // Check to see if this particle is in our priority list
370  if ( ipriority != priorityList.end() ) {
371  // This particle is on our priority list. If it matches,
372  // we break, since we take in order of priority, not deltaR
373  double dist = DeltaR( theJet.p4(), aParticle.p4() );
374  if( dist <= coneSizeToAssociate ) {
375  tempParticle = m;
376  foundPriority = true;
377  }
378  }
379  }
380  // Here we do not want to do priority matching. Ensure the "foundPriority" swtich
381  // is turned off.
382  else {
383  foundPriority = false;
384  }
385 
386  // Default behavior:
387  if ( !foundPriority ) {
388 
389  // skipping all particle but udscbg (is this correct/enough?!?!)
390  bool isAParton = false;
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  wv.theNearest3 = tempNearest;
423 
424  if(nInTheCone != 1) return -1; // rejected --> only one initialParton requested
425  if(theContaminations.empty() ) return tempParticle; //no contamination
426  int initialPartonFlavour = abs( (wv.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) == wv.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  wv.theHeaviest = tempParticle; // Set the heaviest to tempParticle
442  wv.theNearest2 = -1; // Set the rest to -1
443  wv.theNearest3 = -1; // " "
444  wv.theHardest = -1; // " "
445  }
446 
447  return tempParticle;
448 }
449 
450 //define this as a plug-in
452 
T getParameter(std::string const &) const
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
int pdgId() const final
PDG identifier.
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:136
double eta() const final
momentum pseudorapidity
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:519
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
int fillAlgoritDefinition(const Jet &, WorkingVariables &) const
bool exists(std::string const &parameterName) const
checks if a parameter exists
double pt() const final
transverse momentum
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
edm::RefToBaseProd< reco::Jet > JetRefBaseProd
Definition: JetCollection.h:14
edm::EDGetTokenT< GenParticleRefVector > m_ParticleSrcToken
int iEvent
Definition: GenABIO.cc:230
Definition: Jet.py:1
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const LorentzVector & p4() const final
four-momentum Lorentz vector
Definition: LeafCandidate.h:99
JetPartonMatcher(const edm::ParameterSet &)
vector< int > priorityList
size_t numberOfDaughters() const override
number of daughters
~JetPartonMatcher() override
fixed size matrix
HLT enums.
int status() const final
status word
int fillPhysicsDefinition(const Jet &, WorkingVariables &) const
edm::EDGetTokenT< edm::View< reco::Jet > > m_jetsSrcToken
double phi() const final
momentum azimuthal angle
def move(src, dest)
Definition: eostools.py:510
Handle< GenParticleRefVector > particles