Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:31:06

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