CMS 3D CMS Logo

List of all members | Public Types | Public Member Functions | Protected Attributes | Private Member Functions | Private Attributes
HistoryBase Class Reference

Base class to all the history types. More...

#include <HistoryBase.h>

Inheritance diagram for HistoryBase:
TrackHistory VertexHistory

Public Types

typedef std::vector< const HepMC::GenParticle * > GenParticleTrail
 HepMC::GenParticle trail type. More...
 
typedef std::vector< const HepMC::GenVertex * > GenVertexTrail
 GenVertex trail type. More...
 
typedef std::set< const HepMC::GenVertex * > GenVertexTrailHelper
 GenVertex trail helper type. More...
 
typedef std::vector< const reco::GenParticle * > RecoGenParticleTrail
 reco::GenParticle trail type. More...
 
typedef std::set< const reco::GenParticle * > RecoGenParticleTrailHelper
 reco::GenParticle trail helper type. More...
 
typedef std::vector< TrackingParticleRefSimParticleTrail
 SimParticle trail type. More...
 
typedef std::vector< TrackingVertexRefSimVertexTrail
 SimVertex trail type. More...
 

Public Member Functions

void depth (int d)
 Set the depth of the history. More...
 
bool evaluate (TrackingParticleRef tpr)
 Evaluate track history using a TrackingParticleRef. More...
 
bool evaluate (TrackingVertexRef tvr)
 Evaluate track history using a TrackingParticleRef. More...
 
const HepMC::GenParticle * genParticle () const
 Returns a pointer to most primitive status 1 or 2 particle in the genParticleTrail_. More...
 
GenParticleTrail const & genParticleTrail () const
 Return all generated particle (HepMC::GenParticle) in the history. More...
 
GenVertexTrail const & genVertexTrail () const
 Return all generated vertex in the history. More...
 
 HistoryBase ()
 
const reco::GenParticlerecoGenParticle () const
 Returns a pointer to most primitive status 1 or 2 particle in the recoGenParticleTrail_. More...
 
RecoGenParticleTrail const & recoGenParticleTrail () const
 Return all reco::GenParticle in the history. More...
 
const TrackingParticleRefsimParticle () const
 Return the initial tracking particle from the history. More...
 
SimParticleTrail const & simParticleTrail () const
 Return all the simulated particle in the history. More...
 
const TrackingVertexRefsimVertex () const
 Return the initial tracking vertex from the history. More...
 
SimVertexTrail const & simVertexTrail () const
 Return all the simulated vertices in the history. More...
 

Protected Attributes

GenParticleTrail genParticleTrail_
 
GenVertexTrail genVertexTrail_
 
GenVertexTrailHelper genVertexTrailHelper_
 
RecoGenParticleTrail recoGenParticleTrail_
 
RecoGenParticleTrailHelper recoGenParticleTrailHelper_
 
SimParticleTrail simParticleTrail_
 
SimVertexTrail simVertexTrail_
 

Private Member Functions

void resetTrails ()
 Reset trail functions. More...
 
void resetTrails (TrackingParticleRef tpr)
 
void traceGenHistory (HepMC::GenParticle const *)
 Trace all the simulated information for a given pointer to a HepMC::GenParticle. More...
 
void traceGenHistory (HepMC::GenVertex const *)
 Trace all the simulated information for a given pointer to a GenVertex. More...
 
void traceRecoGenHistory (reco::GenParticle const *)
 Trace all the simulated information for a given pointer to a reco::GenParticle. More...
 
bool traceSimHistory (TrackingParticleRef const &, int)
 Trace all the simulated information for a given reference to a TrackingParticle. More...
 
bool traceSimHistory (TrackingVertexRef const &, int)
 Trace all the simulated information for a given reference to a TrackingVertex. More...
 

Private Attributes

int depth_
 

Detailed Description

Base class to all the history types.

Definition at line 12 of file HistoryBase.h.

Member Typedef Documentation

typedef std::vector<const HepMC::GenParticle *> HistoryBase::GenParticleTrail

HepMC::GenParticle trail type.

Definition at line 18 of file HistoryBase.h.

typedef std::vector<const HepMC::GenVertex *> HistoryBase::GenVertexTrail

GenVertex trail type.

Definition at line 27 of file HistoryBase.h.

typedef std::set<const HepMC::GenVertex *> HistoryBase::GenVertexTrailHelper

GenVertex trail helper type.

Definition at line 30 of file HistoryBase.h.

reco::GenParticle trail type.

Definition at line 21 of file HistoryBase.h.

reco::GenParticle trail helper type.

Definition at line 24 of file HistoryBase.h.

SimParticle trail type.

