CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Types | Public Member Functions | Protected 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
 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
< TrackingParticleRef
SimParticleTrail
 SimParticle trail type. More...
 
typedef std::vector
< TrackingVertexRef
SimVertexTrail
 SimVertex trail type. More...
 

Public Member Functions

void depth (int d)
 Set the depth of the history. More...
 
const HepMC::GenParticle * genParticle () const
 Returns a pointer to most primitive status 1 or 2 particle. More...
 
GenParticleTrail const & genParticleTrail () const
 Return all generated particle in the history. More...
 
GenVertexTrail const & genVertexTrail () const
 Return all generated vertex in the history. More...
 
 HistoryBase ()
 
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 Member Functions

bool evaluate (TrackingParticleRef tpr)
 Evaluate track history using a TrackingParticleRef. More...
 
bool evaluate (TrackingVertexRef tvr)
 Evaluate track history using a TrackingParticleRef. More...
 

Protected Attributes

GenParticleTrail genParticleTrail_
 
GenVertexTrail genVertexTrail_
 
GenVertexTrailHelper genVertexTrailHelper_
 
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 GenParticle. More...
 
void traceGenHistory (HepMC::GenVertex const *)
 Trace all the simulated information for a given pointer to a GenVertex. 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

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 21 of file HistoryBase.h.

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

GenVertex trail helper type.

Definition at line 24 of file HistoryBase.h.

SimParticle trail type.

Definition at line 27 of file HistoryBase.h.

SimVertex trail type.

Definition at line 30 of file HistoryBase.h.

Constructor & Destructor Documentation

HistoryBase::HistoryBase ( )
inline

Definition at line 33 of file HistoryBase.h.

References depth_.

34  {
35  // Default depth
36  depth_ = -1;
37  }

Member Function Documentation

void HistoryBase::depth ( int  d)
inline

Set the depth of the history.

Definition at line 47 of file HistoryBase.h.

References ztail::d, and depth_.

Referenced by traceSimHistory(), TrackClassifier::TrackClassifier(), and VertexClassifier::VertexClassifier().

48  {
49  depth_ = d;
50  }
tuple d
Definition: ztail.py:151
bool HistoryBase::evaluate ( TrackingParticleRef  tpr)
inlineprotected

Evaluate track history using a TrackingParticleRef.

Definition at line 114 of file HistoryBase.h.

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

Referenced by TrackHistory::evaluate(), and VertexHistory::evaluate().

115  {
116  resetTrails(tpr);
117  return traceSimHistory(tpr, depth_);
118  }
void resetTrails()
Reset trail functions.
Definition: HistoryBase.h:151
bool traceSimHistory(TrackingParticleRef const &, int)
Trace all the simulated information for a given reference to a TrackingParticle.
Definition: HistoryBase.cc:39
bool HistoryBase::evaluate ( TrackingVertexRef  tvr)
inlineprotected

Evaluate track history using a TrackingParticleRef.

Definition at line 128 of file HistoryBase.h.

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

129  {
130  resetTrails();
131  return traceSimHistory(tvr, depth_);
132  }
void resetTrails()
Reset trail functions.
Definition: HistoryBase.h:151
bool traceSimHistory(TrackingParticleRef const &, int)
Trace all the simulated information for a given reference to a TrackingParticle.
Definition: HistoryBase.cc:39
const HepMC::GenParticle* HistoryBase::genParticle ( ) const
inline

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

Definition at line 89 of file HistoryBase.h.

References genParticleTrail_.

Referenced by TrackClassifier::hadronFlavor(), and GenTrackMatcher::produce().

90  {
91  if ( genParticleTrail_.empty() ) return 0;
92  return genParticleTrail_[genParticleTrail_.size()-1];
93  }
GenParticleTrail genParticleTrail_
Definition: HistoryBase.h:99
GenParticleTrail const& HistoryBase::genParticleTrail ( ) const
inline

Return all generated particle in the history.

Definition at line 71 of file HistoryBase.h.

References genParticleTrail_.

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

72  {
73  return genParticleTrail_;
74  }
GenParticleTrail genParticleTrail_
Definition: HistoryBase.h:99
GenVertexTrail const& HistoryBase::genVertexTrail ( ) const
inline

Return all generated vertex in the history.

Definition at line 65 of file HistoryBase.h.

References genVertexTrail_.

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

66  {
67  return genVertexTrail_;
68  }
GenVertexTrail genVertexTrail_
Definition: HistoryBase.h:98
void HistoryBase::resetTrails ( )
inlineprivate

Reset trail functions.

Definition at line 151 of file HistoryBase.h.

References genParticleTrail_, genVertexTrail_, genVertexTrailHelper_, simParticleTrail_, and simVertexTrail_.

Referenced by evaluate(), and resetTrails().

