CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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
< 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...
 
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
 
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
 
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 *)
 
void traceGenHistory (HepMC::GenVertex const *)
 Trace all the simulated information for a given pointer to a GenVertex. More...
 
void traceRecoGenHistory (reco::GenParticle const *)
 
bool traceSimHistory (TrackingParticleRef const &, int)
 
bool traceSimHistory (TrackingVertexRef const &, int)
 

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

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

GenVertex trail type.

Definition at line 24 of file HistoryBase.h.

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

GenVertex trail helper type.

Definition at line 27 of file HistoryBase.h.

reco::GenParticle trail type.

Definition at line 18 of file HistoryBase.h.

reco::GenParticle trail helper type.

Definition at line 21 of file HistoryBase.h.

SimParticle trail type.

Definition at line 30 of file HistoryBase.h.

SimVertex trail type.

Definition at line 33 of file HistoryBase.h.

Constructor & Destructor Documentation

HistoryBase::HistoryBase ( )
inline

Definition at line 36 of file HistoryBase.h.

References depth_.

36  {
37  // Default depth
38  depth_ = -1;
39  }

Member Function Documentation

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

Evaluate track history using a TrackingParticleRef.

Definition at line 96 of file HistoryBase.h.

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

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

96  {
97  resetTrails(tpr);
98  return traceSimHistory(tpr, depth_);
99  }
void resetTrails()
Reset trail functions.
Definition: HistoryBase.h:149
bool traceSimHistory(TrackingParticleRef const &, int)
Definition: HistoryBase.cc:57
bool HistoryBase::evaluate ( TrackingVertexRef  tvr)
inline

Evaluate track history using a TrackingParticleRef.

Definition at line 109 of file HistoryBase.h.

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

109  {
110  resetTrails();
111  return traceSimHistory(tvr, depth_);
112  }
void resetTrails()
Reset trail functions.
Definition: HistoryBase.h:149
bool traceSimHistory(TrackingParticleRef const &, int)
Definition: HistoryBase.cc:57
const HepMC::GenParticle* HistoryBase::genParticle ( ) const
inline

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

Definition at line 74 of file HistoryBase.h.

References genParticleTrail_.

Referenced by GenTrackMatcher::produce().

74  {
75  if (genParticleTrail_.empty())
76  return nullptr;
77  return genParticleTrail_[genParticleTrail_.size() - 1];
78  }
GenParticleTrail genParticleTrail_
Definition: HistoryBase.h:117
GenParticleTrail const& HistoryBase::genParticleTrail ( ) const
inline

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

Definition at line 61 of file HistoryBase.h.

References genParticleTrail_.

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

61 { return genParticleTrail_; }
GenParticleTrail genParticleTrail_
Definition: HistoryBase.h:117
GenVertexTrail const& HistoryBase::genVertexTrail ( ) const
inline

Return all generated vertex in the history.

Definition at line 58 of file HistoryBase.h.

References genVertexTrail_.

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

58 { return genVertexTrail_; }
GenVertexTrail genVertexTrail_
Definition: HistoryBase.h:116
const reco::GenParticle* HistoryBase::recoGenParticle ( ) const
inline

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

Definition at line 82 of file HistoryBase.h.

References recoGenParticleTrail_.

Referenced by TrackClassifier::hadronFlavor().

82  {
83  if (recoGenParticleTrail_.empty())
84  return nullptr;
86  }
RecoGenParticleTrail recoGenParticleTrail_
Definition: HistoryBase.h:118
RecoGenParticleTrail const& HistoryBase::recoGenParticleTrail ( ) const
inline

Return all reco::GenParticle in the history.

Definition at line 64 of file HistoryBase.h.

References recoGenParticleTrail_.

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

64 { return recoGenParticleTrail_; }
RecoGenParticleTrail recoGenParticleTrail_
Definition: HistoryBase.h:118
void HistoryBase::resetTrails ( )
inlineprivate

Reset trail functions.

Definition at line 149 of file HistoryBase.h.

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

Referenced by evaluate(), and resetTrails().

149  {
150  simParticleTrail_.clear();
151  simVertexTrail_.clear();
152  genVertexTrail_.clear();
153  genParticleTrail_.clear();
154  recoGenParticleTrail_.clear();
156  genVertexTrailHelper_.clear();
157  }
RecoGenParticleTrail recoGenParticleTrail_
Definition: HistoryBase.h:118
SimVertexTrail simVertexTrail_
Definition: HistoryBase.h:119
GenVertexTrail genVertexTrail_
Definition: HistoryBase.h:116
GenVertexTrailHelper genVertexTrailHelper_
Definition: HistoryBase.h:123
SimParticleTrail simParticleTrail_
Definition: HistoryBase.h:120
RecoGenParticleTrailHelper recoGenParticleTrailHelper_
Definition: HistoryBase.h:124
GenParticleTrail genParticleTrail_
Definition: HistoryBase.h:117
void HistoryBase::resetTrails ( TrackingParticleRef  tpr)
inlineprivate