Definition at line 33 of file HistoryBase.h.

SimVertex trail type.

Definition at line 36 of file HistoryBase.h.

Constructor & Destructor Documentation

HistoryBase::HistoryBase ( )
inline

Definition at line 39 of file HistoryBase.h.

References depth_.

40  {
41  // Default depth
42  depth_ = -1;
43  }

Member Function Documentation

void HistoryBase::depth ( int  d)
inline
bool HistoryBase::evaluate ( TrackingParticleRef  tpr)
inline

Evaluate track history using a TrackingParticleRef.

Definition at line 122 of file HistoryBase.h.

References depth_, resetTrails(), and traceSimHistory().

Referenced by TrackHistory::evaluate(), VertexHistory::evaluate(), TrackingNtuple::fillTrackingParticles(), and TrackingParticleBHadronRefSelector::produce().

123  {
124  resetTrails(tpr);
125  return traceSimHistory(tpr, depth_);
126  }
void resetTrails()
Reset trail functions.
Definition: HistoryBase.h:176
bool traceSimHistory(TrackingParticleRef const &, int)
Trace all the simulated information for a given reference to a TrackingParticle.
Definition: HistoryBase.cc:62
bool HistoryBase::evaluate ( TrackingVertexRef  tvr)
inline

Evaluate track history using a TrackingParticleRef.

Definition at line 136 of file HistoryBase.h.

References depth_, resetTrails(), and traceSimHistory().

137  {
138  resetTrails();
139  return traceSimHistory(tvr, depth_);
140  }
void resetTrails()
Reset trail functions.
Definition: HistoryBase.h:176
bool traceSimHistory(TrackingParticleRef const &, int)
Trace all the simulated information for a given reference to a TrackingParticle.
Definition: HistoryBase.cc:62
const HepMC::GenParticle* HistoryBase::genParticle ( ) const
inline

Returns a pointer to most primitive status 1 or 2 particle in the genParticleTrail_.

Definition at line 101 of file HistoryBase.h.

References genParticleTrail_.

Referenced by GenTrackMatcher::produce().

102  {
103  if ( genParticleTrail_.empty() ) return nullptr;
104  return genParticleTrail_[genParticleTrail_.size()-1];
105  }
GenParticleTrail genParticleTrail_
Definition: HistoryBase.h:146
GenParticleTrail const& HistoryBase::genParticleTrail ( ) const
inline

Return all generated particle (HepMC::GenParticle) in the history.

Definition at line 77 of file HistoryBase.h.

References genParticleTrail_.

Referenced by VertexHistoryAnalyzer::analyze(), TrackHistoryAnalyzer::analyze(), and TrackClassifier::vertexInformation().

78  {
79  return genParticleTrail_;
80  }
GenParticleTrail genParticleTrail_
Definition: HistoryBase.h:146
GenVertexTrail const& HistoryBase::genVertexTrail ( ) const
inline

Return all generated vertex in the history.

Definition at line 71 of file HistoryBase.h.

References genVertexTrail_.

Referenced by VertexHistoryAnalyzer::analyze(), TrackHistoryAnalyzer::analyze(), VertexClassifier::processesAtGenerator(), and VertexClassifier::vertexInformation().

72  {
73  return genVertexTrail_;
74  }
GenVertexTrail genVertexTrail_
Definition: HistoryBase.h:145
const reco::GenParticle* HistoryBase::recoGenParticle ( ) const
inline

Returns a pointer to most primitive status 1 or 2 particle in the recoGenParticleTrail_.

Definition at line 108 of file HistoryBase.h.

References recoGenParticleTrail_.

Referenced by TrackClassifier::hadronFlavor().

109  {
110  if ( recoGenParticleTrail_.empty() ) return nullptr;
112  }
RecoGenParticleTrail recoGenParticleTrail_
Definition: HistoryBase.h:147
RecoGenParticleTrail const& HistoryBase::recoGenParticleTrail ( ) const
inline

Return all reco::GenParticle in the history.

Definition at line 83 of file HistoryBase.h.

References recoGenParticleTrail_.

Referenced by TrackingNtuple::fillTrackingParticles(), TrackClassifier::processesAtGenerator(), and TrackingParticleBHadronRefSelector::produce().

84  {
85  return recoGenParticleTrail_;
86  }
RecoGenParticleTrail recoGenParticleTrail_
Definition: HistoryBase.h:147
void HistoryBase::resetTrails ( )
inlineprivate

Reset trail functions.

Definition at line 176 of file HistoryBase.h.

Referenced by evaluate(), and resetTrails().

