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 
75 {
76  trackProducer_ = config.getUntrackedParameter<edm::InputTag> ( "trackProducer" );
77 }
78 
79 
81 
82 
84 {
85  // Track collection
86  edm::Handle<edm::View<reco::Track> > trackCollection;
87  event.getByLabel(trackProducer_, trackCollection);
88 
89  // Set the classifier for a new event
90  classifier_.newEvent(event, setup);
91 
92  // Get a constant reference to the track history associated to the classifier
93  TrackHistory const & tracer = classifier_.history();
94 
95  // Loop over the track collection.
96  for (std::size_t index = 0; index < trackCollection->size(); index++)
97  {
98  std::cout << std::endl << "History for track #" << index << " : " << std::endl;
99 
100  // Classify the track and detect for fakes
102  {
103  // Get the list of TrackingParticles associated to
104  TrackHistory::SimParticleTrail simParticles(tracer.simParticleTrail());
105 
106  // Loop over all simParticles
107  for (std::size_t hindex=0; hindex<simParticles.size(); hindex++)
108  {
109  std::cout << " simParticles [" << hindex << "] : "
110  << particleString(simParticles[hindex]->pdgId())
111  << std::endl;
112  }
113 
114  // Get the list of TrackingVertexes associated to
115  TrackHistory::SimVertexTrail simVertexes(tracer.simVertexTrail());
116 
117  // Loop over all simVertexes
118  if ( !simVertexes.empty() )
119  {
120  for (std::size_t hindex=0; hindex<simVertexes.size(); hindex++)
121  {
122  std::cout << " simVertex [" << hindex << "] : "
123  << vertexString(
124  simVertexes[hindex]->sourceTracks(),
125  simVertexes[hindex]->daughterTracks()
126  )
127  << std::endl;
128  }
129  }
130  else
131  std::cout << " simVertex no found" << std::endl;
132 
133  // Get the list of GenParticles associated to
135 
136  // Loop over all genParticles
137  for (std::size_t hindex=0; hindex<genParticles.size(); hindex++)
138  {
139  std::cout << " genParticles [" << hindex << "] : "
140  << particleString(genParticles[hindex]->pdg_id())
141  << std::endl;
142  }
143 
144  // Get the list of TrackingVertexes associated to
145  TrackHistory::GenVertexTrail genVertexes(tracer.genVertexTrail());
146 
147  // Loop over all simVertexes
148  if ( !genVertexes.empty() )
149  {
150  for (std::size_t hindex=0; hindex<genVertexes.size(); hindex++)
151  {
152  std::cout << " genVertex [" << hindex << "] : "
153  << vertexString(
154  genVertexes[hindex]->particles_in_const_begin(),
155  genVertexes[hindex]->particles_in_const_end(),
156  genVertexes[hindex]->particles_out_const_begin(),
157  genVertexes[hindex]->particles_out_const_end()
158  )
159  << std::endl;
160  }
161  }
162  else
163  std::cout << " genVertex no found" << std::endl;
164  }
165  else
166  std::cout << " fake track" << std::endl;
167 
168  std::cout << " track categories : " << classifier_;
169  std::cout << std::endl;
170  }
171 }
172 
173 
174 void
176 {
177  // Get the particles table.
178  setup.getData( pdt_ );
179 }
180 
181 void
183 {
184  totalTracks_ = 0;
185 }
186 
187 
189 {
190  ParticleData const * pid;
191 
192  std::ostringstream vDescription;
193 
195 
196  if (particleType.isValid())
197  {
198  pid = pdt_->particle(particleType);
199  if (pid)
200  vDescription << pid->name();
201  else
202  vDescription << pdgId;
203  }
204  else
205  vDescription << pdgId;
206 
207  return vDescription.str();
208 }
209 
210 
214 ) const
215 {
216  ParticleData const * pid;
217 
218  std::ostringstream vDescription;
219 
220  for (std::size_t j = 0; j < in.size(); j++)
221  {
222  if (!j) vDescription << "(";
223 
225 
226  if (particleType.isValid())
227  {
228  pid = pdt_->particle(particleType);
229  if (pid)
230  vDescription << pid->name();
231  else
232  vDescription << in[j]->pdgId();
233  }
234  else
235  vDescription << in[j]->pdgId();
236 
237  if (j == in.size() - 1) vDescription << ")";
238  else vDescription << ",";
239  }
240 
241  vDescription << "->";
242 
243  for (std::size_t j = 0; j < out.size(); j++)
244  {
245  if (!j) vDescription << "(";
246 
248 
249  if (particleType.isValid())
250  {
251  pid = pdt_->particle(particleType);
252  if (pid)
253  vDescription << pid->name();
254  else
255  vDescription << out[j]->pdgId();
256  }
257  else
258  vDescription << out[j]->pdgId();
259 
260  if (j == out.size() - 1) vDescription << ")";
261  else vDescription << ",";
262  }
263 
264  return vDescription.str();
265 }
266 
267 
269  HepMC::GenVertex::particles_in_const_iterator in_begin,
270  HepMC::GenVertex::particles_in_const_iterator in_end,
271  HepMC::GenVertex::particles_out_const_iterator out_begin,
272  HepMC::GenVertex::particles_out_const_iterator out_end
273 ) const
274 {
275  ParticleData const * pid;
276 
277  std::ostringstream vDescription;
278 
279  std::size_t j = 0;
280 
281  HepMC::GenVertex::particles_in_const_iterator in, itmp;
282 
283  for (in = in_begin; in != in_end; in++, j++)
284  {
285  if (!j) vDescription << "(";
286 
287  HepPDT::ParticleID particleType((*in)->pdg_id());
288 
289  if (particleType.isValid())
290  {
291  pid = pdt_->particle(particleType);
292  if (pid)
293  vDescription << pid->name();
294  else
295  vDescription << (*in)->pdg_id();
296  }
297  else
298  vDescription << (*in)->pdg_id();
299 
300  itmp = in;
301 
302  if (++itmp == in_end) vDescription << ")";
303  else vDescription << ",";
304  }
305 
306  vDescription << "->";
307  j = 0;
308 
309  HepMC::GenVertex::particles_out_const_iterator out, otmp;
310 
311  for (out = out_begin; out != out_end; out++, j++)
312  {
313  if (!j) vDescription << "(";
314 
315  HepPDT::ParticleID particleType((*out)->pdg_id());
316 
317  if (particleType.isValid())
318  {
319  pid = pdt_->particle(particleType);
320  if (pid)
321  vDescription << pid->name();
322  else
323  vDescription << (*out)->pdg_id();
324  }
325  else
326  vDescription << (*out)->pdg_id();
327 
328  otmp = out;
329 
330  if (++otmp == out_end) vDescription << ")";
331  else vDescription << ",";
332  }
333 
334  return vDescription.str();
335 }
336 
337 
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:78
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.
tuple particleType
Definition: autophobj.py:21
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:16
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:89
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:41
std::vector< TrackingParticleRef > SimParticleTrail
SimParticle trail type.
Definition: HistoryBase.h:27