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