CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch9/src/SimTracker/TrackHistory/plugins/TrackHistoryAnalyzer.cc

Go to the documentation of this file.
00001 /*
00002  *  TrackHistoryAnalyzer.C
00003  *
00004  *  Created by Victor Eduardo Bazterra on 5/31/07.
00005  *
00006  */
00007 
00008 // system include files
00009 #include <iostream>
00010 #include <memory>
00011 #include <string>
00012 #include <sstream>
00013 #include <vector>
00014 
00015 // user include files
00016 
00017 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00018 #include "DataFormats/VertexReco/interface/Vertex.h"
00019 
00020 #include "FWCore/Framework/interface/Frameworkfwd.h"
00021 #include "FWCore/Framework/interface/EDAnalyzer.h"
00022 #include "FWCore/Framework/interface/Event.h"
00023 #include "FWCore/Framework/interface/ESHandle.h"
00024 #include "FWCore/Framework/interface/EventSetup.h"
00025 #include "FWCore/Framework/interface/MakerMacros.h"
00026 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00027 
00028 #include "SimGeneral/HepPDTRecord/interface/ParticleDataTable.h"
00029 #include "SimTracker/TrackHistory/interface/TrackClassifier.h"
00030 
00031 //
00032 // class decleration
00033 //
00034 
00035 class TrackHistoryAnalyzer : public edm::EDAnalyzer
00036 {
00037 public:
00038 
00039     explicit TrackHistoryAnalyzer(const edm::ParameterSet&);
00040     ~TrackHistoryAnalyzer();
00041 
00042 private:
00043 
00044   virtual void beginRun(const edm::Run&, const edm::EventSetup&);
00045     virtual void beginJob() ;
00046     virtual void analyze(const edm::Event&, const edm::EventSetup&);
00047 
00048     // Member data
00049 
00050     edm::InputTag trackProducer_;
00051 
00052     std::size_t totalTracks_;
00053 
00054     edm::ESHandle<ParticleDataTable> pdt_;
00055 
00056     std::string particleString(int) const;
00057 
00058     TrackClassifier classifier_;
00059 
00060     std::string vertexString(
00061         TrackingParticleRefVector,
00062         TrackingParticleRefVector
00063     ) const;
00064 
00065     std::string vertexString(
00066         HepMC::GenVertex::particles_in_const_iterator,
00067         HepMC::GenVertex::particles_in_const_iterator,
00068         HepMC::GenVertex::particles_out_const_iterator,
00069         HepMC::GenVertex::particles_out_const_iterator
00070     ) const;
00071 };
00072 
00073 
00074 TrackHistoryAnalyzer::TrackHistoryAnalyzer(const edm::ParameterSet& config) : classifier_(config)
00075 {
00076     trackProducer_ = config.getUntrackedParameter<edm::InputTag> ( "trackProducer" );
00077 }
00078 
00079 
00080 TrackHistoryAnalyzer::~TrackHistoryAnalyzer() { }
00081 
00082 
00083 void TrackHistoryAnalyzer::analyze(const edm::Event& event, const edm::EventSetup& setup)
00084 {
00085     // Track collection
00086     edm::Handle<edm::View<reco::Track> > trackCollection;
00087     event.getByLabel(trackProducer_, trackCollection);
00088 
00089     // Set the classifier for a new event
00090     classifier_.newEvent(event, setup);
00091 
00092     // Get a constant reference to the track history associated to the classifier
00093     TrackHistory const & tracer = classifier_.history();
00094 
00095     // Loop over the track collection.
00096     for (std::size_t index = 0; index < trackCollection->size(); index++)
00097     {
00098         std::cout << std::endl << "History for track #" << index << " : " << std::endl;
00099 
00100         // Classify the track and detect for fakes
00101         if ( ! classifier_.evaluate( reco::TrackBaseRef(trackCollection, index) ).is(TrackClassifier::Fake) )
00102         {
00103             // Get the list of TrackingParticles associated to
00104             TrackHistory::SimParticleTrail simParticles(tracer.simParticleTrail());
00105 
00106             // Loop over all simParticles
00107             for (std::size_t hindex=0; hindex<simParticles.size(); hindex++)
00108             {
00109                 std::cout << "  simParticles [" << hindex << "] : "
00110                           << particleString(simParticles[hindex]->pdgId())
00111                           << std::endl;
00112             }
00113 
00114             // Get the list of TrackingVertexes associated to
00115             TrackHistory::SimVertexTrail simVertexes(tracer.simVertexTrail());
00116 
00117             // Loop over all simVertexes
00118             if ( !simVertexes.empty() )
00119             {
00120                 for (std::size_t hindex=0; hindex<simVertexes.size(); hindex++)
00121                 {
00122                     std::cout << "  simVertex    [" << hindex << "] : "
00123                               << vertexString(
00124                                   simVertexes[hindex]->sourceTracks(),
00125                                   simVertexes[hindex]->daughterTracks()
00126                               )
00127                               << std::endl;
00128                 }
00129             }
00130             else
00131                 std::cout << "  simVertex no found" << std::endl;
00132 
00133             // Get the list of GenParticles associated to
00134             TrackHistory::GenParticleTrail genParticles(tracer.genParticleTrail());
00135 
00136             // Loop over all genParticles
00137             for (std::size_t hindex=0; hindex<genParticles.size(); hindex++)
00138             {
00139                 std::cout << "  genParticles [" << hindex << "] : "
00140                           << particleString(genParticles[hindex]->pdg_id())
00141                           << std::endl;
00142             }
00143 
00144             // Get the list of TrackingVertexes associated to
00145             TrackHistory::GenVertexTrail genVertexes(tracer.genVertexTrail());
00146 
00147             // Loop over all simVertexes
00148             if ( !genVertexes.empty() )
00149             {
00150                 for (std::size_t hindex=0; hindex<genVertexes.size(); hindex++)
00151                 {
00152                     std::cout << "  genVertex    [" << hindex << "] : "
00153                               << vertexString(
00154                                   genVertexes[hindex]->particles_in_const_begin(),
00155                                   genVertexes[hindex]->particles_in_const_end(),
00156                                   genVertexes[hindex]->particles_out_const_begin(),
00157                                   genVertexes[hindex]->particles_out_const_end()
00158                               )
00159                               << std::endl;
00160                 }
00161             }
00162             else
00163                 std::cout << "  genVertex no found" << std::endl;
00164         }
00165         else
00166             std::cout << "  fake track" << std::endl;
00167 
00168         std::cout << "  track categories : " << classifier_;
00169         std::cout << std::endl;
00170     }
00171 }
00172 
00173 
00174 void
00175 TrackHistoryAnalyzer::beginRun(const edm::Run& run, const edm::EventSetup& setup)
00176 {
00177     // Get the particles table.
00178     setup.getData( pdt_ );
00179 }
00180 
00181 void
00182 TrackHistoryAnalyzer::beginJob()
00183 {
00184     totalTracks_ = 0;
00185 }
00186 
00187 
00188 std::string TrackHistoryAnalyzer::particleString(int pdgId) const
00189 {
00190     ParticleData const * pid;
00191 
00192     std::ostringstream vDescription;
00193 
00194     HepPDT::ParticleID particleType(pdgId);
00195 
00196     if (particleType.isValid())
00197     {
00198         pid = pdt_->particle(particleType);
00199         if (pid)
00200             vDescription << pid->name();
00201         else
00202             vDescription << pdgId;
00203     }
00204     else
00205         vDescription << pdgId;
00206 
00207     return vDescription.str();
00208 }
00209 
00210 
00211 std::string TrackHistoryAnalyzer::vertexString(
00212     TrackingParticleRefVector in,
00213     TrackingParticleRefVector out
00214 ) const
00215 {
00216     ParticleData const * pid;
00217 
00218     std::ostringstream vDescription;
00219 
00220     for (std::size_t j = 0; j < in.size(); j++)
00221     {
00222         if (!j) vDescription << "(";
00223 
00224         HepPDT::ParticleID particleType(in[j]->pdgId());
00225 
00226         if (particleType.isValid())
00227         {
00228             pid = pdt_->particle(particleType);
00229             if (pid)
00230                 vDescription << pid->name();
00231             else
00232                 vDescription << in[j]->pdgId();
00233         }
00234         else
00235             vDescription << in[j]->pdgId();
00236 
00237         if (j == in.size() - 1) vDescription << ")";
00238         else vDescription << ",";
00239     }
00240 
00241     vDescription << "->";
00242 
00243     for (std::size_t j = 0; j < out.size(); j++)
00244     {
00245         if (!j) vDescription << "(";
00246 
00247         HepPDT::ParticleID particleType(out[j]->pdgId());
00248 
00249         if (particleType.isValid())
00250         {
00251             pid = pdt_->particle(particleType);
00252             if (pid)
00253                 vDescription << pid->name();
00254             else
00255                 vDescription << out[j]->pdgId();
00256         }
00257         else
00258             vDescription << out[j]->pdgId();
00259 
00260         if (j == out.size() - 1) vDescription << ")";
00261         else vDescription << ",";
00262     }
00263 
00264     return vDescription.str();
00265 }
00266 
00267 
00268 std::string TrackHistoryAnalyzer::vertexString(
00269     HepMC::GenVertex::particles_in_const_iterator in_begin,
00270     HepMC::GenVertex::particles_in_const_iterator in_end,
00271     HepMC::GenVertex::particles_out_const_iterator out_begin,
00272     HepMC::GenVertex::particles_out_const_iterator out_end
00273 ) const
00274 {
00275     ParticleData const * pid;
00276 
00277     std::ostringstream vDescription;
00278 
00279     std::size_t j = 0;
00280 
00281     HepMC::GenVertex::particles_in_const_iterator in, itmp;
00282 
00283     for (in = in_begin; in != in_end; in++, j++)
00284     {
00285         if (!j) vDescription << "(";
00286 
00287         HepPDT::ParticleID particleType((*in)->pdg_id());
00288 
00289         if (particleType.isValid())
00290         {
00291             pid = pdt_->particle(particleType);
00292             if (pid)
00293                 vDescription << pid->name();
00294             else
00295                 vDescription << (*in)->pdg_id();
00296         }
00297         else
00298             vDescription << (*in)->pdg_id();
00299 
00300         itmp = in;
00301 
00302         if (++itmp == in_end) vDescription << ")";
00303         else vDescription << ",";
00304     }
00305 
00306     vDescription << "->";
00307     j = 0;
00308 
00309     HepMC::GenVertex::particles_out_const_iterator out, otmp;
00310 
00311     for (out = out_begin; out != out_end; out++, j++)
00312     {
00313         if (!j) vDescription << "(";
00314 
00315         HepPDT::ParticleID particleType((*out)->pdg_id());
00316 
00317         if (particleType.isValid())
00318         {
00319             pid = pdt_->particle(particleType);
00320             if (pid)
00321                 vDescription << pid->name();
00322             else
00323                 vDescription << (*out)->pdg_id();
00324         }
00325         else
00326             vDescription << (*out)->pdg_id();
00327 
00328         otmp = out;
00329 
00330         if (++otmp == out_end) vDescription << ")";
00331         else vDescription << ",";
00332     }
00333 
00334     return vDescription.str();
00335 }
00336 
00337 
00338 DEFINE_FWK_MODULE(TrackHistoryAnalyzer);