CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
VertexHistoryAnalyzer Class Reference
Inheritance diagram for VertexHistoryAnalyzer:
edm::EDAnalyzer edm::EDConsumerBase

Public Member Functions

 VertexHistoryAnalyzer (const edm::ParameterSet &)
 
- Public Member Functions inherited from edm::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzer ()
 
ModuleDescription const & moduleDescription () const
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
 EDConsumerBase ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesDependentUpon (const std::string &iProcessName, std::vector< const char * > &oModuleLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Private Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &) override
 
virtual void beginJob () override
 
virtual void beginRun (const edm::Run &, const edm::EventSetup &) override
 
std::string particleString (int) const
 
std::string vertexString (const TrackingParticleRefVector &, const TrackingParticleRefVector &) const
 
std::string vertexString (HepMC::GenVertex::particles_in_const_iterator, HepMC::GenVertex::particles_in_const_iterator, HepMC::GenVertex::particles_out_const_iterator, HepMC::GenVertex::particles_out_const_iterator) const
 

Private Attributes

VertexClassifier classifier_
 
edm::ESHandle< ParticleDataTablepdt_
 
edm::InputTag vertexProducer_
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Definition at line 34 of file VertexHistoryAnalyzer.cc.

Constructor & Destructor Documentation

VertexHistoryAnalyzer::VertexHistoryAnalyzer ( const edm::ParameterSet config)
explicit

Definition at line 70 of file VertexHistoryAnalyzer.cc.

References edm::ParameterSet::getUntrackedParameter(), and vertexProducer_.

70  : classifier_(config,consumesCollector())
71 {
72  vertexProducer_ = config.getUntrackedParameter<edm::InputTag> ("vertexProducer");
73  consumes<edm::View<reco::Vertex>>(vertexProducer_);
74 }
T getUntrackedParameter(std::string const &, T const &) const
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.

Member Function Documentation

void VertexHistoryAnalyzer::analyze ( const edm::Event event,
const edm::EventSetup setup 
)
overrideprivatevirtual

Implements edm::EDAnalyzer.

Definition at line 77 of file VertexHistoryAnalyzer.cc.

References classifier_, gather_cfg::cout, VertexClassifier::evaluate(), VertexCategories::Fake, genParticleCandidates2GenParticles_cfi::genParticles, HistoryBase::genParticleTrail(), HistoryBase::genVertexTrail(), VertexClassifier::history(), cmsHarvester::index, VertexCategories::is(), VertexClassifier::newEvent(), particleString(), benchmark_cfg::pdgId, HistoryBase::simParticleTrail(), HistoryBase::simVertexTrail(), GoodVertex_cfg::vertexCollection, vertexProducer_, and vertexString().

