Public Member Functions | |
TrackHistoryAnalyzer (const edm::ParameterSet &) | |
~TrackHistoryAnalyzer () | |
Private Member Functions | |
virtual void | analyze (const edm::Event &, const edm::EventSetup &) |
virtual void | beginJob () |
virtual void | beginRun (const edm::Run &, const edm::EventSetup &) |
std::string | particleString (int) 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 |
std::string | vertexString (TrackingParticleRefVector, TrackingParticleRefVector) const |
Private Attributes | |
TrackClassifier | classifier_ |
edm::ESHandle< ParticleDataTable > | pdt_ |
std::size_t | totalTracks_ |
edm::InputTag | trackProducer_ |
Definition at line 35 of file TrackHistoryAnalyzer.cc.
TrackHistoryAnalyzer::TrackHistoryAnalyzer | ( | const edm::ParameterSet & | config | ) | [explicit] |
Definition at line 74 of file TrackHistoryAnalyzer.cc.
References edm::ParameterSet::getUntrackedParameter(), and trackProducer_.
: classifier_(config) { trackProducer_ = config.getUntrackedParameter<edm::InputTag> ( "trackProducer" ); }
TrackHistoryAnalyzer::~TrackHistoryAnalyzer | ( | ) |
Definition at line 80 of file TrackHistoryAnalyzer.cc.
{ }
void TrackHistoryAnalyzer::analyze | ( | const edm::Event & | event, |
const edm::EventSetup & | setup | ||
) | [private, virtual] |
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(), getHLTprescales::index, TrackCategories::is(), TrackClassifier::newEvent(), particleString(), benchmark_cfg::pdgId, HistoryBase::simParticleTrail(), HistoryBase::simVertexTrail(), trackProducer_, and vertexString().
{ // Track collection edm::Handle<edm::View<reco::Track> > trackCollection; event.getByLabel(trackProducer_, trackCollection); // Set the classifier for a new event classifier_.newEvent(event, setup); // Get a constant reference to the track history associated to the classifier TrackHistory const & tracer = classifier_.history(); // Loop over the track collection. for (std::size_t index = 0; index < trackCollection->size(); index++) { std::cout << std::endl << "History for track #" << index << " : " << std::endl; // Classify the track and detect for fakes if ( ! classifier_.evaluate( reco::TrackBaseRef(trackCollection, index) ).is(TrackClassifier::Fake) ) { // Get the list of TrackingParticles associated to TrackHistory::SimParticleTrail simParticles(tracer.simParticleTrail()); // Loop over all simParticles for (std::size_t hindex=0; hindex<simParticles.size(); hindex++) { std::cout << " simParticles [" << hindex << "] : " << particleString(simParticles[hindex]->pdgId()) << std::endl; } // Get the list of TrackingVertexes associated to TrackHistory::SimVertexTrail simVertexes(tracer.simVertexTrail()); // Loop over all simVertexes if ( !simVertexes.empty() ) { for (std::size_t hindex=0; hindex<simVertexes.size(); hindex++) { std::cout << " simVertex [" << hindex << "] : " << vertexString( simVertexes[hindex]->sourceTracks(), simVertexes[hindex]->daughterTracks() ) << std::endl; } } else std::cout << " simVertex no found" << std::endl; // Get the list of GenParticles associated to TrackHistory::GenParticleTrail genParticles(tracer.genParticleTrail()); // Loop over all genParticles for (std::size_t hindex=0; hindex<genParticles.size(); hindex++) { std::cout << " genParticles [" << hindex << "] : " << particleString(genParticles[hindex]->pdg_id()) << std::endl; } // Get the list of TrackingVertexes associated to TrackHistory::GenVertexTrail genVertexes(tracer.genVertexTrail()); // Loop over all simVertexes if ( !genVertexes.empty() ) { for (std::size_t hindex=0; hindex<genVertexes.size(); hindex++) { std::cout << " genVertex [" << hindex << "] : " << vertexString( genVertexes[hindex]->particles_in_const_begin(), genVertexes[hindex]->particles_in_const_end(), genVertexes[hindex]->particles_out_const_begin(), genVertexes[hindex]->particles_out_const_end() ) << std::endl; } } else std::cout << " genVertex no found" << std::endl; } else std::cout << " fake track" << std::endl; std::cout << " track categories : " << classifier_; std::cout << std::endl; } }
void TrackHistoryAnalyzer::beginJob | ( | void | ) | [private, virtual] |
Reimplemented from edm::EDAnalyzer.
Definition at line 182 of file TrackHistoryAnalyzer.cc.
References totalTracks_.
{ totalTracks_ = 0; }
void TrackHistoryAnalyzer::beginRun | ( | const edm::Run & | run, |
const edm::EventSetup & | setup | ||
) | [private, virtual] |
Reimplemented from edm::EDAnalyzer.
Definition at line 175 of file TrackHistoryAnalyzer.cc.
References edm::EventSetup::getData(), and 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 evf::utils::pid.
Referenced by analyze().
{ ParticleData const * pid; std::ostringstream vDescription; HepPDT::ParticleID particleType(pdgId); if (particleType.isValid()) { pid = pdt_->particle(particleType); if (pid) vDescription << pid->name(); else vDescription << pdgId; } else vDescription << pdgId; return vDescription.str(); }
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 evf::utils::pid.
{ ParticleData const * pid; std::ostringstream vDescription; std::size_t j = 0; HepMC::GenVertex::particles_in_const_iterator in, itmp; for (in = in_begin; in != in_end; in++, j++) { if (!j) vDescription << "("; HepPDT::ParticleID particleType((*in)->pdg_id()); if (particleType.isValid()) { pid = pdt_->particle(particleType); if (pid) vDescription << pid->name(); else vDescription << (*in)->pdg_id(); } else vDescription << (*in)->pdg_id(); itmp = in; if (++itmp == in_end) vDescription << ")"; else vDescription << ","; } vDescription << "->"; j = 0; HepMC::GenVertex::particles_out_const_iterator out, otmp; for (out = out_begin; out != out_end; out++, j++) { if (!j) vDescription << "("; HepPDT::ParticleID particleType((*out)->pdg_id()); if (particleType.isValid()) { pid = pdt_->particle(particleType); if (pid) vDescription << pid->name(); else vDescription << (*out)->pdg_id(); } else vDescription << (*out)->pdg_id(); otmp = out; if (++otmp == out_end) vDescription << ")"; else vDescription << ","; } return vDescription.str(); }
std::string TrackHistoryAnalyzer::vertexString | ( | TrackingParticleRefVector | in, |
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_, evf::utils::pid, and edm::RefVector< C, T, F >::size().
Referenced by analyze().
{ ParticleData const * pid; std::ostringstream vDescription; for (std::size_t j = 0; j < in.size(); j++) { if (!j) vDescription << "("; HepPDT::ParticleID particleType(in[j]->pdgId()); if (particleType.isValid()) { pid = pdt_->particle(particleType); if (pid) vDescription << pid->name(); else vDescription << in[j]->pdgId(); } else vDescription << in[j]->pdgId(); if (j == in.size() - 1) vDescription << ")"; else vDescription << ","; } vDescription << "->"; for (std::size_t j = 0; j < out.size(); j++) { if (!j) vDescription << "("; HepPDT::ParticleID particleType(out[j]->pdgId()); if (particleType.isValid()) { pid = pdt_->particle(particleType); if (pid) vDescription << pid->name(); else vDescription << out[j]->pdgId(); } else vDescription << out[j]->pdgId(); if (j == out.size() - 1) vDescription << ")"; else vDescription << ","; } return vDescription.str(); }
Definition at line 58 of file TrackHistoryAnalyzer.cc.
Referenced by analyze().
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().
Definition at line 50 of file TrackHistoryAnalyzer.cc.
Referenced by analyze(), and TrackHistoryAnalyzer().