CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
TrackHistoryAnalyzer Class Reference
Inheritance diagram for TrackHistoryAnalyzer:
edm::one::EDAnalyzer<> edm::one::EDAnalyzerBase edm::EDConsumerBase

Public Member Functions

 TrackHistoryAnalyzer (const edm::ParameterSet &)
 
 ~TrackHistoryAnalyzer () override
 
- Public Member Functions inherited from edm::one::EDAnalyzer<>
 EDAnalyzer ()=default
 
SerialTaskQueueglobalLuminosityBlocksQueue () final
 
SerialTaskQueueglobalRunsQueue () final
 
bool wantsGlobalLuminosityBlocks () const final
 
bool wantsGlobalRuns () const final
 
- Public Member Functions inherited from edm::one::EDAnalyzerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzerBase ()
 
ModuleDescription const & moduleDescription () const
 
bool wantsStreamLuminosityBlocks () const
 
bool wantsStreamRuns () const
 
 ~EDAnalyzerBase () override
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
virtual ~EDConsumerBase () noexcept(false)
 

Private Member Functions

void analyze (const edm::Event &, const edm::EventSetup &) override
 
void beginJob () override
 
virtual void beginRun (const edm::Run &, const edm::EventSetup &)
 
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::one::EDAnalyzerBase
typedef EDAnalyzerBase ModuleType
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::one::EDAnalyzerBase
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- 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,consumesCollector())
75 {
76  trackProducer_ = config.getUntrackedParameter<edm::InputTag> ( "trackProducer" );
77  consumes<edm::View<reco::Track>>(trackProducer_);
78 }
T getUntrackedParameter(std::string const &, T const &) const
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
TrackHistoryAnalyzer::~TrackHistoryAnalyzer ( )
override

Definition at line 81 of file TrackHistoryAnalyzer.cc.

81 { }

Member Function Documentation

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

Definition at line 84 of file TrackHistoryAnalyzer.cc.

References classifier_, gather_cfg::cout, TrackClassifier::evaluate(), TrackCategories::Fake, GenHFHadronMatcher_cfi::genParticles, HistoryBase::genParticleTrail(), HistoryBase::genVertexTrail(), TrackClassifier::history(), TrackCategories::is(), TrackClassifier::newEvent(), particleString(), HiggsValidation_cfi::pdg_id, common_cff::pdgId, HistoryBase::simParticleTrail(), HistoryBase::simVertexTrail(), findElectronsInSiStrips_cfi::trackCollection, trackProducer_, and vertexString().

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

Reimplemented from edm::one::EDAnalyzerBase.

Definition at line 183 of file TrackHistoryAnalyzer.cc.

References totalTracks_.

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

Definition at line 176 of file TrackHistoryAnalyzer.cc.

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

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

Definition at line 189 of file TrackHistoryAnalyzer.cc.

References source_particleGun_cfi::ParticleID, objects.autophobj::particleType, common_cff::pdgId, pdt_, and sysUtil::pid.

Referenced by analyze().

190 {
191  ParticleData const * pid;
192 
193  std::ostringstream vDescription;
194 
196 
197  if (particleType.isValid())
198  {
199  pid = pdt_->particle(particleType);
200  if (pid)
201  vDescription << pid->name();
202  else
203  vDescription << pdgId;
204  }
205  else
206  vDescription << pdgId;
207 
208  return vDescription.str();
209 }
HepPDT::ParticleData ParticleData
edm::ESHandle< ParticleDataTable > pdt_
std::string TrackHistoryAnalyzer::vertexString ( const TrackingParticleRefVector in,
const TrackingParticleRefVector out 
) const
private

Definition at line 212 of file TrackHistoryAnalyzer.cc.

References source_particleGun_cfi::ParticleID, objects.autophobj::particleType, common_cff::pdgId, pdt_, sysUtil::pid, and edm::RefVector< C, T, F >::size().

Referenced by analyze().

216 {
217  ParticleData const * pid;
218 
219  std::ostringstream vDescription;
220 
221  for (std::size_t j = 0; j < in.size(); j++)
222  {
223  if (!j) vDescription << "(";
224 
226 
227  if (particleType.isValid())
228  {
229  pid = pdt_->particle(particleType);
230  if (pid)
231  vDescription << pid->name();
232  else
233  vDescription << in[j]->pdgId();
234  }
235  else
236  vDescription << in[j]->pdgId();
237 
238  if (j == in.size() - 1) vDescription << ")";
239  else vDescription << ",";
240  }
241 
242  vDescription << "->";
243 
244  for (std::size_t j = 0; j < out.size(); j++)
245  {
246  if (!j) vDescription << "(";
247 
249 
250  if (particleType.isValid())
251  {
252  pid = pdt_->particle(particleType);
253  if (pid)
254  vDescription << pid->name();
255  else
256  vDescription << out[j]->pdgId();
257  }
258  else
259  vDescription << out[j]->pdgId();
260 
261  if (j == out.size() - 1) vDescription << ")";
262  else vDescription << ",";
263  }
264 
265  return vDescription.str();
266 }
HepPDT::ParticleData ParticleData
edm::ESHandle< ParticleDataTable > pdt_
size_type size() const
Size of the RefVector.
Definition: RefVector.h:107
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 269 of file TrackHistoryAnalyzer.cc.

References DEFINE_FWK_MODULE, recoMuon::in, MillePedeFileConverter_cfg::out, source_particleGun_cfi::ParticleID, objects.autophobj::particleType, pdt_, and sysUtil::pid.

275 {
276  ParticleData const * pid;
277 
278  std::ostringstream vDescription;
279 
280  std::size_t j = 0;
281 
282  HepMC::GenVertex::particles_in_const_iterator in, itmp;
283 
284  for (in = in_begin; in != in_end; in++, j++)
285  {
286  if (!j) vDescription << "(";
287 
288  HepPDT::ParticleID particleType((*in)->pdg_id());
289 
290  if (particleType.isValid())
291  {
292  pid = pdt_->particle(particleType);
293  if (pid)
294  vDescription << pid->name();
295  else
296  vDescription << (*in)->pdg_id();
297  }
298  else
299  vDescription << (*in)->pdg_id();
300 
301  itmp = in;
302 
303  if (++itmp == in_end) vDescription << ")";
304  else vDescription << ",";
305  }
306 
307  vDescription << "->";
308  j = 0;
309 
310  HepMC::GenVertex::particles_out_const_iterator out, otmp;
311 
312  for (out = out_begin; out != out_end; out++, j++)
313  {
314  if (!j) vDescription << "(";
315 
316  HepPDT::ParticleID particleType((*out)->pdg_id());
317 
318  if (particleType.isValid())
319  {
320  pid = pdt_->particle(particleType);
321  if (pid)
322  vDescription << pid->name();
323  else
324  vDescription << (*out)->pdg_id();
325  }
326  else
327  vDescription << (*out)->pdg_id();
328 
329  otmp = out;
330 
331  if (++otmp == out_end) vDescription << ")";
332  else vDescription << ",";
333  }
334 
335  return vDescription.str();
336 }
HepPDT::ParticleData ParticleData
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().