File indexing completed on 2024-09-10 02:59:14
0001
0002
0003
0004
0005
0006
0007
0008
0009 #include <iostream>
0010 #include <memory>
0011 #include <sstream>
0012 #include <string>
0013 #include <vector>
0014
0015
0016
0017 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0018 #include "DataFormats/VertexReco/interface/Vertex.h"
0019
0020 #include "FWCore/Framework/interface/ESHandle.h"
0021 #include "FWCore/Framework/interface/Event.h"
0022 #include "FWCore/Framework/interface/EventSetup.h"
0023 #include "FWCore/Framework/interface/Frameworkfwd.h"
0024 #include "FWCore/Framework/interface/MakerMacros.h"
0025 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0026 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0027
0028 #include "SimGeneral/HepPDTRecord/interface/ParticleDataTable.h"
0029 #include "SimTracker/TrackHistory/interface/TrackClassifier.h"
0030
0031
0032
0033
0034
0035 class TrackHistoryAnalyzer : public edm::one::EDAnalyzer<edm::one::WatchRuns> {
0036 public:
0037 explicit TrackHistoryAnalyzer(const edm::ParameterSet &);
0038 ~TrackHistoryAnalyzer() override = default;
0039
0040 private:
0041 void beginJob() override;
0042 void beginRun(const edm::Run &, const edm::EventSetup &) override;
0043 void analyze(const edm::Event &, const edm::EventSetup &) override;
0044 void endRun(const edm::Run &, const edm::EventSetup &) override {}
0045
0046
0047 const edm::ESGetToken<ParticleDataTable, PDTRecord> pdtToken_;
0048 const edm::EDGetTokenT<edm::View<reco::Track>> trkToken_;
0049
0050 std::size_t totalTracks_;
0051
0052 edm::ESHandle<ParticleDataTable> pdt_;
0053
0054 std::string particleString(int) const;
0055
0056 TrackClassifier classifier_;
0057
0058 std::string vertexString(const TrackingParticleRefVector &, const TrackingParticleRefVector &) const;
0059
0060 std::string vertexString(HepMC::GenVertex::particles_in_const_iterator,
0061 HepMC::GenVertex::particles_in_const_iterator,
0062 HepMC::GenVertex::particles_out_const_iterator,
0063 HepMC::GenVertex::particles_out_const_iterator) const;
0064 };
0065
0066 TrackHistoryAnalyzer::TrackHistoryAnalyzer(const edm::ParameterSet &config)
0067 : pdtToken_(esConsumes<edm::Transition::BeginRun>()),
0068 trkToken_(consumes<edm::View<reco::Track>>(config.getUntrackedParameter<edm::InputTag>("trackProducer"))),
0069 classifier_(config, consumesCollector()) {}
0070
0071 void TrackHistoryAnalyzer::analyze(const edm::Event &event, const edm::EventSetup &setup) {
0072
0073 edm::Handle<edm::View<reco::Track>> trackCollection;
0074 event.getByToken(trkToken_, trackCollection);
0075
0076
0077 classifier_.newEvent(event, setup);
0078
0079
0080 TrackHistory const &tracer = classifier_.history();
0081
0082
0083 for (std::size_t index = 0; index < trackCollection->size(); index++) {
0084 edm::LogPrint("TrackHistoryAnalyzer") << std::endl << "History for track #" << index << " : ";
0085
0086
0087 if (!classifier_.evaluate(reco::TrackBaseRef(trackCollection, index)).is(TrackClassifier::Fake)) {
0088
0089 TrackHistory::SimParticleTrail simParticles(tracer.simParticleTrail());
0090
0091
0092 for (std::size_t hindex = 0; hindex < simParticles.size(); hindex++) {
0093 edm::LogPrint("TrackHistoryAnalyzer")
0094 << " simParticles [" << hindex << "] : " << particleString(simParticles[hindex]->pdgId());
0095 }
0096
0097
0098 TrackHistory::SimVertexTrail simVertexes(tracer.simVertexTrail());
0099
0100
0101 if (!simVertexes.empty()) {
0102 for (std::size_t hindex = 0; hindex < simVertexes.size(); hindex++) {
0103 edm::LogPrint("TrackHistoryAnalyzer")
0104 << " simVertex [" << hindex
0105 << "] : " << vertexString(simVertexes[hindex]->sourceTracks(), simVertexes[hindex]->daughterTracks());
0106 }
0107 } else
0108 edm::LogPrint("TrackHistoryAnalyzer") << " simVertex no found";
0109
0110
0111 TrackHistory::GenParticleTrail genParticles(tracer.genParticleTrail());
0112
0113
0114 for (std::size_t hindex = 0; hindex < genParticles.size(); hindex++) {
0115 edm::LogPrint("TrackHistoryAnalyzer")
0116 << " genParticles [" << hindex << "] : " << particleString(genParticles[hindex]->pdg_id());
0117 }
0118
0119
0120 TrackHistory::GenVertexTrail genVertexes(tracer.genVertexTrail());
0121
0122
0123 if (!genVertexes.empty()) {
0124 for (std::size_t hindex = 0; hindex < genVertexes.size(); hindex++) {
0125 edm::LogPrint("TrackHistoryAnalyzer") << " genVertex [" << hindex << "] : "
0126 << vertexString(genVertexes[hindex]->particles_in_const_begin(),
0127 genVertexes[hindex]->particles_in_const_end(),
0128 genVertexes[hindex]->particles_out_const_begin(),
0129 genVertexes[hindex]->particles_out_const_end());
0130 }
0131 } else
0132 edm::LogPrint("TrackHistoryAnalyzer") << " genVertex no found";
0133 } else
0134 edm::LogPrint("TrackHistoryAnalyzer") << " fake track";
0135
0136 edm::LogPrint("TrackHistoryAnalyzer") << " track categories : " << classifier_;
0137 }
0138 }
0139
0140 void TrackHistoryAnalyzer::beginRun(const edm::Run &run, const edm::EventSetup &setup) {
0141
0142 pdt_ = setup.getHandle(pdtToken_);
0143 }
0144
0145 void TrackHistoryAnalyzer::beginJob() { totalTracks_ = 0; }
0146
0147 std::string TrackHistoryAnalyzer::particleString(int pdgId) const {
0148 ParticleData const *pid;
0149
0150 std::ostringstream vDescription;
0151
0152 HepPDT::ParticleID particleType(pdgId);
0153
0154 if (particleType.isValid()) {
0155 pid = pdt_->particle(particleType);
0156 if (pid)
0157 vDescription << pid->name();
0158 else
0159 vDescription << pdgId;
0160 } else
0161 vDescription << pdgId;
0162
0163 return vDescription.str();
0164 }
0165
0166 std::string TrackHistoryAnalyzer::vertexString(const TrackingParticleRefVector &in,
0167 const TrackingParticleRefVector &out) const {
0168 ParticleData const *pid;
0169
0170 std::ostringstream vDescription;
0171
0172 for (std::size_t j = 0; j < in.size(); j++) {
0173 if (!j)
0174 vDescription << "(";
0175
0176 HepPDT::ParticleID particleType(in[j]->pdgId());
0177
0178 if (particleType.isValid()) {
0179 pid = pdt_->particle(particleType);
0180 if (pid)
0181 vDescription << pid->name();
0182 else
0183 vDescription << in[j]->pdgId();
0184 } else
0185 vDescription << in[j]->pdgId();
0186
0187 if (j == in.size() - 1)
0188 vDescription << ")";
0189 else
0190 vDescription << ",";
0191 }
0192
0193 vDescription << "->";
0194
0195 for (std::size_t j = 0; j < out.size(); j++) {
0196 if (!j)
0197 vDescription << "(";
0198
0199 HepPDT::ParticleID particleType(out[j]->pdgId());
0200
0201 if (particleType.isValid()) {
0202 pid = pdt_->particle(particleType);
0203 if (pid)
0204 vDescription << pid->name();
0205 else
0206 vDescription << out[j]->pdgId();
0207 } else
0208 vDescription << out[j]->pdgId();
0209
0210 if (j == out.size() - 1)
0211 vDescription << ")";
0212 else
0213 vDescription << ",";
0214 }
0215
0216 return vDescription.str();
0217 }
0218
0219 std::string TrackHistoryAnalyzer::vertexString(HepMC::GenVertex::particles_in_const_iterator in_begin,
0220 HepMC::GenVertex::particles_in_const_iterator in_end,
0221 HepMC::GenVertex::particles_out_const_iterator out_begin,
0222 HepMC::GenVertex::particles_out_const_iterator out_end) const {
0223 ParticleData const *pid;
0224
0225 std::ostringstream vDescription;
0226
0227 std::size_t j = 0;
0228
0229 HepMC::GenVertex::particles_in_const_iterator in, itmp;
0230
0231 for (in = in_begin; in != in_end; in++, j++) {
0232 if (!j)
0233 vDescription << "(";
0234
0235 HepPDT::ParticleID particleType((*in)->pdg_id());
0236
0237 if (particleType.isValid()) {
0238 pid = pdt_->particle(particleType);
0239 if (pid)
0240 vDescription << pid->name();
0241 else
0242 vDescription << (*in)->pdg_id();
0243 } else
0244 vDescription << (*in)->pdg_id();
0245
0246 itmp = in;
0247
0248 if (++itmp == in_end)
0249 vDescription << ")";
0250 else
0251 vDescription << ",";
0252 }
0253
0254 vDescription << "->";
0255 j = 0;
0256
0257 HepMC::GenVertex::particles_out_const_iterator out, otmp;
0258
0259 for (out = out_begin; out != out_end; out++, j++) {
0260 if (!j)
0261 vDescription << "(";
0262
0263 HepPDT::ParticleID particleType((*out)->pdg_id());
0264
0265 if (particleType.isValid()) {
0266 pid = pdt_->particle(particleType);
0267 if (pid)
0268 vDescription << pid->name();
0269 else
0270 vDescription << (*out)->pdg_id();
0271 } else
0272 vDescription << (*out)->pdg_id();
0273
0274 otmp = out;
0275
0276 if (++otmp == out_end)
0277 vDescription << ")";
0278 else
0279 vDescription << ",";
0280 }
0281
0282 return vDescription.str();
0283 }
0284
0285 DEFINE_FWK_MODULE(TrackHistoryAnalyzer);