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 
27 
30 
31 //
32 // class decleration
33 //
34 
35 class VertexHistoryAnalyzer : public edm::one::EDAnalyzer<edm::one::WatchRuns> {
36 public:
37  explicit VertexHistoryAnalyzer(const edm::ParameterSet &);
38  ~VertexHistoryAnalyzer() override = default;
39 
40 private:
41  void beginRun(const edm::Run &, const edm::EventSetup &) override;
42  void endRun(const edm::Run &, const edm::EventSetup &) override{};
43  void analyze(const edm::Event &, const edm::EventSetup &) override;
44 
45  // Member data
48 
50 
52 
53  std::string particleString(int) const;
54 
56 
57  std::string vertexString(HepMC::GenVertex::particles_in_const_iterator,
58  HepMC::GenVertex::particles_in_const_iterator,
59  HepMC::GenVertex::particles_out_const_iterator,
60  HepMC::GenVertex::particles_out_const_iterator) const;
61 };
62 
64  : pdtToken_(esConsumes<edm::Transition::BeginRun>()),
65  vtxToken_(consumes<edm::View<reco::Vertex>>(config.getUntrackedParameter<edm::InputTag>("vertexProducer"))),
66  classifier_(config, consumesCollector()) {}
67 
69  // Set the classifier for a new event
70  classifier_.newEvent(event, setup);
71 
72  // Vertex collection
74  event.getByToken(vtxToken_, vertexCollection);
75 
76  // Get a constant reference to the track history associated to the classifier
77  VertexHistory const &tracer = classifier_.history();
78 
79  // Loop over the track collection.
80  for (std::size_t index = 0; index < vertexCollection->size(); index++) {
81  edm::LogPrint("VertexHistoryAnalyzer") << std::endl << "History for vertex #" << index << " : ";
82 
83  // Classify the track and detect for fakes
84  if (!classifier_.evaluate(reco::VertexBaseRef(vertexCollection, index)).is(VertexClassifier::Fake)) {
85  // Get the list of TrackingParticles associated to
87 
88  // Loop over all simParticles
89  for (std::size_t hindex = 0; hindex < simParticles.size(); hindex++) {
90  edm::LogPrint("VertexHistoryAnalyzer")
91  << " simParticles [" << hindex << "] : " << particleString(simParticles[hindex]->pdgId());
92  }
93 
94  // Get the list of TrackingVertexes associated to
95  VertexHistory::SimVertexTrail simVertexes(tracer.simVertexTrail());
96 
97  // Loop over all simVertexes
98  if (!simVertexes.empty()) {
99  for (std::size_t hindex = 0; hindex < simVertexes.size(); hindex++) {
100  edm::LogPrint("VertexHistoryAnalyzer")
101  << " simVertex [" << hindex
102  << "] : " << vertexString(simVertexes[hindex]->sourceTracks(), simVertexes[hindex]->daughterTracks());
103  }
104  } else
105  edm::LogPrint("VertexHistoryAnalyzer") << " simVertex no found";
106 
107  // Get the list of GenParticles associated to
109 
110  // Loop over all genParticles
111  for (std::size_t hindex = 0; hindex < genParticles.size(); hindex++) {
112  edm::LogPrint("VertexHistoryAnalyzer")
113  << " genParticles [" << hindex << "] : " << particleString(genParticles[hindex]->pdg_id());
114  }
115 
116  // Get the list of TrackingVertexes associated to
117  VertexHistory::GenVertexTrail genVertexes(tracer.genVertexTrail());
118 
119  // Loop over all simVertexes
120  if (!genVertexes.empty()) {
121  for (std::size_t hindex = 0; hindex < genVertexes.size(); hindex++) {
122  edm::LogPrint("VertexHistoryAnalyzer") << " genVertex [" << hindex << "] : "
123  << vertexString(genVertexes[hindex]->particles_in_const_begin(),
124  genVertexes[hindex]->particles_in_const_end(),
125  genVertexes[hindex]->particles_out_const_begin(),
126  genVertexes[hindex]->particles_out_const_end());
127  }
128  } else
129  edm::LogPrint("VertexHistoryAnalyzer") << " genVertex no found";
130  } else
131  edm::LogPrint("VertexHistoryAnalyzer") << " fake vertex";
132 
133  edm::LogPrint("VertexHistoryAnalyzer") << " vertex categories : " << classifier_;
134  }
135 }
136 
138  // Get the particles table.
139  pdt_ = setup.getHandle(pdtToken_);
140 }
141 
143  ParticleData const *pid;
144 
145  std::ostringstream vDescription;
146 
147  HepPDT::ParticleID particleType(pdgId);
148 
149  if (particleType.isValid()) {
150  pid = pdt_->particle(particleType);
151  if (pid)
152  vDescription << pid->name();
153  else
154  vDescription << pdgId;
155  } else
156  vDescription << pdgId;
157 
158  return vDescription.str();
159 }
160 
162  const TrackingParticleRefVector &out) const {
163  ParticleData const *pid;
164 
165  std::ostringstream vDescription;
166 
167  for (std::size_t j = 0; j < in.size(); j++) {
168  if (!j)
169  vDescription << "(";
170 
171  HepPDT::ParticleID particleType(in[j]->pdgId());
172 
173  if (particleType.isValid()) {
174  pid = pdt_->particle(particleType);
175  if (pid)
176  vDescription << pid->name();
177  else
178  vDescription << in[j]->pdgId();
179  } else
180  vDescription << in[j]->pdgId();
181 
182  if (j == in.size() - 1)
183  vDescription << ")";
184  else
185  vDescription << ",";
186  }
187 
188  vDescription << "->";
189 
190  for (std::size_t j = 0; j < out.size(); j++) {
191  if (!j)
192  vDescription << "(";
193 
194  HepPDT::ParticleID particleType(out[j]->pdgId());
195 
196  if (particleType.isValid()) {
197  pid = pdt_->particle(particleType);
198  if (pid)
199  vDescription << pid->name();
200  else
201  vDescription << out[j]->pdgId();
202  } else
203  vDescription << out[j]->pdgId();
204 
205  if (j == out.size() - 1)
206  vDescription << ")";
207  else
208  vDescription << ",";
209  }
210 
211  return vDescription.str();
212 }
213 
214 std::string VertexHistoryAnalyzer::vertexString(HepMC::GenVertex::particles_in_const_iterator in_begin,
215  HepMC::GenVertex::particles_in_const_iterator in_end,
216  HepMC::GenVertex::particles_out_const_iterator out_begin,
217  HepMC::GenVertex::particles_out_const_iterator out_end) const {
218  ParticleData const *pid;
219 
220  std::ostringstream vDescription;
221 
222  std::size_t j = 0;
223 
224  HepMC::GenVertex::particles_in_const_iterator in, itmp;
225 
226  for (in = in_begin; in != in_end; in++, j++) {
227  if (!j)
228  vDescription << "(";
229 
230  HepPDT::ParticleID particleType((*in)->pdg_id());
231 
232  if (particleType.isValid()) {
233  pid = pdt_->particle(particleType);
234  if (pid)
235  vDescription << pid->name();
236  else
237  vDescription << (*in)->pdg_id();
238  } else
239  vDescription << (*in)->pdg_id();
240 
241  itmp = in;
242 
243  if (++itmp == in_end)
244  vDescription << ")";
245  else
246  vDescription << ",";
247  }
248 
249  vDescription << "->";
250  j = 0;
251 
252  HepMC::GenVertex::particles_out_const_iterator out, otmp;
253 
254  for (out = out_begin; out != out_end; out++, j++) {
255  if (!j)
256  vDescription << "(";
257 
258  HepPDT::ParticleID particleType((*out)->pdg_id());
259 
260  if (particleType.isValid()) {
261  pid = pdt_->particle(particleType);
262  if (pid)
263  vDescription << pid->name();
264  else
265  vDescription << (*out)->pdg_id();
266  } else
267  vDescription << (*out)->pdg_id();
268 
269  otmp = out;
270 
271  if (++otmp == out_end)
272  vDescription << ")";
273  else
274  vDescription << ",";
275  }
276 
277  return vDescription.str();
278 }
279 
This class traces the simulated and generated history of a given track.
Definition: VertexHistory.h:18
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 .
void beginRun(const edm::Run &, const edm::EventSetup &) override
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
edm::ESHandle< ParticleDataTable > pdt_
std::string particleString(int) const
Transition
Definition: Transition.h:12
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
Log< level::Warning, true > LogPrint
void endRun(const edm::Run &, const edm::EventSetup &) override
VertexHistoryAnalyzer(const edm::ParameterSet &)
SimParticleTrail const & simParticleTrail() const
Return all the simulated particle in the history.
Definition: HistoryBase.h:55
const edm::ESGetToken< ParticleDataTable, PDTRecord > pdtToken_
tuple config
parse the configuration file
VertexClassifier const & evaluate(reco::VertexBaseRef const &)
Classify the RecoVertex in categories.
const edm::EDGetTokenT< edm::View< reco::Vertex > > vtxToken_
size_type size() const
Size of the RefVector.
Definition: RefVector.h:102
virtual void newEvent(edm::Event const &, edm::EventSetup const &)
Pre-process event information (for accessing reconstraction information)
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:151
GenVertexTrail const & genVertexTrail() const
Return all generated vertex in the history.
Definition: HistoryBase.h:58
~VertexHistoryAnalyzer() override=default
std::vector< TrackingVertexRef > SimVertexTrail
SimVertex trail type.
Definition: HistoryBase.h:33
ESGetTokenH3DDVariant esConsumes(std::string const &Reccord, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
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