177  {
178  simParticleTrail_.clear();
179  simVertexTrail_.clear();
180  genVertexTrail_.clear();
181  genParticleTrail_.clear();
182  recoGenParticleTrail_.clear();
184  genVertexTrailHelper_.clear();
185  }
RecoGenParticleTrail recoGenParticleTrail_
Definition: HistoryBase.h:147
SimVertexTrail simVertexTrail_
Definition: HistoryBase.h:148
GenVertexTrail genVertexTrail_
Definition: HistoryBase.h:145
GenVertexTrailHelper genVertexTrailHelper_
Definition: HistoryBase.h:152
SimParticleTrail simParticleTrail_
Definition: HistoryBase.h:149
RecoGenParticleTrailHelper recoGenParticleTrailHelper_
Definition: HistoryBase.h:153
GenParticleTrail genParticleTrail_
Definition: HistoryBase.h:146
void HistoryBase::resetTrails ( TrackingParticleRef  tpr)
inlineprivate

Definition at line 187 of file HistoryBase.h.

References resetTrails().

188  {
189  resetTrails();
190  simParticleTrail_.push_back(tpr);
191  }
void resetTrails()
Reset trail functions.
Definition: HistoryBase.h:176
SimParticleTrail simParticleTrail_
Definition: HistoryBase.h:149
const TrackingParticleRef& HistoryBase::simParticle ( ) const
inline

Return the initial tracking particle from the history.

Definition at line 89 of file HistoryBase.h.

References simParticleTrail_.

Referenced by TrackClassifier::reconstructionInformation(), and TrackClassifier::simulationInformation().

90  {
91  return simParticleTrail_[0];
92  }
SimParticleTrail simParticleTrail_
Definition: HistoryBase.h:149
SimParticleTrail const& HistoryBase::simParticleTrail ( ) const
inline

Return all the simulated particle in the history.

Definition at line 65 of file HistoryBase.h.

References simParticleTrail_.

Referenced by VertexHistoryAnalyzer::analyze(), TrackHistoryAnalyzer::analyze(), TrackClassifier::processesAtSimulation(), TrackClassifier::qualityInformation(), and TrackClassifier::vertexInformation().

66  {
67  return simParticleTrail_;
68  }
SimParticleTrail simParticleTrail_
Definition: HistoryBase.h:149
const TrackingVertexRef& HistoryBase::simVertex ( ) const
inline

Return the initial tracking vertex from the history.

Definition at line 95 of file HistoryBase.h.

References simVertexTrail_.

Referenced by SVTagInfoValidationAnalyzer::analyze(), recoBSVTagInfoValidationAnalyzer::analyze(), and VertexClassifier::simulationInformation().

96  {
97  return simVertexTrail_[0];
98  }
SimVertexTrail simVertexTrail_
Definition: HistoryBase.h:148
SimVertexTrail const& HistoryBase::simVertexTrail ( ) const
inline

Return all the simulated vertices in the history.

Definition at line 59 of file HistoryBase.h.

References simVertexTrail_.

Referenced by VertexHistoryAnalyzer::analyze(), TrackHistoryAnalyzer::analyze(), VertexClassifier::processesAtSimulation(), and VertexClassifier::vertexInformation().

60  {
61  return simVertexTrail_;
62  }
SimVertexTrail simVertexTrail_
Definition: HistoryBase.h:148
void HistoryBase::traceGenHistory ( HepMC::GenParticle const *  genParticle)
private

Trace all the simulated information for a given pointer to a HepMC::GenParticle.

Definition at line 31 of file HistoryBase.cc.

References funct::abs(), depth_, and genParticleTrail_.

Referenced by traceGenHistory(), and traceSimHistory().

32 {
33  // Third stop criteria: status abs(depth_) particles after the hadronization.
34  // The after hadronization is done by detecting the pdg_id pythia code from 88 to 99
35  if ( genParticle->status() <= abs(depth_) && (genParticle->pdg_id() < 88 || genParticle->pdg_id() > 99) )
36  {
37  genParticleTrail_.push_back(genParticle);
38  // Get the producer vertex and trace it history
39  traceGenHistory( genParticle->production_vertex() );
40  }
41 }
void traceGenHistory(HepMC::GenParticle const *)
Trace all the simulated information for a given pointer to a HepMC::GenParticle.
Definition: HistoryBase.cc:31
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
GenParticleTrail genParticleTrail_
Definition: HistoryBase.h:146
const HepMC::GenParticle * genParticle() const
Returns a pointer to most primitive status 1 or 2 particle in the genParticleTrail_.
Definition: HistoryBase.h:101
void HistoryBase::traceGenHistory ( HepMC::GenVertex const *  genVertex)
private