152  {
153  simParticleTrail_.clear();
154  simVertexTrail_.clear();
155  genVertexTrail_.clear();
156  genParticleTrail_.clear();
157  genVertexTrailHelper_.clear();
158  }
SimVertexTrail simVertexTrail_
Definition: HistoryBase.h:100
GenVertexTrail genVertexTrail_
Definition: HistoryBase.h:98
GenVertexTrailHelper genVertexTrailHelper_
Definition: HistoryBase.h:104
SimParticleTrail simParticleTrail_
Definition: HistoryBase.h:101
GenParticleTrail genParticleTrail_
Definition: HistoryBase.h:99
void HistoryBase::resetTrails ( TrackingParticleRef  tpr)
inlineprivate

Definition at line 160 of file HistoryBase.h.

References resetTrails(), and simParticleTrail_.

161  {
162  resetTrails();
163  simParticleTrail_.push_back(tpr);
164  }
void resetTrails()
Reset trail functions.
Definition: HistoryBase.h:151
SimParticleTrail simParticleTrail_
Definition: HistoryBase.h:101
const TrackingParticleRef& HistoryBase::simParticle ( ) const
inline

Return the initial tracking particle from the history.

Definition at line 77 of file HistoryBase.h.

References simParticleTrail_.

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

78  {
79  return simParticleTrail_[0];
80  }
SimParticleTrail simParticleTrail_
Definition: HistoryBase.h:101
SimParticleTrail const& HistoryBase::simParticleTrail ( ) const
inline

Return all the simulated particle in the history.

Definition at line 59 of file HistoryBase.h.

References simParticleTrail_.

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

60  {
61  return simParticleTrail_;
62  }
SimParticleTrail simParticleTrail_
Definition: HistoryBase.h:101
const TrackingVertexRef& HistoryBase::simVertex ( ) const
inline

Return the initial tracking vertex from the history.

Definition at line 83 of file HistoryBase.h.

References simVertexTrail_.

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

84  {
85  return simVertexTrail_[0];
86  }
SimVertexTrail simVertexTrail_
Definition: HistoryBase.h:100
SimVertexTrail const& HistoryBase::simVertexTrail ( ) const
inline

Return all the simulated vertices in the history.

Definition at line 53 of file HistoryBase.h.

References simVertexTrail_.

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

54  {
55  return simVertexTrail_;
56  }
SimVertexTrail simVertexTrail_
Definition: HistoryBase.h:100
void HistoryBase::traceGenHistory ( HepMC::GenParticle const *  genParticle)
private

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

Definition at line 8 of file HistoryBase.cc.

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

Referenced by traceGenHistory(), and traceSimHistory().

9 {
10  // Third stop criteria: status abs(depth_) particles after the hadronization.
11  // The after hadronization is done by detecting the pdg_id pythia code from 88 to 99
12  if ( genParticle->status() <= abs(depth_) && (genParticle->pdg_id() < 88 || genParticle->pdg_id() > 99) )
13  {
14  genParticleTrail_.push_back(genParticle);
15  // Get the producer vertex and trace it history
16  traceGenHistory( genParticle->production_vertex() );
17  }
18 }
void traceGenHistory(HepMC::GenParticle const *)
Trace all the simulated information for a given pointer to a GenParticle.
Definition: HistoryBase.cc:8
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
GenParticleTrail genParticleTrail_
Definition: HistoryBase.h:99
const HepMC::GenParticle * genParticle() const
Returns a pointer to most primitive status 1 or 2 particle.
Definition: HistoryBase.h:89
void HistoryBase::traceGenHistory ( HepMC::GenVertex const *  genVertex)
private

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

Definition at line 21 of file HistoryBase.cc.

References genVertexTrail_, genVertexTrailHelper_, and traceGenHistory().

22 {
23  // Verify if has a vertex associated
24  if (genVertex)
25  {
26  // Skip if already exist in the collection
27  if ( genVertexTrailHelper_.find(genVertex) != genVertexTrailHelper_.end() )
28  return;
29  // Add vertex to the history
30  genVertexTrail_.push_back(genVertex);
31  genVertexTrailHelper_.insert(genVertex);
32  // Verify if the vertex has incoming particles
33  if ( genVertex->particles_in_size() )
34  traceGenHistory( *(genVertex->particles_in_const_begin()) );
35  }
36 }
GenVertexTrail genVertexTrail_
Definition: HistoryBase.h:98
void traceGenHistory(HepMC::GenParticle const *)
Trace all the simulated information for a given pointer to a GenParticle.
Definition: HistoryBase.cc:8
GenVertexTrailHelper genVertexTrailHelper_
Definition: HistoryBase.h:104
bool HistoryBase::traceSimHistory ( TrackingParticleRef const &  trackingParticle,
int  depth 
)
private

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

Definition at line 39 of file HistoryBase.cc.

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

Referenced by evaluate(), and traceSimHistory().

