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
TrackHistoryAnalyzer Class Reference
Inheritance diagram for TrackHistoryAnalyzer:
edm::EDAnalyzer edm::EDConsumerBase

Public Member Functions

 TrackHistoryAnalyzer (const edm::ParameterSet &)
 
 ~TrackHistoryAnalyzer ()
 
- 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
 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
 
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

TrackClassifier classifier_
 
edm::ESHandle< ParticleDataTablepdt_
 
std::size_t totalTracks_
 
edm::InputTag trackProducer_
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
- 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 35 of file TrackHistoryAnalyzer.cc.

Constructor & Destructor Documentation

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

Definition at line 74 of file TrackHistoryAnalyzer.cc.

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

74  : classifier_(config)
75 {
76  trackProducer_ = config.getUntrackedParameter<edm::InputTag> ( "trackProducer" );
77 }
T getUntrackedParameter(std::string const &, T const &) const
TrackHistoryAnalyzer::~TrackHistoryAnalyzer ( )

Definition at line 80 of file TrackHistoryAnalyzer.cc.

80 { }

Member Function Documentation

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

Implements edm::EDAnalyzer.

Definition at line 83 of file TrackHistoryAnalyzer.cc.

References classifier_, gather_cfg::cout, TrackClassifier::evaluate(), TrackCategories::Fake, genParticleCandidates2GenParticles_cfi::genParticles, HistoryBase::genParticleTrail(), HistoryBase::genVertexTrail(), TrackClassifier::history(), cmsHarvester::index, TrackCategories::is(), TrackClassifier::newEvent(), particleString(), benchmark_cfg::pdgId, HistoryBase::simParticleTrail(), HistoryBase::simVertexTrail(), trackProducer_, and vertexString().

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 }
std::string vertexString(const TrackingParticleRefVector &, const TrackingParticleRefVector &) const
SimVertexTrail const & simVertexTrail() const
Return all the simulated vertices in the history.
Definition: HistoryBase.h:53
std::string particleString(int) const
TrackHistory const & history() const
Returns a reference to the track history used in the classification.
std::vector< const HepMC::GenParticle * > GenParticleTrail
GenParticle trail type.
Definition: HistoryBase.h:18
TrackClassifier const & evaluate(reco::TrackBaseRef const &)
Classify the RecoTrack in categories.
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
void newEvent(edm::Event const &, edm::EventSetup const &)
Pre-process event information (for accessing reconstraction information)
bool is(Category category) const
Returns track flag for a given category.
SimParticleTrail const & simParticleTrail() const
Return all the simulated particle in the history.
Definition: HistoryBase.h:59
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
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 TrackHistoryAnalyzer::beginJob ( void  )
overrideprivatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 182 of file TrackHistoryAnalyzer.cc.

References totalTracks_.

183 {
184  totalTracks_ = 0;
185 }
void TrackHistoryAnalyzer::beginRun ( const edm::Run run,
const edm::EventSetup setup 
)
overrideprivatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 175 of file TrackHistoryAnalyzer.cc.

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

176 {
177  // Get the particles table.
178  setup.getData( pdt_ );
179 }
void getData(T &iHolder) const
Definition: EventSetup.h:67
edm::ESHandle< ParticleDataTable > pdt_
std::string TrackHistoryAnalyzer::particleString ( int  pdgId) const
private

Definition at line 188 of file TrackHistoryAnalyzer.cc.

References RecoTau_DiTaus_pt_20-420_cfg::ParticleID, benchmark_cfg::pdgId, pdt_, and sysUtil::pid.

Referenced by analyze().

189 {
190  ParticleData const * pid;
191 
192  std::ostringstream vDescription;
193 
194  HepPDT::ParticleID particleType(pdgId);
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 }
HepPDT::ParticleData ParticleData
tuple pid
Definition: sysUtil.py:22
edm::ESHandle< ParticleDataTable > pdt_
std::string TrackHistoryAnalyzer::vertexString ( const TrackingParticleRefVector in,
const TrackingParticleRefVector out 
) const
private

Definition at line 211 of file TrackHistoryAnalyzer.cc.

References j, RecoTau_DiTaus_pt_20-420_cfg::ParticleID, benchmark_cfg::pdgId, pdt_, sysUtil::pid, and edm::RefVector< C, T, F >::size().

Referenced by analyze().

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 
224  HepPDT::ParticleID particleType(in[j]->pdgId());
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 
247  HepPDT::ParticleID particleType(out[j]->pdgId());
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 }
int j
Definition: DBlmapReader.cc:9
HepPDT::ParticleData ParticleData
tuple pid
Definition: sysUtil.py:22
edm::ESHandle< ParticleDataTable > pdt_
size_type size() const
Size of the RefVector.
Definition: RefVector.h:89
std::string TrackHistoryAnalyzer::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 268 of file TrackHistoryAnalyzer.cc.

References recoMuon::in, j, dbtoconf::out, RecoTau_DiTaus_pt_20-420_cfg::ParticleID, pdt_, and sysUtil::pid.

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 }
int j
Definition: DBlmapReader.cc:9
HepPDT::ParticleData ParticleData
tuple out
Definition: dbtoconf.py:99
tuple pid
Definition: sysUtil.py:22
edm::ESHandle< ParticleDataTable > pdt_

Member Data Documentation

TrackClassifier TrackHistoryAnalyzer::classifier_
private

Definition at line 58 of file TrackHistoryAnalyzer.cc.

Referenced by analyze().

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

Definition at line 54 of file TrackHistoryAnalyzer.cc.

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

std::size_t TrackHistoryAnalyzer::totalTracks_
private

Definition at line 52 of file TrackHistoryAnalyzer.cc.

Referenced by beginJob().

edm::InputTag TrackHistoryAnalyzer::trackProducer_
private

Definition at line 50 of file TrackHistoryAnalyzer.cc.

Referenced by analyze(), and TrackHistoryAnalyzer().