Definition at line 159 of file HistoryBase.h.

References resetTrails(), and simParticleTrail_.

159  {
160  resetTrails();
161  simParticleTrail_.push_back(tpr);
162  }
void resetTrails()
Reset trail functions.
Definition: HistoryBase.h:149
SimParticleTrail simParticleTrail_
Definition: HistoryBase.h:120
const TrackingParticleRef& HistoryBase::simParticle ( ) const
inline

Return the initial tracking particle from the history.

Definition at line 67 of file HistoryBase.h.

References simParticleTrail_.

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

67 { return simParticleTrail_[0]; }
SimParticleTrail simParticleTrail_
Definition: HistoryBase.h:120
SimParticleTrail const& HistoryBase::simParticleTrail ( ) const
inline

Return all the simulated particle in the history.

Definition at line 55 of file HistoryBase.h.

References simParticleTrail_.

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

55 { return simParticleTrail_; }
SimParticleTrail simParticleTrail_
Definition: HistoryBase.h:120
const TrackingVertexRef& HistoryBase::simVertex ( ) const
inline

Return the initial tracking vertex from the history.

Definition at line 70 of file HistoryBase.h.

References simVertexTrail_.

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

70 { return simVertexTrail_[0]; }
SimVertexTrail simVertexTrail_
Definition: HistoryBase.h:119
SimVertexTrail const& HistoryBase::simVertexTrail ( ) const
inline

Return all the simulated vertices in the history.

Definition at line 52 of file HistoryBase.h.

References simVertexTrail_.

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

52 { return simVertexTrail_; }
SimVertexTrail simVertexTrail_
Definition: HistoryBase.h:119
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().

31  {
32  // Third stop criteria: status abs(depth_) particles after the hadronization.
33  // The after hadronization is done by detecting the pdg_id pythia code from 88
34  // to 99
35  if (genParticle->status() <= abs(depth_) && (genParticle->pdg_id() < 88 || genParticle->pdg_id() > 99)) {
36  genParticleTrail_.push_back(genParticle);
37  // Get the producer vertex and trace it history
38  traceGenHistory(genParticle->production_vertex());
39  }
40 }
void traceGenHistory(HepMC::GenParticle const *)
Definition: HistoryBase.cc:31
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
GenParticleTrail genParticleTrail_
Definition: HistoryBase.h:117
const HepMC::GenParticle * genParticle() const
Definition: HistoryBase.h:74
void HistoryBase::traceGenHistory ( HepMC::GenVertex const *  genVertex)
private

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

Definition at line 42 of file HistoryBase.cc.

References genVertexTrail_, genVertexTrailHelper_, and traceGenHistory().

42  {
43  // Verify if has a vertex associated
44  if (genVertex) {
45  // Skip if already exist in the collection
46  if (genVertexTrailHelper_.find(genVertex) != genVertexTrailHelper_.end())
47  return;
48  // Add vertex to the history
49  genVertexTrail_.push_back(genVertex);
50  genVertexTrailHelper_.insert(genVertex);
51  // Verify if the vertex has incoming particles
52  if (genVertex->particles_in_size())
53  traceGenHistory(*(genVertex->particles_in_const_begin()));
54  }
55 }
GenVertexTrail genVertexTrail_
Definition: HistoryBase.h:116
void traceGenHistory(HepMC::GenParticle const *)
Definition: HistoryBase.cc:31
GenVertexTrailHelper genVertexTrailHelper_
Definition: HistoryBase.h:123
void HistoryBase::traceRecoGenHistory ( reco::GenParticle const *  genParticle)
private

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

Definition at line 7 of file HistoryBase.cc.

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

Referenced by traceSimHistory().

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

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

Definition at line 57 of file HistoryBase.cc.

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

Referenced by evaluate(), and traceSimHistory().

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

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

Definition at line 76 of file HistoryBase.cc.

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

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

Member Data Documentation

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

Definition at line 117 of file HistoryBase.h.

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

GenVertexTrail HistoryBase::genVertexTrail_
protected

Definition at line 116 of file HistoryBase.h.

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

GenVertexTrailHelper HistoryBase::genVertexTrailHelper_
protected

Definition at line 123 of file HistoryBase.h.

Referenced by resetTrails(), and traceGenHistory().

RecoGenParticleTrail HistoryBase::recoGenParticleTrail_
protected
RecoGenParticleTrailHelper HistoryBase::recoGenParticleTrailHelper_
protected

Definition at line 124 of file HistoryBase.h.

Referenced by resetTrails(), and traceRecoGenHistory().

SimParticleTrail HistoryBase::simParticleTrail_
protected

Definition at line 120 of file HistoryBase.h.

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

SimVertexTrail HistoryBase::simVertexTrail_
protected

Definition at line 119 of file HistoryBase.h.

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