Trace all the simulated information for a given pointer to a GenVertex.

Definition at line 44 of file HistoryBase.cc.

References genVertexTrail_, genVertexTrailHelper_, and traceGenHistory().

45 {
46  // Verify if has a vertex associated
47  if (genVertex)
48  {
49  // Skip if already exist in the collection
50  if ( genVertexTrailHelper_.find(genVertex) != genVertexTrailHelper_.end() )
51  return;
52  // Add vertex to the history
53  genVertexTrail_.push_back(genVertex);
54  genVertexTrailHelper_.insert(genVertex);
55  // Verify if the vertex has incoming particles
56  if ( genVertex->particles_in_size() )
57  traceGenHistory( *(genVertex->particles_in_const_begin()) );
58  }
59 }
GenVertexTrail genVertexTrail_
Definition: HistoryBase.h:145
void traceGenHistory(HepMC::GenParticle const *)
Trace all the simulated information for a given pointer to a HepMC::GenParticle.
Definition: HistoryBase.cc:31
GenVertexTrailHelper genVertexTrailHelper_
Definition: HistoryBase.h:152
void HistoryBase::traceRecoGenHistory ( reco::GenParticle const *  genParticle)
private

Trace all the simulated information for a given pointer to a reco::GenParticle.

Definition at line 8 of file HistoryBase.cc.

References funct::abs(), depth_, reco::CompositeRefCandidateT< D >::mother(), reco::LeafCandidate::pdgId(), recoGenParticleTrail_, recoGenParticleTrailHelper_, and reco::LeafCandidate::status().

Referenced by traceSimHistory().

9 {
10  // fill the trace of reco::GenParticle for correct reading of the History flags etc.
11  // This is called from the TrackingParticle->genParticle_begin() which is a reco::GenParticle
12 
13  // Take onlt genParticles with a status() smaller than "depth_". Typically depth is 2 and so you trace back untill you reach the hard partons
14  // see for status(): https://twiki.cern.ch/twiki/bin/view/CMSPublic/WorkBookGenParticleCandidate#GenPCand
15  if ( genParticle->status() <= abs(depth_) && (genParticle->pdgId() < 88 || genParticle->pdgId() > 99) )
16  {
17  //If the particle is already in the history, it is looping and you should stop
19  return;
20  }
23  // Get the genParticle's mother and trace its history
24  if (genParticle->mother() != nullptr){
26  }
27  }
28 }
RecoGenParticleTrail recoGenParticleTrail_
Definition: HistoryBase.h:147
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void traceRecoGenHistory(reco::GenParticle const *)
Trace all the simulated information for a given pointer to a reco::GenParticle.
Definition: HistoryBase.cc:8
RecoGenParticleTrailHelper recoGenParticleTrailHelper_
Definition: HistoryBase.h:153
const HepMC::GenParticle * genParticle() const
Returns a pointer to most primitive status 1 or 2 particle in the genParticleTrail_.
Definition: HistoryBase.h:101
bool HistoryBase::traceSimHistory ( TrackingParticleRef const &  trackingParticle,
int  depth 
)
private

Trace all the simulated information for a given reference to a TrackingParticle.

Definition at line 62 of file HistoryBase.cc.

References depth(), depth_, LogDebug, and traceRecoGenHistory().

Referenced by evaluate(), and traceSimHistory().

63 {
64  // first stop condition: if the required depth is reached
65  if ( depth == depth_ && depth_ >= 0 ) return true;
66 
67  // second stop condition: if a gen particle is associated to the TP
68  if ( !trackingParticle->genParticles().empty() )
69  {
70  LogDebug("TrackHistory") << "Particle " << trackingParticle->pdgId() << " has a GenParicle image."
71  << std::endl;
72 
73  traceRecoGenHistory(&(**(trackingParticle->genParticle_begin())));
74 
75 
76  }
77 
78  LogDebug("TrackHistory") << "No GenParticle image for " << trackingParticle->pdgId() << std::endl;
79 
80  // if no association between TP and genParticles is found: get a reference to the TP's parent vertex and trace its history
81  return traceSimHistory( trackingParticle->parentVertex(), depth );
82 }
#define LogDebug(id)
void traceRecoGenHistory(reco::GenParticle const *)
Trace all the simulated information for a given pointer to a reco::GenParticle.
Definition: HistoryBase.cc:8
void depth(int d)
Set the depth of the history.
Definition: HistoryBase.h:53
bool traceSimHistory(TrackingParticleRef const &, int)
Trace all the simulated information for a given reference to a TrackingParticle.
Definition: HistoryBase.cc:62
bool HistoryBase::traceSimHistory ( TrackingVertexRef const &  trackingVertex,
int  depth 
)
private

