CMS 3D CMS Logo

VertexClassifier.cc
Go to the documentation of this file.
1 /*
2  * VertexClassifier.C
3  */
4 
5 #include <cmath>
6 #include <cstdlib>
7 #include <iostream>
8 
9 #include "HepPDT/ParticleID.hh"
10 
12 
13 
14 #define update(a, b) do { (a) = (a) | (b); } while(0)
15 
16 
18  edm::ConsumesCollector&& collector) :
20  tracer_(config,std::move(collector)),
21  hepMCLabel_( config.getUntrackedParameter<edm::InputTag>("hepMC") )
22 {
23  collector.consumes<edm::HepMCProduct>(hepMCLabel_);
24  // Set the history depth after hadronization
25  tracer_.depth(-2);
26 
27  // Set the minimum decay length for detecting long decays
28  longLivedDecayLength_ = config.getUntrackedParameter<double>("longLivedDecayLength");
29 
30  // Set the distance for clustering vertices
31  vertexClusteringDistance_ = config.getUntrackedParameter<double>("vertexClusteringDistance");
32 }
33 
34 
36 {
37  // Get the new event information for the tracer
38  tracer_.newEvent(event, setup);
39 
40  // Get hepmc of the event
41  event.getByLabel(hepMCLabel_, mcInformation_);
42 
43  // Get the partivle data table
45 
46  // Create the list of primary vertices associated to the event
48 }
49 
50 
52 {
53  // Initializing the category vector
54  reset();
55 
56  // Associate and evaluate the vertex history (check for fakes)
57  if ( tracer_.evaluate(vertex) )
58  {
59  // Get all the information related to the simulation details
61 
62  // Get all the information related to decay process
64 
65  // Get information about conversion and other interactions
67 
68  // Get geometrical information about the vertices
70 
71  // Check for unkown classification
72  unknownVertex();
73  }
74  else
75  flags_[Fake] = true;
76 
77  return *this;
78 }
79 
80 
82 {
83  // Initializing the category vector
84  reset();
85 
86  // Trace the history for the given TP
87  tracer_.evaluate(vertex);
88 
89  // Check for a reconstructed track
90  if ( tracer_.recoVertex().isNonnull() )
91  flags_[Reconstructed] = true;
92  else
93  flags_[Reconstructed] = false;
94 
95  // Get all the information related to the simulation details
97 
98  // Get all the information related to decay process
100 
101  // Get information about conversion and other interactions
103 
104  // Get geometrical information about the vertices
106 
107  // Check for unkown classification
108  unknownVertex();
109 
110  return *this;
111 }
112 
113 
115 {
116  // Get the event id for the initial TP.
117  EncodedEventId eventId = tracer_.simVertex()->eventId();
118  // Check for signal events
119  flags_[SignalEvent] = !eventId.bunchCrossing() && !eventId.event();
120 }
121 
122 
124 {
125  // Get the generated vetices from track history
126  VertexHistory::GenVertexTrail const & genVertexTrail = tracer_.genVertexTrail();
127 
128  // Loop over the generated vertices
129  for (
130  VertexHistory::GenVertexTrail::const_iterator ivertex = genVertexTrail.begin();
131  ivertex != genVertexTrail.end();
132  ++ivertex
133  )
134  {
135  // Get the pointer to the vertex by removing the const-ness (no const methos in HepMC::GenVertex)
136  HepMC::GenVertex * vertex = const_cast<HepMC::GenVertex *>(*ivertex);
137 
138  // Loop over the sources looking for specific decays
139  for (
140  HepMC::GenVertex::particle_iterator iparent = vertex->particles_begin(HepMC::parents);
141  iparent != vertex->particles_end(HepMC::parents);
142  ++iparent
143  )
144  {
145  // Collect the pdgid of the parent
146  int pdgid = std::abs((*iparent)->pdg_id());
147  // Get particle type
149 
150  // Check if the particle type is valid one
151  if (particleID.isValid())
152  {
153  // Get particle data
154  ParticleData const * particleData = particleDataTable_->particle(particleID);
155  // Check if the particle exist in the table
156  if (particleData)
157  {
158  // Check if their life time is bigger than longLivedDecayLength_
159  if ( particleData->lifetime() > longLivedDecayLength_ )
160  {
161  // Check for B, C weak decays and long lived decays
162  update(flags_[BWeakDecay], particleID.hasBottom());
163  update(flags_[CWeakDecay], particleID.hasCharm());
164  update(flags_[LongLivedDecay], true);
165  }
166  // Check Tau, Ks and Lambda decay
167  update(flags_[TauDecay], pdgid == 15);
168  update(flags_[KsDecay], pdgid == 310);
169  update(flags_[LambdaDecay], pdgid == 3122);
170  update(flags_[JpsiDecay], pdgid == 443);
171  update(flags_[XiDecay], pdgid == 3312);
172  update(flags_[OmegaDecay], pdgid == 3334);
173  update(flags_[SigmaPlusDecay], pdgid == 3222);
174  update(flags_[SigmaMinusDecay], pdgid == 3112);
175  }
176  }
177  }
178  }
179 }
180 
181 
183 {
184  VertexHistory::SimVertexTrail const & simVertexTrail = tracer_.simVertexTrail();
185 
186  for (
187  VertexHistory::SimVertexTrail::const_iterator ivertex = simVertexTrail.begin();
188  ivertex != simVertexTrail.end();
189  ++ivertex
190  )
191  {
192  // pdgid of the real source parent vertex
193  int pdgid = 0;
194 
195  // select the original source in case of combined vertices
196  bool flag = false;
198 
199  for (its = (*ivertex)->sourceTracks_begin(); its != (*ivertex)->sourceTracks_end(); ++its)
200  {
201  for (itd = (*ivertex)->daughterTracks_begin(); itd != (*ivertex)->daughterTracks_end(); ++itd)
202  if (itd != its)
203  {
204  flag = true;
205  break;
206  }
207  if (flag)
208  break;
209  }
210  // Collect the pdgid of the original source track
211  if ( its != (*ivertex)->sourceTracks_end() )
212  pdgid = std::abs((*its)->pdgId());
213  else
214  pdgid = 0;
215 
216  // Geant4 process type is selected using first Geant4 vertex assigned to
217  // the TrackingVertex
218  unsigned int processG4 = 0;
219 
220  if((*ivertex)->nG4Vertices() > 0) {
221  processG4 = (*(*ivertex)->g4Vertices_begin()).processType();
222  }
223 
224  unsigned int process = g4toCMSProcMap_.processId(processG4);
225 
226  // Flagging all the different processes
227  update(
229  process != CMS::Undefined &&
230  process != CMS::Unknown &&
231  process != CMS::Primary
232  );
233 
238  update(flags_[DecayProcess], process == CMS::Decay);
241  update(flags_[EIoniProcess], process == CMS::EIoni);
242  update(flags_[HIoniProcess], process == CMS::HIoni);
243  update(flags_[MuIoniProcess], process == CMS::MuIoni);
244  update(flags_[PhotonProcess], process == CMS::Photon);
247  update(flags_[EBremProcess], process == CMS::EBrem);
249  update(flags_[MuBremProcess], process == CMS::MuBrem);
250  update(flags_[MuNuclProcess], process == CMS::MuNucl);
251 
252 
253  // Loop over the simulated particles
254  for (
255  TrackingVertex::tp_iterator iparticle = (*ivertex)->daughterTracks_begin();
256  iparticle != (*ivertex)->daughterTracks_end();
257  ++iparticle
258  )
259  {
260 
261  if ( (*iparticle)->numberOfTrackerLayers() )
262  {
263 
264  // Special treatment for decays
265  if (process == CMS::Decay)
266  {
267  // Get particle type
269  // Check if the particle type is valid one
270  if (particleID.isValid())
271  {
272  // Get particle data
273  ParticleData const * particleData = particleDataTable_->particle(particleID);
274  // Check if the particle exist in the table
275  if (particleData)
276  {
277  // Check if their life time is bigger than 1e-14
278  if ( particleDataTable_->particle(particleID)->lifetime() > longLivedDecayLength_ )
279  {
280  // Check for B, C weak decays and long lived decays
281  update(flags_[BWeakDecay], particleID.hasBottom());
282  update(flags_[CWeakDecay], particleID.hasCharm());
283  update(flags_[LongLivedDecay], true);
284  }
285  // Check Tau, Ks and Lambda decay
286  update(flags_[TauDecay], pdgid == 15);
287  update(flags_[KsDecay], pdgid == 310);
288  update(flags_[LambdaDecay], pdgid == 3122);
289  update(flags_[JpsiDecay], pdgid == 443);
290  update(flags_[XiDecay], pdgid == 3312);
291  update(flags_[OmegaDecay], pdgid == 3334);
292  update(flags_[SigmaPlusDecay], pdgid == 3222);
293  update(flags_[SigmaMinusDecay], pdgid == 3112);
294  }
295  }
296  }
297  }
298  }
299  }
300 }
301 
302 
304 {
305  // Helper class for clusterization
306  typedef std::multimap<double, HepMC::ThreeVector> Clusters;
307  typedef std::pair<double, HepMC::ThreeVector> ClusterPair;
308 
309  Clusters clusters;
310 
311  // Get the main primary vertex from the list
312  GeneratedPrimaryVertex const & genpv = genpvs_.back();
313 
314  // Get the generated history of the tracks
315  const VertexHistory::GenVertexTrail & genVertexTrail = tracer_.genVertexTrail();
316 
317  // Unit transformation from mm to cm
318  double const mm = 0.1;
319 
320  // Loop over the generated vertexes
321  for (
322  VertexHistory::GenVertexTrail::const_iterator ivertex = genVertexTrail.begin();
323  ivertex != genVertexTrail.end();
324  ++ivertex
325  )
326  {
327  // Check vertex exist
328  if (*ivertex)
329  {
330  // Measure the distance2 respecto the primary vertex
331  HepMC::ThreeVector p = (*ivertex)->point3d();
332  double distance = sqrt( pow(p.x() * mm - genpv.x, 2) + pow(p.y() * mm - genpv.y, 2) + pow(p.z() * mm - genpv.z, 2) );
333 
334  // If there is not any clusters add the first vertex.
335  if ( clusters.empty() )
336  {
337  clusters.insert( ClusterPair(distance, HepMC::ThreeVector(p.x() * mm, p.y() * mm, p.z() * mm)) );
338  continue;
339  }
340 
341  // Check if there is already a cluster in the given distance from primary vertex
342  Clusters::const_iterator icluster = clusters.lower_bound(distance - vertexClusteringDistance_);
343 
344  if ( icluster == clusters.upper_bound(distance + vertexClusteringDistance_) )
345  {
346  clusters.insert ( ClusterPair(distance, HepMC::ThreeVector(p.x() * mm, p.y() * mm, p.z() * mm)) );
347  continue;
348  }
349 
350  bool cluster = false;
351 
352  // Looping over the vertex clusters of a given distance from primary vertex
353  for (;
354  icluster != clusters.upper_bound(distance + vertexClusteringDistance_);
355  ++icluster
356  )
357  {
358  double difference = sqrt (
359  pow(p.x() * mm - icluster->second.x(), 2) +
360  pow(p.y() * mm - icluster->second.y(), 2) +
361  pow(p.z() * mm - icluster->second.z(), 2)
362  );
363 
364  if ( difference < vertexClusteringDistance_ )
365  {
366  cluster = true;
367  break;
368  }
369  }
370 
371  if (!cluster) clusters.insert ( ClusterPair(distance, HepMC::ThreeVector(p.x() * mm, p.y() * mm, p.z() * mm)) );
372  }
373  }
374 
375  const VertexHistory::SimVertexTrail & simVertexTrail = tracer_.simVertexTrail();
376 
377  // Loop over the generated particles
378  for (
379  VertexHistory::SimVertexTrail::const_reverse_iterator ivertex = simVertexTrail.rbegin();
380  ivertex != simVertexTrail.rend();
381  ++ivertex
382  )
383  {
384  // Look for those with production vertex
385  TrackingVertex::LorentzVector p = (*ivertex)->position();
386 
387  double distance = sqrt( pow(p.x() - genpv.x, 2) + pow(p.y() - genpv.y, 2) + pow(p.z() - genpv.z, 2) );
388 
389  // If there is not any clusters add the first vertex.
390  if ( clusters.empty() )
391  {
392  clusters.insert( ClusterPair(distance, HepMC::ThreeVector(p.x(), p.y(), p.z())) );
393  continue;
394  }
395 
396  // Check if there is already a cluster in the given distance from primary vertex
397  Clusters::const_iterator icluster = clusters.lower_bound(distance - vertexClusteringDistance_);
398 
399  if ( icluster == clusters.upper_bound(distance + vertexClusteringDistance_) )
400  {
401  clusters.insert ( ClusterPair(distance, HepMC::ThreeVector(p.x(), p.y(), p.z())) );
402  continue;
403  }
404 
405  bool cluster = false;
406 
407  // Looping over the vertex clusters of a given distance from primary vertex
408  for (;
409  icluster != clusters.upper_bound(distance + vertexClusteringDistance_);
410  ++icluster
411  )
412  {
413  double difference = sqrt (
414  pow(p.x() - icluster->second.x(), 2) +
415  pow(p.y() - icluster->second.y(), 2) +
416  pow(p.z() - icluster->second.z(), 2)
417  );
418 
419  if ( difference < vertexClusteringDistance_ )
420  {
421  cluster = true;
422  break;
423  }
424  }
425 
426  if (!cluster) clusters.insert ( ClusterPair(distance, HepMC::ThreeVector(p.x(), p.y(), p.z())) );
427  }
428 
429  if ( clusters.size() == 1 )
430  flags_[PrimaryVertex] = true;
431  else if ( clusters.size() == 2 )
432  flags_[SecondaryVertex] = true;
433  else
434  flags_[TertiaryVertex] = true;
435 }
436 
437 
439 {
440  return !p->end_vertex() && p->status() == 1;
441 }
442 
443 
445 {
446  const ParticleData * part = particleDataTable_->particle( p->pdg_id() );
447  if (part)
448  return part->charge()!=0;
449  else
450  {
451  // the new/improved particle table doesn't know anti-particles
452  return particleDataTable_->particle( -p->pdg_id() ) != nullptr;
453  }
454 }
455 
456 
458 {
459  genpvs_.clear();
460 
461  const HepMC::GenEvent * event = mcInformation_->GetEvent();
462 
463  if (event)
464  {
465  int idx = 0;
466 
467  // Loop over the different GenVertex
468  for ( HepMC::GenEvent::vertex_const_iterator ivertex = event->vertices_begin(); ivertex != event->vertices_end(); ++ivertex )
469  {
470  bool hasParentVertex = false;
471 
472  // Loop over the parents looking to see if they are coming from a production vertex
473  for (
474  HepMC::GenVertex::particle_iterator iparent = (*ivertex)->particles_begin(HepMC::parents);
475  iparent != (*ivertex)->particles_end(HepMC::parents);
476  ++iparent
477  )
478  if ( (*iparent)->production_vertex() )
479  {
480  hasParentVertex = true;
481  break;
482  }
483 
484  // Reject those vertices with parent vertices
485  if (hasParentVertex) continue;
486 
487  // Get the position of the vertex
488  HepMC::FourVector pos = (*ivertex)->position();
489 
490  double const mm = 0.1;
491 
492  GeneratedPrimaryVertex pv(pos.x()*mm, pos.y()*mm, pos.z()*mm);
493 
494  std::vector<GeneratedPrimaryVertex>::iterator ientry = genpvs_.begin();
495 
496  // Search for a VERY close vertex in the list
497  for (; ientry != genpvs_.end(); ++ientry)
498  {
499  double distance = sqrt( pow(pv.x - ientry->x, 2) + pow(pv.y - ientry->y, 2) + pow(pv.z - ientry->z, 2) );
500  if ( distance < vertexClusteringDistance_ )
501  break;
502  }
503 
504  // Check if there is not a VERY close vertex and added to the list
505  if (ientry == genpvs_.end())
506  ientry = genpvs_.insert(ientry,pv);
507 
508  // Add the vertex barcodes to the new or existent vertices
509  ientry->genVertex.push_back((*ivertex)->barcode());
510 
511  // Collect final state descendants
512  for (
513  HepMC::GenVertex::particle_iterator idecendants = (*ivertex)->particles_begin(HepMC::descendants);
514  idecendants != (*ivertex)->particles_end(HepMC::descendants);
515  ++idecendants
516  )
517  {
518  if (isFinalstateParticle(*idecendants))
519  if ( find(ientry->finalstateParticles.begin(), ientry->finalstateParticles.end(), (*idecendants)->barcode()) == ientry->finalstateParticles.end() )
520  {
521  ientry->finalstateParticles.push_back((*idecendants)->barcode());
522  HepMC::FourVector m = (*idecendants)->momentum();
523 
524  ientry->ptot.setPx(ientry->ptot.px() + m.px());
525  ientry->ptot.setPy(ientry->ptot.py() + m.py());
526  ientry->ptot.setPz(ientry->ptot.pz() + m.pz());
527  ientry->ptot.setE(ientry->ptot.e() + m.e());
528  ientry->ptsq += m.perp() * m.perp();
529 
530  if ( m.perp() > 0.8 && std::abs(m.pseudoRapidity()) < 2.5 && isCharged(*idecendants) ) ientry->nGenTrk++;
531  }
532  }
533  idx++;
534  }
535  }
536 
537  std::sort(genpvs_.begin(), genpvs_.end());
538 }
edm::Handle< edm::HepMCProduct > mcInformation_
T getUntrackedParameter(std::string const &, T const &) const
VertexHistory tracer_
int event() const
get the contents of the subdetector field (should be protected?)
TPRegexp parents
Definition: eve_filter.cc:21
void newEvent(const edm::Event &, const edm::EventSetup &)
Pre-process event information (for accessing reconstruction information)
void vertexInformation()
Get geometrical information about the vertices.
std::vector< GeneratedPrimaryVertex > genpvs_
const edm::InputTag hepMCLabel_
SimVertexTrail const & simVertexTrail() const
Return all the simulated vertices in the history.
Definition: HistoryBase.h:59
bool isCharged(const HepMC::GenParticle *)
def setup(process, global_tag, zero_tesla=False)
Definition: GeneralSetup.py:1
Get track history and classify it in function of their .
void reset()
Reset the categories flags.
bool isNonnull() const
Checks for non-null.
Definition: RefToBase.h:337
Definition: config.py:1
const reco::VertexBaseRef & recoVertex() const
Return a reference to the reconstructed track.
Definition: VertexHistory.h:64
bool isFinalstateParticle(const HepMC::GenParticle *)
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
void getData(T &iHolder) const
Definition: EventSetup.h:78
math::XYZTLorentzVectorD LorentzVector
bool evaluate(TrackingVertexRef tvr)
Evaluate track history using a TrackingParticleRef.
Definition: VertexHistory.h:42
Flags flags_
Flag containers.
T sqrt(T t)
Definition: SSEVec.h:18
int bunchCrossing() const
get the detector field from this detid
double longLivedDecayLength_
def pv(vc)
Definition: MetAnalyzer.py:6
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< const HepMC::GenVertex * > GenVertexTrail
GenVertex trail type.
Definition: HistoryBase.h:27
edm::ESHandle< ParticleDataTable > particleDataTable_
const G4toCMSLegacyProcTypeMap g4toCMSProcMap_
const unsigned int processId(unsigned int g4ProcessId) const
Definition: Utils.cc:59
HepPDT::ParticleData ParticleData
const TrackingVertexRef & simVertex() const
Return the initial tracking vertex from the history.
Definition: HistoryBase.h:95
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:38
VertexClassifier(edm::ParameterSet const &pset, edm::ConsumesCollector &&)
Constructor by ParameterSet.
part
Definition: HCALResponse.h:20
void processesAtGenerator()
Get all the information related to decay process.
VertexClassifier const & evaluate(reco::VertexBaseRef const &)
Classify the RecoVertex in categories.
HLT enums.
void depth(int d)
Set the depth of the history.
Definition: HistoryBase.h:53
double vertexClusteringDistance_
virtual void newEvent(edm::Event const &, edm::EventSetup const &)
Pre-process event information (for accessing reconstraction information)
void processesAtSimulation()
Get information about conversion and other interactions.
GenVertexTrail const & genVertexTrail() const
Return all generated vertex in the history.
Definition: HistoryBase.h:71
#define update(a, b)
std::vector< TrackingVertexRef > SimVertexTrail
SimVertex trail type.
Definition: HistoryBase.h:36
Auxiliary class holding simulated primary vertices.
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
def move(src, dest)
Definition: eostools.py:510
Definition: event.py:1
void simulationInformation()
Get all the information related to the simulation details.