78 {
79  // Set the classifier for a new event
80  classifier_.newEvent(event, setup);
81 
82  // Vertex collection
84  event.getByLabel(vertexProducer_, vertexCollection);
85 
86  // Get a constant reference to the track history associated to the classifier
87  VertexHistory const & tracer = classifier_.history();
88 
89  // Loop over the track collection.
90  for (std::size_t index = 0; index < vertexCollection->size(); index++)
91  {
92  std::cout << std::endl << "History for vertex #" << index << " : " << std::endl;
93 
94  // Classify the track and detect for fakes
96  {
97  // Get the list of TrackingParticles associated to
99 
100  // Loop over all simParticles
101  for (std::size_t hindex=0; hindex<simParticles.size(); hindex++)
102  {
103  std::cout << " simParticles [" << hindex << "] : "
104  << particleString(simParticles[hindex]->pdgId())
105  << std::endl;
106  }
107 
108  // Get the list of TrackingVertexes associated to
109  VertexHistory::SimVertexTrail simVertexes(tracer.simVertexTrail());
110 
111  // Loop over all simVertexes
112  if ( !simVertexes.empty() )
113  {
114  for (std::size_t hindex=0; hindex<simVertexes.size(); hindex++)
115  {
116  std::cout << " simVertex [" << hindex << "] : "
117  << vertexString(
118  simVertexes[hindex]->sourceTracks(),
119  simVertexes[hindex]->daughterTracks()
120  )
121  << std::endl;
122  }
123  }
124  else
125  std::cout << " simVertex no found" << std::endl;
126 
127  // Get the list of GenParticles associated to
129 
130  // Loop over all genParticles
131  for (std::size_t hindex=0; hindex<genParticles.size(); hindex++)
132  {
133  std::cout << " genParticles [" << hindex << "] : "
134  << particleString(genParticles[hindex]->pdg_id())
135  << std::endl;
136  }
137 
138  // Get the list of TrackingVertexes associated to
139  VertexHistory::GenVertexTrail genVertexes(tracer.genVertexTrail());
140 
141  // Loop over all simVertexes
142  if ( !genVertexes.empty() )
143  {
144  for (std::size_t hindex=0; hindex<genVertexes.size(); hindex++)
145  {
146  std::cout << " genVertex [" << hindex << "] : "
147  << vertexString(
148  genVertexes[hindex]->particles_in_const_begin(),
149  genVertexes[hindex]->particles_in_const_end(),
150  genVertexes[hindex]->particles_out_const_begin(),
151  genVertexes[hindex]->particles_out_const_end()
152  )
153  << std::endl;
154  }
155  }
156  else
157  std::cout << " genVertex no found" << std::endl;
158  }
159  else
160  std::cout << " fake vertex" << std::endl;
161 
162  std::cout << " vertex categories : " << classifier_;
163  std::cout << std::endl;
164  }
165 }
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.
std::string vertexString(const TrackingParticleRefVector &, const TrackingParticleRefVector &) const
SimVertexTrail const & simVertexTrail() const
Return all the simulated vertices in the history.
Definition: HistoryBase.h:53
tuple vertexCollection
std::vector< const HepMC::GenParticle * > GenParticleTrail
GenParticle trail type.
Definition: HistoryBase.h:18
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:21
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.
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
GenParticleTrail const & genParticleTrail() const
Return all generated particle in the history.
Definition: HistoryBase.h:71
std::vector< TrackingParticleRef > SimParticleTrail
SimParticle trail type.
Definition: HistoryBase.h:27
void VertexHistoryAnalyzer::beginJob ( void  )
overrideprivatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 176 of file VertexHistoryAnalyzer.cc.

177 {
178 
179 }
void VertexHistoryAnalyzer::beginRun ( const edm::Run run,
const edm::EventSetup setup 
)
overrideprivatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 169 of file VertexHistoryAnalyzer.cc.

References edm::EventSetup::getData(), and pdt_.

170 {
171  // Get the particles table.
172  setup.getData( pdt_ );
173 }
void getData(T &iHolder) const
Definition: EventSetup.h:79
edm::ESHandle< ParticleDataTable > pdt_
std::string VertexHistoryAnalyzer::particleString ( int  pdgId) const
private

Definition at line 182 of file VertexHistoryAnalyzer.cc.

References HLT_25ns14e33_v1_cff::particleType, benchmark_cfg::pdgId, pdt_, and sysUtil::pid.

Referenced by analyze().

183 {
184  ParticleData const * pid;
185 
186  std::ostringstream vDescription;
187 
188  HepPDT::ParticleID particleType(pdgId);
189 
190  if (particleType.isValid())
191  {
192  pid = pdt_->particle(particleType);
193  if (pid)
194  vDescription << pid->name();
195  else
196  vDescription << pdgId;
197  }
198  else
199  vDescription << pdgId;
200 
201  return vDescription.str();
202 }
edm::ESHandle< ParticleDataTable > pdt_
HepPDT::ParticleData ParticleData
tuple pid
Definition: sysUtil.py:22
std::string VertexHistoryAnalyzer::vertexString ( const TrackingParticleRefVector in,
const TrackingParticleRefVector out 
) const
private

Definition at line 205 of file VertexHistoryAnalyzer.cc.

References j, HLT_25ns14e33_v1_cff::particleType, benchmark_cfg::pdgId, pdt_, sysUtil::pid, and edm::RefVector< C, T, F >::size().

Referenced by analyze().

