CMS 3D CMS Logo

TrackHistoryAnalyzer.cc
Go to the documentation of this file.
1 /*
2  * TrackHistoryAnalyzer.C
3  *
4  * Created by Victor Eduardo Bazterra on 5/31/07.
5  *
6  */
7 
8 // system include files
9 #include <iostream>
10 #include <memory>
11 #include <sstream>
12 #include <string>
13 #include <vector>
14 
15 // user include files
16 
19 
27 
30 
31 //
32 // class decleration
33 //
34 
35 class TrackHistoryAnalyzer : public edm::one::EDAnalyzer<edm::one::WatchRuns> {
36 public:
37  explicit TrackHistoryAnalyzer(const edm::ParameterSet &);
38  ~TrackHistoryAnalyzer() override = default;
39 
40 private:
41  void beginJob() override;
42  void beginRun(const edm::Run &, const edm::EventSetup &) override;
43  void analyze(const edm::Event &, const edm::EventSetup &) override;
44  void endRun(const edm::Run &, const edm::EventSetup &) override{};
45 
46  // Member data
49 
50  std::size_t totalTracks_;
51 
53 
54  std::string particleString(int) const;
55 
57 
59 
60  std::string vertexString(HepMC::GenVertex::particles_in_const_iterator,
61  HepMC::GenVertex::particles_in_const_iterator,
62  HepMC::GenVertex::particles_out_const_iterator,
63  HepMC::GenVertex::particles_out_const_iterator) const;
64 };
65 
67  : pdtToken_(esConsumes<edm::Transition::BeginRun>()),
68  trkToken_(consumes<edm::View<reco::Track>>(config.getUntrackedParameter<edm::InputTag>("trackProducer"))),
69  classifier_(config, consumesCollector()) {}
70 
72  // Track collection
74  event.getByToken(trkToken_, trackCollection);
75 
76  // Set the classifier for a new event
78 
79  // Get a constant reference to the track history associated to the classifier
80  TrackHistory const &tracer = classifier_.history();
81 
82  // Loop over the track collection.
83  for (std::size_t index = 0; index < trackCollection->size(); index++) {
84  edm::LogPrint("TrackHistoryAnalyzer") << std::endl << "History for track #" << index << " : ";
85 
86  // Classify the track and detect for fakes
88  // Get the list of TrackingParticles associated to
89  TrackHistory::SimParticleTrail simParticles(tracer.simParticleTrail());
90 
91  // Loop over all simParticles
92  for (std::size_t hindex = 0; hindex < simParticles.size(); hindex++) {
93  edm::LogPrint("TrackHistoryAnalyzer")
94  << " simParticles [" << hindex << "] : " << particleString(simParticles[hindex]->pdgId());
95  }
96 
97  // Get the list of TrackingVertexes associated to
98  TrackHistory::SimVertexTrail simVertexes(tracer.simVertexTrail());
99 
100  // Loop over all simVertexes
101  if (!simVertexes.empty()) {
102  for (std::size_t hindex = 0; hindex < simVertexes.size(); hindex++) {
103  edm::LogPrint("TrackHistoryAnalyzer")
104  << " simVertex [" << hindex
105  << "] : " << vertexString(simVertexes[hindex]->sourceTracks(), simVertexes[hindex]->daughterTracks());
106  }
107  } else
108  edm::LogPrint("TrackHistoryAnalyzer") << " simVertex no found";
109 
110  // Get the list of GenParticles associated to
112 
113  // Loop over all genParticles
114  for (std::size_t hindex = 0; hindex < genParticles.size(); hindex++) {
115  edm::LogPrint("TrackHistoryAnalyzer")
116  << " genParticles [" << hindex << "] : " << particleString(genParticles[hindex]->pdg_id());
117  }
118 
119  // Get the list of TrackingVertexes associated to
120  TrackHistory::GenVertexTrail genVertexes(tracer.genVertexTrail());
121 
122  // Loop over all simVertexes
123  if (!genVertexes.empty()) {
124  for (std::size_t hindex = 0; hindex < genVertexes.size(); hindex++) {
125  edm::LogPrint("TrackHistoryAnalyzer") << " genVertex [" << hindex << "] : "
126  << vertexString(genVertexes[hindex]->particles_in_const_begin(),
127  genVertexes[hindex]->particles_in_const_end(),
128  genVertexes[hindex]->particles_out_const_begin(),
129  genVertexes[hindex]->particles_out_const_end());
130  }
131  } else
132  edm::LogPrint("TrackHistoryAnalyzer") << " genVertex no found";
133  } else
134  edm::LogPrint("TrackHistoryAnalyzer") << " fake track";
135 
136  edm::LogPrint("TrackHistoryAnalyzer") << " track categories : " << classifier_;
137  }
138 }
139 
141  // Get the particles table.
142  pdt_ = setup.getHandle(pdtToken_);
143 }
144 
146 
148  ParticleData const *pid;
149 
150  std::ostringstream vDescription;
151 
153 
154  if (particleType.isValid()) {
155  pid = pdt_->particle(particleType);
156  if (pid)
157  vDescription << pid->name();
158  else
159  vDescription << pdgId;
160  } else
161  vDescription << pdgId;
162 
163  return vDescription.str();
164 }
165 
167  const TrackingParticleRefVector &out) const {
168  ParticleData const *pid;
169 
170  std::ostringstream vDescription;
171 
172  for (std::size_t j = 0; j < in.size(); j++) {
173  if (!j)
174  vDescription << "(";
175 
177 
178  if (particleType.isValid()) {
179  pid = pdt_->particle(particleType);
180  if (pid)
181  vDescription << pid->name();
182  else
183  vDescription << in[j]->pdgId();
184  } else
185  vDescription << in[j]->pdgId();
186 
187  if (j == in.size() - 1)
188  vDescription << ")";
189  else
190  vDescription << ",";
191  }
192 
193  vDescription << "->";
194 
195  for (std::size_t j = 0; j < out.size(); j++) {
196  if (!j)
197  vDescription << "(";
198 
200 
201  if (particleType.isValid()) {
202  pid = pdt_->particle(particleType);
203  if (pid)
204  vDescription << pid->name();
205  else
206  vDescription << out[j]->pdgId();
207  } else
208  vDescription << out[j]->pdgId();
209 
210  if (j == out.size() - 1)
211  vDescription << ")";
212  else
213  vDescription << ",";
214  }
215 
216  return vDescription.str();
217 }
218 
219 std::string TrackHistoryAnalyzer::vertexString(HepMC::GenVertex::particles_in_const_iterator in_begin,
220  HepMC::GenVertex::particles_in_const_iterator in_end,
221  HepMC::GenVertex::particles_out_const_iterator out_begin,
222  HepMC::GenVertex::particles_out_const_iterator out_end) const {
223  ParticleData const *pid;
224 
225  std::ostringstream vDescription;
226 
227  std::size_t j = 0;
228 
229  HepMC::GenVertex::particles_in_const_iterator in, itmp;
230 
231  for (in = in_begin; in != in_end; in++, j++) {
232  if (!j)
233  vDescription << "(";
234 
235  HepPDT::ParticleID particleType((*in)->pdg_id());
236 
237  if (particleType.isValid()) {
238  pid = pdt_->particle(particleType);
239  if (pid)
240  vDescription << pid->name();
241  else
242  vDescription << (*in)->pdg_id();
243  } else
244  vDescription << (*in)->pdg_id();
245 
246  itmp = in;
247 
248  if (++itmp == in_end)
249  vDescription << ")";
250  else
251  vDescription << ",";
252  }
253 
254  vDescription << "->";
255  j = 0;
256 
257  HepMC::GenVertex::particles_out_const_iterator out, otmp;
258 
259  for (out = out_begin; out != out_end; out++, j++) {
260  if (!j)
261  vDescription << "(";
262 
263  HepPDT::ParticleID particleType((*out)->pdg_id());
264 
265  if (particleType.isValid()) {
266  pid = pdt_->particle(particleType);
267  if (pid)
268  vDescription << pid->name();
269  else
270  vDescription << (*out)->pdg_id();
271  } else
272  vDescription << (*out)->pdg_id();
273 
274  otmp = out;
275 
276  if (++otmp == out_end)
277  vDescription << ")";
278  else
279  vDescription << ",";
280  }
281 
282  return vDescription.str();
283 }
284 
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
const edm::ESGetToken< ParticleDataTable, PDTRecord > pdtToken_
void endRun(const edm::Run &, const edm::EventSetup &) override
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
GenVertexTrail const & genVertexTrail() const
Return all generated vertex in the history.
Definition: HistoryBase.h:58
void analyze(const edm::Event &, const edm::EventSetup &) override
SimVertexTrail const & simVertexTrail() const
Return all the simulated vertices in the history.
Definition: HistoryBase.h:52
Definition: config.py:1
std::string vertexString(const TrackingParticleRefVector &, const TrackingParticleRefVector &) const
SimParticleTrail const & simParticleTrail() const
Return all the simulated particle in the history.
Definition: HistoryBase.h:55
void beginRun(const edm::Run &, const edm::EventSetup &) override
std::vector< const HepMC::GenParticle * > GenParticleTrail
HepMC::GenParticle trail type.
Definition: HistoryBase.h:15
TrackHistory const & history() const
Returns a reference to the track history used in the classification.
~TrackHistoryAnalyzer() override=default
bool is(Category category) const
Returns track flag for a given category.
TrackHistoryAnalyzer(const edm::ParameterSet &)
Transition
Definition: Transition.h:12
This class traces the simulated and generated history of a given track.
Definition: TrackHistory.h:17
std::vector< const HepMC::GenVertex * > GenVertexTrail
GenVertex trail type.
Definition: HistoryBase.h:24
void newEvent(edm::Event const &, edm::EventSetup const &)
Pre-process event information (for accessing reconstraction information)
const edm::EDGetTokenT< edm::View< reco::Track > > trkToken_
HepPDT::ParticleData ParticleData
Log< level::Warning, true > LogPrint
Get track history and classify it in function of their .
TrackClassifier const & evaluate(reco::TrackBaseRef const &)
Classify the RecoTrack in categories.
GenParticleTrail const & genParticleTrail() const
Return all generated particle (HepMC::GenParticle) in the history.
Definition: HistoryBase.h:61
fixed size matrix
HLT enums.
edm::ESHandle< ParticleDataTable > pdt_
std::string particleString(int) const
std::vector< TrackingVertexRef > SimVertexTrail
SimVertex trail type.
Definition: HistoryBase.h:33
Definition: event.py:1
Definition: Run.h:45
std::vector< TrackingParticleRef > SimParticleTrail
SimParticle trail type.
Definition: HistoryBase.h:30