CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups 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 <sstream>
12 #include <string>
13 #include <vector>
14 
15 // user include files
18 
26 
29 
30 //
31 // class decleration
32 //
33 
35 public:
36  explicit VertexHistoryAnalyzer(const edm::ParameterSet &);
37 
38 private:
39  virtual void beginRun(const edm::Run &, const edm::EventSetup &);
40  void beginJob() override;
41  void analyze(const edm::Event &, const edm::EventSetup &) override;
42 
43  // Member data
44 
46 
48 
50 
51  std::string particleString(int) const;
52 
54 
55  std::string vertexString(HepMC::GenVertex::particles_in_const_iterator,
56  HepMC::GenVertex::particles_in_const_iterator,
57  HepMC::GenVertex::particles_out_const_iterator,
58  HepMC::GenVertex::particles_out_const_iterator) const;
59 };
60 
62  : classifier_(config, consumesCollector()) {
63  vertexProducer_ = config.getUntrackedParameter<edm::InputTag>("vertexProducer");
64  consumes<edm::View<reco::Vertex>>(vertexProducer_);
65 }
66 
68  // Set the classifier for a new event
69  classifier_.newEvent(event, setup);
70 
71  // Vertex collection
73  event.getByLabel(vertexProducer_, vertexCollection);
74 
75  // Get a constant reference to the track history associated to the classifier
76  VertexHistory const &tracer = classifier_.history();
77 
78  // Loop over the track collection.
79  for (std::size_t index = 0; index < vertexCollection->size(); index++) {
80  std::cout << std::endl << "History for vertex #" << index << " : " << std::endl;
81 
82  // Classify the track and detect for fakes
84  // Get the list of TrackingParticles associated to
86 
87  // Loop over all simParticles
88  for (std::size_t hindex = 0; hindex < simParticles.size(); hindex++) {
89  std::cout << " simParticles [" << hindex << "] : " << particleString(simParticles[hindex]->pdgId())
90  << std::endl;
91  }
92 
93  // Get the list of TrackingVertexes associated to
94  VertexHistory::SimVertexTrail simVertexes(tracer.simVertexTrail());
95 
96  // Loop over all simVertexes
97  if (!simVertexes.empty()) {
98  for (std::size_t hindex = 0; hindex < simVertexes.size(); hindex++) {
99  std::cout << " simVertex [" << hindex << "] : "
100  << vertexString(simVertexes[hindex]->sourceTracks(), simVertexes[hindex]->daughterTracks())
101  << std::endl;
102  }
103  } else
104  std::cout << " simVertex no found" << std::endl;
105 
106  // Get the list of GenParticles associated to
108 
109  // Loop over all genParticles
110  for (std::size_t hindex = 0; hindex < genParticles.size(); hindex++) {
111  std::cout << " genParticles [" << hindex << "] : " << particleString(genParticles[hindex]->pdg_id())
112  << std::endl;
113  }
114 
115  // Get the list of TrackingVertexes associated to
116  VertexHistory::GenVertexTrail genVertexes(tracer.genVertexTrail());
117 
118  // Loop over all simVertexes
119  if (!genVertexes.empty()) {
120  for (std::size_t hindex = 0; hindex < genVertexes.size(); hindex++) {
121  std::cout << " genVertex [" << hindex << "] : "
122  << vertexString(genVertexes[hindex]->particles_in_const_begin(),
123  genVertexes[hindex]->particles_in_const_end(),
124  genVertexes[hindex]->particles_out_const_begin(),
125  genVertexes[hindex]->particles_out_const_end())
126  << std::endl;
127  }
128  } else
129  std::cout << " genVertex no found" << std::endl;
130  } else
131  std::cout << " fake vertex" << std::endl;
132 
133  std::cout << " vertex categories : " << classifier_;
134  std::cout << std::endl;
135  }
136 }
137 
139  // Get the particles table.
140  setup.getData(pdt_);
141 }
142 
144 
146  ParticleData const *pid;
147 
148  std::ostringstream vDescription;
149 
150  HepPDT::ParticleID particleType(pdgId);
151 
152  if (particleType.isValid()) {
153  pid = pdt_->particle(particleType);
154  if (pid)
155  vDescription << pid->name();
156  else
157  vDescription << pdgId;
158  } else
159  vDescription << pdgId;
160 
161  return vDescription.str();
162 }
163 
165  const TrackingParticleRefVector &out) const {
166  ParticleData const *pid;
167 
168  std::ostringstream vDescription;
169 
170  for (std::size_t j = 0; j < in.size(); j++) {
171  if (!j)
172  vDescription << "(";
173 
174  HepPDT::ParticleID particleType(in[j]->pdgId());
175 
176  if (particleType.isValid()) {
177  pid = pdt_->particle(particleType);
178  if (pid)
179  vDescription << pid->name();
180  else
181  vDescription << in[j]->pdgId();
182  } else
183  vDescription << in[j]->pdgId();
184 
185  if (j == in.size() - 1)
186  vDescription << ")";
187  else
188  vDescription << ",";
189  }
190 
191  vDescription << "->";
192 
193  for (std::size_t j = 0; j < out.size(); j++) {
194  if (!j)
195  vDescription << "(";
196 
197  HepPDT::ParticleID particleType(out[j]->pdgId());
198 
199  if (particleType.isValid()) {
200  pid = pdt_->particle(particleType);
201  if (pid)
202  vDescription << pid->name();
203  else
204  vDescription << out[j]->pdgId();
205  } else
206  vDescription << out[j]->pdgId();
207 
208  if (j == out.size() - 1)
209  vDescription << ")";
210  else
211  vDescription << ",";
212  }
213 
214  return vDescription.str();
215 }
216 
217 std::string VertexHistoryAnalyzer::vertexString(HepMC::GenVertex::particles_in_const_iterator in_begin,
218  HepMC::GenVertex::particles_in_const_iterator in_end,
219  HepMC::GenVertex::particles_out_const_iterator out_begin,
220  HepMC::GenVertex::particles_out_const_iterator out_end) const {
221  ParticleData const *pid;
222 
223  std::ostringstream vDescription;
224 
225  std::size_t j = 0;
226 
227  HepMC::GenVertex::particles_in_const_iterator in, itmp;
228 
229  for (in = in_begin; in != in_end; in++, j++) {
230  if (!j)
231  vDescription << "(";
232 
233  HepPDT::ParticleID particleType((*in)->pdg_id());
234 
235  if (particleType.isValid()) {
236  pid = pdt_->particle(particleType);
237  if (pid)
238  vDescription << pid->name();
239  else
240  vDescription << (*in)->pdg_id();
241  } else
242  vDescription << (*in)->pdg_id();
243 
244  itmp = in;
245 
246  if (++itmp == in_end)
247  vDescription << ")";
248  else
249  vDescription << ",";
250  }
251 
252  vDescription << "->";
253  j = 0;
254 
255  HepMC::GenVertex::particles_out_const_iterator out, otmp;
256 
257  for (out = out_begin; out != out_end; out++, j++) {
258  if (!j)
259  vDescription << "(";
260 
261  HepPDT::ParticleID particleType((*out)->pdg_id());
262 
263  if (particleType.isValid()) {
264  pid = pdt_->particle(particleType);
265  if (pid)
266  vDescription << pid->name();
267  else
268  vDescription << (*out)->pdg_id();
269  } else
270  vDescription << (*out)->pdg_id();
271 
272  otmp = out;
273 
274  if (++otmp == out_end)
275  vDescription << ")";
276  else
277  vDescription << ",";
278  }
279 
280  return vDescription.str();
281 }
282 
This class traces the simulated and generated history of a given track.
Definition: VertexHistory.h:18
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:16
std::string vertexString(const TrackingParticleRefVector &, const TrackingParticleRefVector &) const
SimVertexTrail const & simVertexTrail() const
Return all the simulated vertices in the history.
Definition: HistoryBase.h:52
Get track history and classify it in function of their .
virtual void beginRun(const edm::Run &, const edm::EventSetup &)
void analyze(const edm::Event &, const edm::EventSetup &) override
tuple vertexCollection
std::vector< const HepMC::GenParticle * > GenParticleTrail
HepMC::GenParticle trail type.
Definition: HistoryBase.h:15
bool getData(T &iHolder) const
Definition: EventSetup.h:128
edm::ESHandle< ParticleDataTable > pdt_
std::string particleString(int) const
bool is(Category category) const
Returns track flag for a given category.
std::vector< const HepMC::GenVertex * > GenVertexTrail
GenVertex trail type.
Definition: HistoryBase.h:24
HepPDT::ParticleData ParticleData
VertexHistoryAnalyzer(const edm::ParameterSet &)
SimParticleTrail const & simParticleTrail() const
Return all the simulated particle in the history.
Definition: HistoryBase.h:55
tuple config
parse the configuration file
VertexClassifier const & evaluate(reco::VertexBaseRef const &)
Classify the RecoVertex in categories.
size_type size() const
Size of the RefVector.
Definition: RefVector.h:102
tuple cout
Definition: gather_cfg.py:144
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:58
std::vector< TrackingVertexRef > SimVertexTrail
SimVertex trail type.
Definition: HistoryBase.h:33
GenParticleTrail const & genParticleTrail() const
Return all generated particle (HepMC::GenParticle) in the history.
Definition: HistoryBase.h:61
Definition: Run.h:45
std::vector< TrackingParticleRef > SimParticleTrail
SimParticle trail type.
Definition: HistoryBase.h:30