CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
VertexHistoryAnalyzer.cc
Go to the documentation of this file.
1 /*
2  * VertexHistoryAnalyzer.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
18 
26 
29 
30 //
31 // class decleration
32 //
33 
35 {
36 public:
37 
39 
40 private:
41 
42  virtual void beginRun(const edm::Run&,const edm::EventSetup&);
43  virtual void beginJob() ;
44  virtual void analyze(const edm::Event&, const edm::EventSetup&);
45 
46  // Member data
47 
49 
51 
53 
54  std::string particleString(int) const;
55 
59  ) const;
60 
62  HepMC::GenVertex::particles_in_const_iterator,
63  HepMC::GenVertex::particles_in_const_iterator,
64  HepMC::GenVertex::particles_out_const_iterator,
65  HepMC::GenVertex::particles_out_const_iterator
66  ) const;
67 };
68 
69 
71 {
72  vertexProducer_ = config.getUntrackedParameter<edm::InputTag> ("vertexProducer");
73 }
74 
75 
77 {
78  // Set the classifier for a new event
79  classifier_.newEvent(event, setup);
80 
81  // Vertex collection
83  event.getByLabel(vertexProducer_, vertexCollection);
84 
85  // Get a constant reference to the track history associated to the classifier
86  VertexHistory const & tracer = classifier_.history();
87 
88  // Loop over the track collection.
89  for (std::size_t index = 0; index < vertexCollection->size(); index++)
90  {
91  std::cout << std::endl << "History for vertex #" << index << " : " << std::endl;
92 
93  // Classify the track and detect for fakes
95  {
96  // Get the list of TrackingParticles associated to
98 
99  // Loop over all simParticles
100  for (std::size_t hindex=0; hindex<simParticles.size(); hindex++)
101  {
102  std::cout << " simParticles [" << hindex << "] : "
103  << particleString(simParticles[hindex]->pdgId())
104  << std::endl;
105  }
106 
107  // Get the list of TrackingVertexes associated to
108  VertexHistory::SimVertexTrail simVertexes(tracer.simVertexTrail());
109 
110  // Loop over all simVertexes
111  if ( !simVertexes.empty() )
112  {
113  for (std::size_t hindex=0; hindex<simVertexes.size(); hindex++)
114  {
115  std::cout << " simVertex [" << hindex << "] : "
116  << vertexString(
117  simVertexes[hindex]->sourceTracks(),
118  simVertexes[hindex]->daughterTracks()
119  )
120  << std::endl;
121  }
122  }
123  else
124  std::cout << " simVertex no found" << std::endl;
125 
126  // Get the list of GenParticles associated to
128 
129  // Loop over all genParticles
130  for (std::size_t hindex=0; hindex<genParticles.size(); hindex++)
131  {
132  std::cout << " genParticles [" << hindex << "] : "
133  << particleString(genParticles[hindex]->pdg_id())
134  << std::endl;
135  }
136 
137  // Get the list of TrackingVertexes associated to
138  VertexHistory::GenVertexTrail genVertexes(tracer.genVertexTrail());
139 
140  // Loop over all simVertexes
141  if ( !genVertexes.empty() )
142  {
143  for (std::size_t hindex=0; hindex<genVertexes.size(); hindex++)
144  {
145  std::cout << " genVertex [" << hindex << "] : "
146  << vertexString(
147  genVertexes[hindex]->particles_in_const_begin(),
148  genVertexes[hindex]->particles_in_const_end(),
149  genVertexes[hindex]->particles_out_const_begin(),
150  genVertexes[hindex]->particles_out_const_end()
151  )
152  << std::endl;
153  }
154  }
155  else
156  std::cout << " genVertex no found" << std::endl;
157  }
158  else
159  std::cout << " fake vertex" << std::endl;
160 
161  std::cout << " vertex categories : " << classifier_;
162  std::cout << std::endl;
163  }
164 }
165 
166 
167 void
169 {
170  // Get the particles table.
171  setup.getData( pdt_ );
172 }
173 
174 void
176 {
177 
178 }
179 
180 
182 {
183  ParticleData const * pid;
184 
185  std::ostringstream vDescription;
186 
187  HepPDT::ParticleID particleType(pdgId);
188 
189  if (particleType.isValid())
190  {
191  pid = pdt_->particle(particleType);
192  if (pid)
193  vDescription << pid->name();
194  else
195  vDescription << pdgId;
196  }
197  else
198  vDescription << pdgId;
199 
200  return vDescription.str();
201 }
202 
203 
207 ) const
208 {
209  ParticleData const * pid;
210 
211  std::ostringstream vDescription;
212 
213  for (std::size_t j = 0; j < in.size(); j++)
214  {
215  if (!j) vDescription << "(";
216 
217  HepPDT::ParticleID particleType(in[j]->pdgId());
218 
219  if (particleType.isValid())
220  {
221  pid = pdt_->particle(particleType);
222  if (pid)
223  vDescription << pid->name();
224  else
225  vDescription << in[j]->pdgId();
226  }
227  else
228  vDescription << in[j]->pdgId();
229 
230  if (j == in.size() - 1) vDescription << ")";
231  else vDescription << ",";
232  }
233 
234  vDescription << "->";
235 
236  for (std::size_t j = 0; j < out.size(); j++)
237  {
238  if (!j) vDescription << "(";
239 
240  HepPDT::ParticleID particleType(out[j]->pdgId());
241 
242  if (particleType.isValid())
243  {
244  pid = pdt_->particle(particleType);
245  if (pid)
246  vDescription << pid->name();
247  else
248  vDescription << out[j]->pdgId();
249  }
250  else
251  vDescription << out[j]->pdgId();
252 
253  if (j == out.size() - 1) vDescription << ")";
254  else vDescription << ",";
255  }
256 
257  return vDescription.str();
258 }
259 
260 
262  HepMC::GenVertex::particles_in_const_iterator in_begin,
263  HepMC::GenVertex::particles_in_const_iterator in_end,
264  HepMC::GenVertex::particles_out_const_iterator out_begin,
265  HepMC::GenVertex::particles_out_const_iterator out_end
266 ) const
267 {
268  ParticleData const * pid;
269 
270  std::ostringstream vDescription;
271 
272  std::size_t j = 0;
273 
274  HepMC::GenVertex::particles_in_const_iterator in, itmp;
275 
276  for (in = in_begin; in != in_end; in++, j++)
277  {
278  if (!j) vDescription << "(";
279 
280  HepPDT::ParticleID particleType((*in)->pdg_id());
281 
282  if (particleType.isValid())
283  {
284  pid = pdt_->particle(particleType);
285  if (pid)
286  vDescription << pid->name();
287  else
288  vDescription << (*in)->pdg_id();
289  }
290  else
291  vDescription << (*in)->pdg_id();
292 
293  itmp = in;
294 
295  if (++itmp == in_end) vDescription << ")";
296  else vDescription << ",";
297  }
298 
299  vDescription << "->";
300  j = 0;
301 
302  HepMC::GenVertex::particles_out_const_iterator out, otmp;
303 
304  for (out = out_begin; out != out_end; out++, j++)
305  {
306  if (!j) vDescription << "(";
307 
308  HepPDT::ParticleID particleType((*out)->pdg_id());
309 
310  if (particleType.isValid())
311  {
312  pid = pdt_->particle(particleType);
313  if (pid)
314  vDescription << pid->name();
315  else
316  vDescription << (*out)->pdg_id();
317  }
318  else
319  vDescription << (*out)->pdg_id();
320 
321  otmp = out;
322 
323  if (++otmp == out_end) vDescription << ")";
324  else vDescription << ",";
325  }
326 
327  return vDescription.str();
328 }
329 
330 
This class traces the simulated and generated history of a given track.
Definition: VertexHistory.h:17
T getUntrackedParameter(std::string const &, T const &) const
VertexHistory const & history() const
Returns a reference to the vertex history used in the classification.
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
std::string vertexString(const TrackingParticleRefVector &, const TrackingParticleRefVector &) const
SimVertexTrail const & simVertexTrail() const
Return all the simulated vertices in the history.
Definition: HistoryBase.h:53
Get track history and classify it in function of their .
virtual void beginRun(const edm::Run &, const edm::EventSetup &)
tuple vertexCollection
void getData(T &iHolder) const
Definition: EventSetup.h:67
std::vector< const HepMC::GenParticle * > GenParticleTrail
GenParticle trail type.
Definition: HistoryBase.h:18
virtual void analyze(const edm::Event &, const edm::EventSetup &)
edm::ESHandle< ParticleDataTable > pdt_
std::string particleString(int) const
int j
Definition: DBlmapReader.cc:9
bool is(Category category) const
Returns track flag for a given category.
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
HepPDT::ParticleData ParticleData
tuple out
Definition: dbtoconf.py:99
VertexHistoryAnalyzer(const edm::ParameterSet &)
SimParticleTrail const & simParticleTrail() const
Return all the simulated particle in the history.
Definition: HistoryBase.h:59
VertexClassifier const & evaluate(reco::VertexBaseRef const &)
Classify the RecoVertex in categories.
size_type size() const
Size of the RefVector.
Definition: RefVector.h:89
tuple cout
Definition: gather_cfg.py:121
virtual void newEvent(edm::Event const &, edm::EventSetup const &)
Pre-process event information (for accessing reconstraction information)
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:36
std::vector< TrackingParticleRef > SimParticleTrail
SimParticle trail type.
Definition: HistoryBase.h:27