40 {
41  // first stop condition: if the required depth is reached
42  if ( depth == depth_ && depth_ >= 0 ) return true;
43 
44  // second stop condition: if a gen particle is associated to the TP
45  if ( !trackingParticle->genParticles().empty() )
46  {
47  LogDebug("TrackHistory") << "Particle " << trackingParticle->pdgId() << " has a GenParicle image."
48  << std::endl;
49  //#warning "This file has been modified just to get it to compile without any regard as to whether it still functions as intended"
50  // this code does not compile likely due to the wrong type of iterator
51  #ifdef REMOVED_JUST_TO_GET_IT_TO_COMPILE__THIS_CODE_NEEDS_TO_BE_CHECKED
52  traceGenHistory(&(**(trackingParticle->genParticle_begin())));
53  #endif
54 
55  }
56 
57  LogDebug("TrackHistory") << "No GenParticle image for " << trackingParticle->pdgId() << std::endl;
58 
59  // get a reference to the TP's parent vertex and trace it history
60  return traceSimHistory( trackingParticle->parentVertex(), depth );
61 }
#define LogDebug(id)
void traceGenHistory(HepMC::GenParticle const *)
Trace all the simulated information for a given pointer to a GenParticle.
Definition: HistoryBase.cc:8
void depth(int d)
Set the depth of the history.
Definition: HistoryBase.h:47
bool traceSimHistory(TrackingParticleRef const &, int)
Trace all the simulated information for a given reference to a TrackingParticle.
Definition: HistoryBase.cc:39
bool HistoryBase::traceSimHistory ( TrackingVertexRef const &  trackingVertex,
int  depth 
)
private

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

Definition at line 64 of file HistoryBase.cc.

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

65 {
66  // verify if the parent vertex exists
67  if ( trackingVertex.isNonnull() )
68  {
69  // save the vertex in the trail
70  simVertexTrail_.push_back(trackingVertex);
71 
72  if ( !trackingVertex->sourceTracks().empty() )
73  {
74  LogDebug("TrackHistory") << "Moving on to the parent particle." << std::endl;
75 
76  // select the original source in case of combined vertices
77  bool flag = false;
79 
80  for (its = trackingVertex->sourceTracks_begin(); its != trackingVertex->sourceTracks_end(); its++)
81  {
82  for (itd = trackingVertex->daughterTracks_begin(); itd != trackingVertex->daughterTracks_end(); itd++)
83  if (itd != its)
84  {
85  flag = true;
86  break;
87  }
88  if (flag)
89  break;
90  }
91 
92  // verify if the new particle is not in the trail (looping partiles)
93  if (
94  std::find(
95  simParticleTrail_.begin(),
96  simParticleTrail_.end(),
97  *its
98  ) != simParticleTrail_.end()
99  )
100  {
101  LogDebug("TrackHistory") << "WARNING: Looping track found." << std::endl;
102  return false;
103  }
104 
105  // save particle in the trail
106  simParticleTrail_.push_back(*its);
107  return traceSimHistory (*its, --depth);
108  }
109  else if ( !trackingVertex->genVertices().empty() )
110  {
111  // navigate over all the associated generated vertexes
112  LogDebug("TrackHistory") << "Vertex has " << trackingVertex->genVertices().size() << "GenVertex image." << std::endl;
113  for (
114  TrackingVertex::genv_iterator ivertex = trackingVertex->genVertices_begin();
115  ivertex != trackingVertex->genVertices_end();
116  ++ivertex
117  )
118  traceGenHistory(&(**(ivertex)));
119  return true;
120  }
121  else
122  {
123  LogDebug("TrackHistory") << "WARNING: Source track for tracking vertex cannot be found." << std::endl;
124  }
125  }
126  else
127  {
128  LogDebug("TrackHistory") << " WARNING: Vertex cannot be found.";
129  }
130 
131  return false;
132 }
#define LogDebug(id)
SimVertexTrail simVertexTrail_
Definition: HistoryBase.h:100
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
void traceGenHistory(HepMC::GenParticle const *)
Trace all the simulated information for a given pointer to a GenParticle.
Definition: HistoryBase.cc:8
SimParticleTrail simParticleTrail_
Definition: HistoryBase.h:101
void depth(int d)
Set the depth of the history.
Definition: HistoryBase.h:47
bool traceSimHistory(TrackingParticleRef const &, int)
Trace all the simulated information for a given reference to a TrackingParticle.
Definition: HistoryBase.cc:39

Member Data Documentation

int HistoryBase::depth_
private

Definition at line 136 of file HistoryBase.h.

Referenced by depth(), evaluate(), HistoryBase(), traceGenHistory(), and traceSimHistory().

GenParticleTrail HistoryBase::genParticleTrail_
protected

Definition at line 99 of file HistoryBase.h.

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

GenVertexTrail HistoryBase::genVertexTrail_
protected

Definition at line 98 of file HistoryBase.h.

Referenced by genVertexTrail(), resetTrails(), and traceGenHistory().

GenVertexTrailHelper HistoryBase::genVertexTrailHelper_
protected

Definition at line 104 of file HistoryBase.h.

Referenced by resetTrails(), and traceGenHistory().

SimParticleTrail HistoryBase::simParticleTrail_
protected

Definition at line 101 of file HistoryBase.h.

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

SimVertexTrail HistoryBase::simVertexTrail_
protected

Definition at line 100 of file HistoryBase.h.

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