209 {
210  ParticleData const * pid;
211 
212  std::ostringstream vDescription;
213 
214  for (std::size_t j = 0; j < in.size(); j++)
215  {
216  if (!j) vDescription << "(";
217 
218  HepPDT::ParticleID particleType(in[j]->pdgId());
219 
220  if (particleType.isValid())
221  {
222  pid = pdt_->particle(particleType);
223  if (pid)
224  vDescription << pid->name();
225  else
226  vDescription << in[j]->pdgId();
227  }
228  else
229  vDescription << in[j]->pdgId();
230 
231  if (j == in.size() - 1) vDescription << ")";
232  else vDescription << ",";
233  }
234 
235  vDescription << "->";
236 
237  for (std::size_t j = 0; j < out.size(); j++)
238  {
239  if (!j) vDescription << "(";
240 
241  HepPDT::ParticleID particleType(out[j]->pdgId());
242 
243  if (particleType.isValid())
244  {
245  pid = pdt_->particle(particleType);
246  if (pid)
247  vDescription << pid->name();
248  else
249  vDescription << out[j]->pdgId();
250  }
251  else
252  vDescription << out[j]->pdgId();
253 
254  if (j == out.size() - 1) vDescription << ")";
255  else vDescription << ",";
256  }
257 
258  return vDescription.str();
259 }
edm::ESHandle< ParticleDataTable > pdt_
int j
Definition: DBlmapReader.cc:9
HepPDT::ParticleData ParticleData
tuple pid
Definition: sysUtil.py:22
size_type size() const
Size of the RefVector.
Definition: RefVector.h:99
std::string VertexHistoryAnalyzer::vertexString ( HepMC::GenVertex::particles_in_const_iterator  in_begin,
HepMC::GenVertex::particles_in_const_iterator  in_end,
HepMC::GenVertex::particles_out_const_iterator  out_begin,
HepMC::GenVertex::particles_out_const_iterator  out_end 
) const
private

Definition at line 262 of file VertexHistoryAnalyzer.cc.

References recoMuon::in, j, dbtoconf::out, HLT_25ns14e33_v1_cff::particleType, pdt_, and sysUtil::pid.

268 {
269  ParticleData const * pid;
270 
271  std::ostringstream vDescription;
272 
273  std::size_t j = 0;
274 
275  HepMC::GenVertex::particles_in_const_iterator in, itmp;
276 
277  for (in = in_begin; in != in_end; in++, j++)
278  {
279  if (!j) vDescription << "(";
280 
281  HepPDT::ParticleID particleType((*in)->pdg_id());
282 
283  if (particleType.isValid())
284  {
285  pid = pdt_->particle(particleType);
286  if (pid)
287  vDescription << pid->name();
288  else
289  vDescription << (*in)->pdg_id();
290  }
291  else
292  vDescription << (*in)->pdg_id();
293 
294  itmp = in;
295 
296  if (++itmp == in_end) vDescription << ")";
297  else vDescription << ",";
298  }
299 
300  vDescription << "->";
301  j = 0;
302 
303  HepMC::GenVertex::particles_out_const_iterator out, otmp;
304 
305  for (out = out_begin; out != out_end; out++, j++)
306  {
307  if (!j) vDescription << "(";
308 
309  HepPDT::ParticleID particleType((*out)->pdg_id());
310 
311  if (particleType.isValid())
312  {
313  pid = pdt_->particle(particleType);
314  if (pid)
315  vDescription << pid->name();
316  else
317  vDescription << (*out)->pdg_id();
318  }
319  else
320  vDescription << (*out)->pdg_id();
321 
322  otmp = out;
323 
324  if (++otmp == out_end) vDescription << ")";
325  else vDescription << ",";
326  }
327 
328  return vDescription.str();
329 }
edm::ESHandle< ParticleDataTable > pdt_
int j
Definition: DBlmapReader.cc:9
HepPDT::ParticleData ParticleData
tuple out
Definition: dbtoconf.py:99
tuple pid
Definition: sysUtil.py:22

Member Data Documentation

VertexClassifier VertexHistoryAnalyzer::classifier_
private

Definition at line 48 of file VertexHistoryAnalyzer.cc.

Referenced by analyze().

edm::ESHandle<ParticleDataTable> VertexHistoryAnalyzer::pdt_
private

Definition at line 52 of file VertexHistoryAnalyzer.cc.

Referenced by beginRun(), particleString(), and vertexString().

edm::InputTag VertexHistoryAnalyzer::vertexProducer_
private

Definition at line 50 of file VertexHistoryAnalyzer.cc.

Referenced by analyze(), and VertexHistoryAnalyzer().