Trace all the simulated information for a given reference to a TrackingVertex.

Definition at line 85 of file HistoryBase.cc.

References spr::find(), RemoveAddSevLevel::flag, edm::Ref< C, T, F >::isNonnull(), LogDebug, simParticleTrail_, simVertexTrail_, traceGenHistory(), and traceSimHistory().

86 {
87  // verify if the parent vertex exists
88  if ( trackingVertex.isNonnull() )
89  {
90  // save the vertex in the trail
91  simVertexTrail_.push_back(trackingVertex);
92 
93  if ( !trackingVertex->sourceTracks().empty() )
94  {
95  LogDebug("TrackHistory") << "Moving on to the parent particle." << std::endl;
96 
97  // select the original source in case of combined vertices
98  bool flag = false;
100 
101  for (its = trackingVertex->sourceTracks_begin(); its != trackingVertex->sourceTracks_end(); its++)
102  {
103  for (itd = trackingVertex->daughterTracks_begin(); itd != trackingVertex->daughterTracks_end(); itd++)
104  if (itd != its)
105  {
106  flag = true;
107  break;
108  }
109  if (flag)
110  break;
111  }
112 
113  if(!flag)
114  return false;
115 
116  // verify if the new particle is not in the trail (looping partiles)
117  if (
118  std::find(
119  simParticleTrail_.begin(),
120  simParticleTrail_.end(),
121  *its
122  ) != simParticleTrail_.end()
123  )
124  {
125  LogDebug("TrackHistory") << "WARNING: Looping track found." << std::endl;
126  return false;
127  }
128 
129  // save particle in the trail
130  simParticleTrail_.push_back(*its);
131  return traceSimHistory (*its, --depth);
132  }
133  else if ( !trackingVertex->genVertices().empty() )
134  {
135  // navigate over all the associated generated vertexes
136  LogDebug("TrackHistory") << "Vertex has " << trackingVertex->genVertices().size() << "GenVertex image." << std::endl;
137  for (
138  TrackingVertex::genv_iterator ivertex = trackingVertex->genVertices_begin();
139  ivertex != trackingVertex->genVertices_end();
140  ++ivertex
141  )
142  traceGenHistory(&(**(ivertex)));
143  return true;
144  }
145  else
146  {
147  LogDebug("TrackHistory") << "WARNING: Source track for tracking vertex cannot be found." << std::endl;
148  }
149  }
150  else
151  {
152  LogDebug("TrackHistory") << " WARNING: Vertex cannot be found.";
153  }
154 
155  return false;
156 }
#define LogDebug(id)
SimVertexTrail simVertexTrail_
Definition: HistoryBase.h:148
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
void traceGenHistory(HepMC::GenParticle const *)
Trace all the simulated information for a given pointer to a HepMC::GenParticle.
Definition: HistoryBase.cc:31
SimParticleTrail simParticleTrail_
Definition: HistoryBase.h:149
void depth(int d)
Set the depth of the history.
Definition: HistoryBase.h:53
bool traceSimHistory(TrackingParticleRef const &, int)
Trace all the simulated information for a given reference to a TrackingParticle.
Definition: HistoryBase.cc:62

Member Data Documentation

int HistoryBase::depth_
private
GenParticleTrail HistoryBase::genParticleTrail_
protected

Definition at line 146 of file HistoryBase.h.

Referenced by genParticle(), genParticleTrail(), and traceGenHistory().

GenVertexTrail HistoryBase::genVertexTrail_
protected

Definition at line 145 of file HistoryBase.h.

Referenced by genVertexTrail(), and traceGenHistory().

GenVertexTrailHelper HistoryBase::genVertexTrailHelper_
protected

Definition at line 152 of file HistoryBase.h.

Referenced by traceGenHistory().

RecoGenParticleTrail HistoryBase::recoGenParticleTrail_
protected

Definition at line 147 of file HistoryBase.h.

Referenced by recoGenParticle(), recoGenParticleTrail(), and traceRecoGenHistory().

RecoGenParticleTrailHelper HistoryBase::recoGenParticleTrailHelper_
protected

Definition at line 153 of file HistoryBase.h.

Referenced by traceRecoGenHistory().

SimParticleTrail HistoryBase::simParticleTrail_
protected

Definition at line 149 of file HistoryBase.h.

Referenced by simParticle(), simParticleTrail(), and traceSimHistory().

SimVertexTrail HistoryBase::simVertexTrail_
protected

Definition at line 148 of file HistoryBase.h.

Referenced by simVertex(), simVertexTrail(), and traceSimHistory().