CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10/src/PhysicsTools/HepMCCandAlgos/plugins/ParticleListDrawer.cc

Go to the documentation of this file.
00001 #include <memory>
00002 #include <string>
00003 #include <iostream>
00004 #include <sstream>
00005 
00006 #include "FWCore/Framework/interface/Frameworkfwd.h"
00007 #include "FWCore/Framework/interface/EDAnalyzer.h"
00008 #include "FWCore/Framework/interface/ESHandle.h"
00009 #include "FWCore/Framework/interface/Event.h"
00010 #include "FWCore/Framework/interface/EventSetup.h"
00011 #include "FWCore/Framework/interface/MakerMacros.h"
00012 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00013 #include "FWCore/Utilities/interface/InputTag.h"
00014 #include "SimGeneral/HepPDTRecord/interface/ParticleDataTable.h"
00015 #include "DataFormats/Candidate/interface/Candidate.h"
00016 #include "DataFormats/Candidate/interface/CandidateFwd.h"
00017 #include "DataFormats/Common/interface/Ref.h"
00018 
00037 using namespace std;
00038 using namespace reco;
00039 using namespace edm;
00040 
00041 class ParticleListDrawer : public edm::EDAnalyzer {
00042   public:
00043     explicit ParticleListDrawer(const edm::ParameterSet & );
00044     ~ParticleListDrawer() {};
00045     void analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup);
00046 
00047   private:
00048     std::string getParticleName( int id ) const;
00049 
00050     edm::InputTag src_;
00051     edm::ESHandle<ParticleDataTable> pdt_;
00052     int maxEventsToPrint_; // Must be signed, because -1 is used for no limit
00053     unsigned int nEventAnalyzed_;
00054     bool printOnlyHardInteraction_;
00055     bool printVertex_;
00056     bool useMessageLogger_;
00057 };
00058 
00059 ParticleListDrawer::ParticleListDrawer(const edm::ParameterSet & pset) :
00060   src_(pset.getParameter<InputTag>("src")),
00061   maxEventsToPrint_ (pset.getUntrackedParameter<int>("maxEventsToPrint",1)),
00062   nEventAnalyzed_(0),
00063   printOnlyHardInteraction_(pset.getUntrackedParameter<bool>("printOnlyHardInteraction", false)),
00064   printVertex_(pset.getUntrackedParameter<bool>("printVertex", false)),
00065   useMessageLogger_(pset.getUntrackedParameter<bool>("useMessageLogger", false)) {
00066 }
00067 
00068 std::string ParticleListDrawer::getParticleName(int id) const
00069 {
00070   const ParticleData * pd = pdt_->particle( id );
00071   if (!pd) {
00072     std::ostringstream ss;
00073     ss << "P" << id;
00074     return ss.str();
00075   } else
00076     return pd->name();
00077 }
00078 
00079 void ParticleListDrawer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {  
00080   Handle<reco::CandidateView> particles;
00081   iEvent.getByLabel (src_, particles );
00082   iSetup.getData( pdt_ );
00083 
00084   if(maxEventsToPrint_ < 0 || nEventAnalyzed_ < static_cast<unsigned int>(maxEventsToPrint_)) {
00085     ostringstream out;
00086     char buf[256];
00087 
00088     out << endl
00089         << "[ParticleListDrawer] analysing particle collection " << src_.label() << endl;
00090 
00091     snprintf(buf, 256, " idx  |    ID -       Name |Stat|  Mo1  Mo2  Da1  Da2 |nMo nDa|    pt       eta     phi   |     px         py         pz        m     |"); 
00092     out << buf;
00093     if (printVertex_) {
00094       snprintf(buf, 256, "        vx       vy        vz     |");
00095       out << buf;
00096     }
00097     out << endl; 
00098 
00099     int idx  = -1;
00100     int iMo1 = -1;
00101     int iMo2 = -1;
00102     int iDa1 = -1;
00103     int iDa2 = -1;
00104     vector<const reco::Candidate *> cands;
00105     vector<const Candidate *>::const_iterator found = cands.begin();
00106     for(CandidateView::const_iterator p = particles->begin();
00107         p != particles->end(); ++ p) {
00108       cands.push_back(&*p);
00109     }
00110 
00111     for(CandidateView::const_iterator p  = particles->begin();
00112         p != particles->end(); 
00113         p ++) {
00114       if (printOnlyHardInteraction_ && p->status() != 3) continue;
00115 
00116       // Particle Name
00117       int id = p->pdgId();
00118       string particleName = getParticleName(id);
00119       
00120       // Particle Index
00121       idx =  p - particles->begin();
00122 
00123       // Particles Mothers and Daighters
00124       iMo1 = -1;
00125       iMo2 = -1;
00126       iDa1 = -1;
00127       iDa2 = -1;
00128       int nMo = p->numberOfMothers();
00129       int nDa = p->numberOfDaughters();
00130 
00131       found = find(cands.begin(), cands.end(), p->mother(0));
00132       if(found != cands.end()) iMo1 = found - cands.begin() ;
00133 
00134       found = find(cands.begin(), cands.end(), p->mother(nMo-1));
00135       if(found != cands.end()) iMo2 = found - cands.begin() ;
00136      
00137       found = find(cands.begin(), cands.end(), p->daughter(0));
00138       if(found != cands.end()) iDa1 = found - cands.begin() ;
00139 
00140       found = find(cands.begin(), cands.end(), p->daughter(nDa-1));
00141       if(found != cands.end()) iDa2 = found - cands.begin() ;
00142 
00143       char buf[256];
00144       snprintf(buf, 256,
00145              " %4d | %5d - %10s | %2d | %4d %4d %4d %4d | %2d %2d | %7.3f %10.3f %6.3f | %10.3f %10.3f %10.3f %8.3f |",
00146              idx,
00147              p->pdgId(),
00148              particleName.c_str(),
00149              p->status(),
00150              iMo1,iMo2,iDa1,iDa2,nMo,nDa,
00151              p->pt(),
00152              p->eta(),
00153              p->phi(),
00154              p->px(),
00155              p->py(),
00156              p->pz(),
00157              p->mass()
00158             );
00159       out << buf;
00160 
00161       if (printVertex_) {
00162         snprintf(buf, 256, " %10.3f %10.3f %10.3f |",
00163                  p->vertex().x(),
00164                  p->vertex().y(),
00165                  p->vertex().z());
00166         out << buf;
00167       }
00168 
00169       out << endl;
00170     }
00171     nEventAnalyzed_++;
00172 
00173     if (useMessageLogger_)
00174       LogVerbatim("ParticleListDrawer") << out.str();
00175     else
00176       cout << out.str();
00177   }
00178 }
00179 
00180 DEFINE_FWK_MODULE(ParticleListDrawer);
00181