CMS 3D CMS Logo

HistoryBase.h
Go to the documentation of this file.
1 #ifndef HistoryBase_h
2 #define HistoryBase_h
3 
4 #include <set>
5 
10 
13 {
14 
15 public:
16 
18  typedef std::vector<const HepMC::GenParticle *> GenParticleTrail;
19 
21  typedef std::vector<const reco::GenParticle *> RecoGenParticleTrail;
22 
24  typedef std::set<const reco::GenParticle *> RecoGenParticleTrailHelper;
25 
27  typedef std::vector<const HepMC::GenVertex *> GenVertexTrail;
28 
30  typedef std::set<const HepMC::GenVertex *> GenVertexTrailHelper;
31 
33  typedef std::vector<TrackingParticleRef> SimParticleTrail;
34 
36  typedef std::vector<TrackingVertexRef> SimVertexTrail;
37 
38  // Default constructor
40  {
41  // Default depth
42  depth_ = -1;
43  }
44 
46  /* Set TrackHistory to given depth. Positive values
47  constrain the number of TrackingVertex visit in the history.
48  Negatives values set the limit of the iteration over generated
49  information i.e. (-1 -> status 1 or -2 -> status 2 particles).
50 
51  /param[in] depth the history
52  */
53  void depth(int d)
54  {
55  depth_ = d;
56  }
57 
59  SimVertexTrail const & simVertexTrail() const
60  {
61  return simVertexTrail_;
62  }
63 
65  SimParticleTrail const & simParticleTrail() const
66  {
67  return simParticleTrail_;
68  }
69 
71  GenVertexTrail const & genVertexTrail() const
72  {
73  return genVertexTrail_;
74  }
75 
77  GenParticleTrail const & genParticleTrail() const
78  {
79  return genParticleTrail_;
80  }
81 
83  RecoGenParticleTrail const & recoGenParticleTrail() const
84  {
85  return recoGenParticleTrail_;
86  }
87 
90  {
91  return simParticleTrail_[0];
92  }
93 
95  const TrackingVertexRef & simVertex() const
96  {
97  return simVertexTrail_[0];
98  }
99 
102  {
103  if ( genParticleTrail_.empty() ) return 0;
104  return genParticleTrail_[genParticleTrail_.size()-1];
105  }
106 
109  {
110  if ( recoGenParticleTrail_.empty() ) return 0;
112  }
113 
115  /* Return false when the history cannot be determined upto a given depth.
116  If not depth is pass to the function no restriction are apply to it.
117 
118  /param[in] TrackingParticleRef of a simulated track
119  /param[in] depth of the track history
120  /param[out] boolean that is true when history can be determined
121  */
123  {
124  resetTrails(tpr);
125  return traceSimHistory(tpr, depth_);
126  }
127 
129  /* Return false when the history cannot be determined upto a given depth.
130  If not depth is pass to the function no restriction are apply to it.
131 
132  /param[in] TrackingVertexRef of a simulated vertex
133  /param[in] depth of the track history
134  /param[out] boolean that is true when history can be determined
135  */
137  {
138  resetTrails();
139  return traceSimHistory(tvr, depth_);
140  }
141 
142 protected:
143 
144  // History cointainers
145  GenVertexTrail genVertexTrail_;
146  GenParticleTrail genParticleTrail_;
147  RecoGenParticleTrail recoGenParticleTrail_;
148  SimVertexTrail simVertexTrail_;
149  SimParticleTrail simParticleTrail_;
150 
151  // Helper function to speedup search
152  GenVertexTrailHelper genVertexTrailHelper_;
153  RecoGenParticleTrailHelper recoGenParticleTrailHelper_;
154 
155 private:
156 
157  int depth_;
158 
160  bool traceSimHistory (TrackingParticleRef const &, int);
161 
163  bool traceSimHistory (TrackingVertexRef const &, int);
164 
166  void traceGenHistory (HepMC::GenParticle const *);
167 
169  void traceRecoGenHistory (reco::GenParticle const *);
170 
171 
173  void traceGenHistory (HepMC::GenVertex const *);
174 
176  void resetTrails()
177  {
178  simParticleTrail_.clear();
179  simVertexTrail_.clear();
180  genVertexTrail_.clear();
181  genParticleTrail_.clear();
182  recoGenParticleTrail_.clear();
183  recoGenParticleTrailHelper_.clear();
184  genVertexTrailHelper_.clear();
185  }
186 
188  {
189  resetTrails();
190  simParticleTrail_.push_back(tpr);
191  }
192 };
193 
194 #endif
std::set< const HepMC::GenVertex * > GenVertexTrailHelper
GenVertex trail helper type.
Definition: HistoryBase.h:30
RecoGenParticleTrail recoGenParticleTrail_
Definition: HistoryBase.h:147
void resetTrails()
Reset trail functions.
Definition: HistoryBase.h:176
const reco::GenParticle * recoGenParticle() const
Returns a pointer to most primitive status 1 or 2 particle in the recoGenParticleTrail_.
Definition: HistoryBase.h:108
SimVertexTrail const & simVertexTrail() const
Return all the simulated vertices in the history.
Definition: HistoryBase.h:59
const TrackingParticleRef & simParticle() const
Return the initial tracking particle from the history.
Definition: HistoryBase.h:89
SimVertexTrail simVertexTrail_
Definition: HistoryBase.h:148
std::vector< const HepMC::GenParticle * > GenParticleTrail
HepMC::GenParticle trail type.
Definition: HistoryBase.h:18
std::vector< const reco::GenParticle * > RecoGenParticleTrail
reco::GenParticle trail type.
Definition: HistoryBase.h:21
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
std::set< const reco::GenParticle * > RecoGenParticleTrailHelper
reco::GenParticle trail helper type.
Definition: HistoryBase.h:24
void resetTrails(TrackingParticleRef tpr)
Definition: HistoryBase.h:187
void traceRecoGenHistory(reco::GenParticle const *)
Trace all the simulated information for a given pointer to a reco::GenParticle.
Definition: HistoryBase.cc:8
std::vector< const HepMC::GenVertex * > GenVertexTrail
GenVertex trail type.
Definition: HistoryBase.h:27
const TrackingVertexRef & simVertex() const
Return the initial tracking vertex from the history.
Definition: HistoryBase.h:95
SimParticleTrail simParticleTrail_
Definition: HistoryBase.h:149
Base class to all the history types.
Definition: HistoryBase.h:12
SimParticleTrail const & simParticleTrail() const
Return all the simulated particle in the history.
Definition: HistoryBase.h:65
bool evaluate(TrackingParticleRef tpr)
Evaluate track history using a TrackingParticleRef.
Definition: HistoryBase.h:122
void depth(int d)
Set the depth of the history.
Definition: HistoryBase.h:53
bool evaluate(TrackingVertexRef tvr)
Evaluate track history using a TrackingParticleRef.
Definition: HistoryBase.h:136
bool traceSimHistory(TrackingParticleRef const &, int)
Trace all the simulated information for a given reference to a TrackingParticle.
Definition: HistoryBase.cc:62
RecoGenParticleTrailHelper recoGenParticleTrailHelper_
Definition: HistoryBase.h:153
GenVertexTrail const & genVertexTrail() const
Return all generated vertex in the history.
Definition: HistoryBase.h:71
RecoGenParticleTrail const & recoGenParticleTrail() const
Return all reco::GenParticle in the history.
Definition: HistoryBase.h:83
std::vector< TrackingVertexRef > SimVertexTrail
SimVertex trail type.
Definition: HistoryBase.h:36
GenParticleTrail genParticleTrail_
Definition: HistoryBase.h:146
GenParticleTrail const & genParticleTrail() const
Return all generated particle (HepMC::GenParticle) in the history.
Definition: HistoryBase.h:77
const HepMC::GenParticle * genParticle() const
Returns a pointer to most primitive status 1 or 2 particle in the genParticleTrail_.
Definition: HistoryBase.h:101
std::vector< TrackingParticleRef > SimParticleTrail
SimParticle trail type.
Definition: HistoryBase.h:33