CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 <string>
12 #include <sstream>
13 #include <vector>
14 
15 // user include files
16 
19 
27 
30 
31 //
32 // class decleration
33 //
34 
36 {
37 public:
38 
39  explicit TrackHistoryAnalyzer(const edm::ParameterSet&);
41 
42 private:
43 
44  virtual void beginRun(const edm::Run&, const edm::EventSetup&) override;
45  virtual void beginJob() override ;
46  virtual void analyze(const edm::Event&, const edm::EventSetup&) override;
47 
48  // Member data
49 
51 
52  std::size_t totalTracks_;
53 
55 
56  std::string particleString(int) const;
57 
59 
63  ) const;
64 
66  HepMC::GenVertex::particles_in_const_iterator,
67  HepMC::GenVertex::particles_in_const_iterator,
68  HepMC::GenVertex::particles_out_const_iterator,
69  HepMC::GenVertex::particles_out_const_iterator
70  ) const;
71 };
72 
73 
74 TrackHistoryAnalyzer::TrackHistoryAnalyzer(const edm::ParameterSet& config) : classifier_(config,consumesCollector())
75 {
76  trackProducer_ = config.getUntrackedParameter<edm::InputTag> ( "trackProducer" );
77  consumes<edm::View<reco::Track>>(trackProducer_);
78 }
79 
80 
82 
83 
85 {
86  // Track collection
88  event.getByLabel(trackProducer_, trackCollection);
89 
90  // Set the classifier for a new event
91  classifier_.newEvent(event, setup);
92 
93  // Get a constant reference to the track history associated to the classifier
94  TrackHistory const & tracer = classifier_.history();
95 
96  // Loop over the track collection.
97  for (std::size_t index = 0; index < trackCollection->size(); index++)
98  {
99  std::cout << std::endl << "History for track #" << index << " : " << std::endl;
100 
101  // Classify the track and detect for fakes
103  {
104  // Get the list of TrackingParticles associated to
105  TrackHistory::SimParticleTrail simParticles(tracer.simParticleTrail());
106 
107  // Loop over all simParticles
108  for (std::size_t hindex=0; hindex<simParticles.size(); hindex++)
109  {
110  std::cout << " simParticles [" << hindex << "] : "
111  << particleString(simParticles[hindex]->pdgId())
112  << std::endl;
113  }
114 
115  // Get the list of TrackingVertexes associated to
116  TrackHistory::SimVertexTrail simVertexes(tracer.simVertexTrail());
117 
118  // Loop over all simVertexes
119  if ( !simVertexes.empty() )
120  {
121  for (std::size_t hindex=0; hindex<simVertexes.size(); hindex++)
122  {
123  std::cout << " simVertex [" << hindex << "] : "
124  << vertexString(
125  simVertexes[hindex]->sourceTracks(),
126  simVertexes[hindex]->daughterTracks()
127  )
128  << std::endl;
129  }
130  }
131  else
132  std::cout << " simVertex no found" << std::endl;
133 
134  // Get the list of GenParticles associated to
136 
137  // Loop over all genParticles
138  for (std::size_t hindex=0; hindex<genParticles.size(); hindex++)
139  {
140  std::cout << " genParticles [" << hindex << "] : "
141  << particleString(genParticles[hindex]->pdg_id())
142  << std::endl;
143  }
144 
145  // Get the list of TrackingVertexes associated to
146  TrackHistory::GenVertexTrail genVertexes(tracer.genVertexTrail());
147 
148  // Loop over all simVertexes
149  if ( !genVertexes.empty() )
150  {
151  for (std::size_t hindex=0; hindex<genVertexes.size(); hindex++)
152  {
153  std::cout << " genVertex [" << hindex << "] : "
154  << vertexString(
155  genVertexes[hindex]->particles_in_const_begin(),
156  genVertexes[hindex]->particles_in_const_end(),
157  genVertexes[hindex]->particles_out_const_begin(),
158  genVertexes[hindex]->particles_out_const_end()
159  )
160  << std::endl;
161  }
162  }
163  else
164  std::cout << " genVertex no found" << std::endl;
165  }
166  else
167  std::cout << " fake track" << std::endl;
168 
169  std::cout << " track categories : " << classifier_;
170  std::cout << std::endl;
171  }
172 }
173 
174 
175 void
177 {
178  // Get the particles table.
179  setup.getData( pdt_ );
180 }
181 
182 void
184 {
185  totalTracks_ = 0;
186 }
187 
188 
190 {
191  ParticleData const * pid;
192 
193  std::ostringstream vDescription;
194 
195  HepPDT::ParticleID particleType(pdgId);
196 
197  if (particleType.isValid())
198  {
199  pid = pdt_->particle(particleType);
200  if (pid)
201  vDescription << pid->name();
202  else
203  vDescription << pdgId;
204  }
205  else
206  vDescription << pdgId;
207 
208  return vDescription.str();
209 }
210 
211 
215 ) const
216 {
217  ParticleData const * pid;
218 
219  std::ostringstream vDescription;
220 
221  for (std::size_t j = 0; j < in.size(); j++)
222  {
223  if (!j) vDescription << "(";
224 
225  HepPDT::ParticleID particleType(in[j]->pdgId());
226 
227  if (particleType.isValid())
228  {
229  pid = pdt_->particle(particleType);
230  if (pid)
231  vDescription << pid->name();
232  else
233  vDescription << in[j]->pdgId();
234  }
235  else
236  vDescription << in[j]->pdgId();
237 
238  if (j == in.size() - 1) vDescription << ")";
239  else vDescription << ",";
240  }
241 
242  vDescription << "->";
243 
244  for (std::size_t j = 0; j < out.size(); j++)
245  {
246  if (!j) vDescription << "(";
247 
248  HepPDT::ParticleID particleType(out[j]->pdgId());
249 
250  if (particleType.isValid())
251  {
252  pid = pdt_->particle(particleType);
253  if (pid)
254  vDescription << pid->name();
255  else
256  vDescription << out[j]->pdgId();
257  }
258  else
259  vDescription << out[j]->pdgId();
260 
261  if (j == out.size() - 1) vDescription << ")";
262  else vDescription << ",";
263  }
264 
265  return vDescription.str();
266 }
267 
268 
270  HepMC::GenVertex::particles_in_const_iterator in_begin,
271  HepMC::GenVertex::particles_in_const_iterator in_end,
272  HepMC::GenVertex::particles_out_const_iterator out_begin,
273  HepMC::GenVertex::particles_out_const_iterator out_end
274 ) const
275 {
276  ParticleData const * pid;
277 
278  std::ostringstream vDescription;
279 
280  std::size_t j = 0;
281 
282  HepMC::GenVertex::particles_in_const_iterator in, itmp;
283 
284  for (in = in_begin; in != in_end; in++, j++)
285  {
286  if (!j) vDescription << "(";
287 
288  HepPDT::ParticleID particleType((*in)->pdg_id());
289 
290  if (particleType.isValid())
291  {
292  pid = pdt_->particle(particleType);
293  if (pid)
294  vDescription << pid->name();
295  else
296  vDescription << (*in)->pdg_id();
297  }
298  else
299  vDescription << (*in)->pdg_id();
300 
301  itmp = in;
302 
303  if (++itmp == in_end) vDescription << ")";
304  else vDescription << ",";
305  }
306 
307  vDescription << "->";
308  j = 0;
309 
310  HepMC::GenVertex::particles_out_const_iterator out, otmp;
311 
312  for (out = out_begin; out != out_end; out++, j++)
313  {
314  if (!j) vDescription << "(";
315 
316  HepPDT::ParticleID particleType((*out)->pdg_id());
317 
318  if (particleType.isValid())
319  {
320  pid = pdt_->particle(particleType);
321  if (pid)
322  vDescription << pid->name();
323  else
324  vDescription << (*out)->pdg_id();
325  }
326  else
327  vDescription << (*out)->pdg_id();
328 
329  otmp = out;
330 
331  if (++otmp == out_end) vDescription << ")";
332  else vDescription << ",";
333  }
334 
335  return vDescription.str();
336 }
337 
338 
T getUntrackedParameter(std::string const &, T const &) const
std::string vertexString(const TrackingParticleRefVector &, const TrackingParticleRefVector &) const
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
SimVertexTrail const & simVertexTrail() const
Return all the simulated vertices in the history.
Definition: HistoryBase.h:53
virtual void analyze(const edm::Event &, const edm::EventSetup &) override
std::string particleString(int) const
virtual void beginJob() override
TrackHistory const & history() const
Returns a reference to the track history used in the classification.
void getData(T &iHolder) const
Definition: EventSetup.h:79
virtual void beginRun(const edm::Run &, const edm::EventSetup &) override
std::vector< const HepMC::GenParticle * > GenParticleTrail
GenParticle trail type.
Definition: HistoryBase.h:18
TrackClassifier const & evaluate(reco::TrackBaseRef const &)
Classify the RecoTrack in categories.
TrackHistoryAnalyzer(const edm::ParameterSet &)
int j
Definition: DBlmapReader.cc:9
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:21
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
void newEvent(edm::Event const &, edm::EventSetup const &)
Pre-process event information (for accessing reconstraction information)
HepPDT::ParticleData ParticleData
tuple out
Definition: dbtoconf.py:99
bool is(Category category) const
Returns track flag for a given category.
Get track history and classify it in function of their .
tuple pid
Definition: sysUtil.py:22
SimParticleTrail const & simParticleTrail() const
Return all the simulated particle in the history.
Definition: HistoryBase.h:59
edm::ESHandle< ParticleDataTable > pdt_
size_type size() const
Size of the RefVector.
Definition: RefVector.h:99
tuple cout
Definition: gather_cfg.py:121
GenVertexTrail const & genVertexTrail() const
Return all generated vertex in the history.
Definition: HistoryBase.h:65
std::vector< TrackingVertexRef > SimVertexTrail
SimVertex trail type.
Definition: HistoryBase.h:30
void setup(std::vector< TH2F > &depth, std::string name, std::string units="")
GenParticleTrail const & genParticleTrail() const
Return all generated particle in the history.
Definition: HistoryBase.h:71
Definition: Run.h:43
std::vector< TrackingParticleRef > SimParticleTrail
SimParticle trail type.
Definition: HistoryBase.h:27