CMS 3D CMS Logo

CMSSW_4_4_3_patch1/src/PhysicsTools/HepMCCandAlgos/plugins/ParticleTreeDrawer.cc

Go to the documentation of this file.
00001 /* class ParticleTreeDrawer
00002  *
00003  * \author Luca Lista, INFN
00004  */
00005 #include <sstream>
00006 
00007 #include "FWCore/Framework/interface/EDAnalyzer.h"
00008 #include "FWCore/Utilities/interface/InputTag.h"
00009 #include "FWCore/Framework/interface/ESHandle.h"
00010 #include "SimGeneral/HepPDTRecord/interface/ParticleDataTable.h"
00011 #include "DataFormats/Candidate/interface/Candidate.h"
00012 
00013 class ParticleTreeDrawer : public edm::EDAnalyzer {
00014 public:
00015   ParticleTreeDrawer( const edm::ParameterSet & );
00016 private:
00017   std::string getParticleName( int id ) const;
00018   void analyze( const edm::Event &, const edm::EventSetup & );
00019   edm::InputTag src_;
00020   void printDecay( const reco::Candidate &, const std::string & pre ) const;
00021   edm::ESHandle<ParticleDataTable> pdt_;
00023   bool printP4_, printPtEtaPhi_, printVertex_, printStatus_, printIndex_;
00025   typedef std::vector<int> vint;
00026   vint status_;
00028   void printInfo( const reco::Candidate & ) const;
00030   bool accept( const reco::Candidate & ) const;
00032   bool hasValidDaughters( const reco::Candidate & ) const;
00034   std::vector<const reco::Candidate *> cands_;
00035 };
00036 
00037 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00038 #include "FWCore/Framework/interface/Event.h"
00039 #include "DataFormats/Common/interface/Handle.h"
00040 #include "DataFormats/Common/interface/View.h"
00041 #include "FWCore/Framework/interface/EventSetup.h"
00042 #include "FWCore/Utilities/interface/EDMException.h"
00043 #include <iostream>
00044 #include <algorithm>
00045 using namespace std;
00046 using namespace edm;
00047 using namespace reco;
00048 
00049 ParticleTreeDrawer::ParticleTreeDrawer( const ParameterSet & cfg ) :
00050   src_( cfg.getParameter<InputTag>( "src" ) ),
00051   printP4_( cfg.getUntrackedParameter<bool>( "printP4", false ) ),
00052   printPtEtaPhi_( cfg.getUntrackedParameter<bool>( "printPtEtaPhi", false ) ),
00053   printVertex_( cfg.getUntrackedParameter<bool>( "printVertex", false ) ),
00054   printStatus_( cfg.getUntrackedParameter<bool>( "printStatus", false ) ),
00055   printIndex_( cfg.getUntrackedParameter<bool>( "printIndex", false ) ),
00056   status_( cfg.getUntrackedParameter<vint>( "status", vint() ) ) {
00057 }
00058 
00059 bool ParticleTreeDrawer::accept( const reco::Candidate & c ) const {
00060   if ( status_.size() == 0 ) return true;
00061   return find( status_.begin(), status_.end(), c.status() ) != status_.end();
00062 }
00063 
00064 bool ParticleTreeDrawer::hasValidDaughters( const reco::Candidate & c ) const {
00065   size_t ndau = c.numberOfDaughters();
00066   for( size_t i = 0; i < ndau; ++ i )
00067     if ( accept( * c.daughter( i ) ) )
00068       return true;
00069   return false;
00070 }
00071 
00072 std::string ParticleTreeDrawer::getParticleName(int id) const
00073 {
00074   const ParticleData * pd = pdt_->particle( id );
00075   if (!pd) {
00076     std::ostringstream ss;
00077     ss << "P" << id;
00078     return ss.str();
00079   } else
00080     return pd->name();
00081 }
00082 
00083 void ParticleTreeDrawer::analyze( const Event & event, const EventSetup & es ) {  
00084   es.getData( pdt_ );
00085   Handle<View<Candidate> > particles;
00086   event.getByLabel( src_, particles );
00087   cands_.clear();
00088   for( View<Candidate>::const_iterator p = particles->begin();
00089        p != particles->end(); ++ p ) {
00090     cands_.push_back( & * p );
00091   }
00092   for( View<Candidate>::const_iterator p = particles->begin();
00093        p != particles->end(); ++ p ) {
00094     if ( accept( * p ) ) {
00095       if ( p->mother() == 0 ) {
00096         cout << "-- decay tree: --" << endl;
00097         printDecay( * p, "" );
00098       }
00099     }
00100   }
00101 }
00102 
00103 void ParticleTreeDrawer::printInfo( const Candidate & c ) const {
00104   if ( printP4_ ) cout << " (" << c.px() << ", " << c.py() << ", " << c.pz() << "; " << c.energy() << ")"; 
00105   if ( printPtEtaPhi_ ) cout << " [" << c.pt() << ": " << c.eta() << ", " << c.phi() << "]";
00106   if ( printVertex_ ) cout << " {" << c.vx() << ", " << c.vy() << ", " << c.vz() << "}";
00107   if ( printStatus_ ) cout << "{status: " << c.status() << "}";
00108   if ( printIndex_ ) {
00109     int idx = -1;
00110     vector<const Candidate *>::const_iterator found = find( cands_.begin(), cands_.end(), & c );
00111     if ( found != cands_.end() ) {
00112       idx = found - cands_.begin();
00113       cout << " <idx: " << idx << ">";
00114     }
00115   }
00116 }
00117 
00118 void ParticleTreeDrawer::printDecay( const Candidate & c, const string & pre ) const {
00119   cout << getParticleName( c.pdgId() );
00120   printInfo( c );
00121   cout << endl;
00122 
00123   size_t ndau = c.numberOfDaughters(), validDau = 0;
00124   for( size_t i = 0; i < ndau; ++ i )
00125     if ( accept( * c.daughter( i ) ) )
00126       ++ validDau;
00127   if ( validDau == 0 ) return;
00128   
00129   bool lastLevel = true;
00130   for( size_t i = 0; i < ndau; ++ i ) {
00131     if ( hasValidDaughters( * c.daughter( i ) ) ) {
00132       lastLevel = false;
00133       break;
00134     }      
00135   }
00136   
00137   if ( lastLevel ) {
00138     cout << pre << "+-> ";
00139     size_t vd = 0;
00140     for( size_t i = 0; i < ndau; ++ i ) {
00141       const Candidate * d = c.daughter( i );
00142       if ( accept( * d ) ) {
00143         cout << getParticleName( d->pdgId() );
00144         printInfo( * d );
00145         if ( vd != validDau - 1 )
00146           cout << " ";
00147         vd ++;
00148       }
00149     }
00150     cout << endl;
00151     return;
00152   }
00153 
00154   for( size_t i = 0; i < ndau; ++i ) {
00155     const Candidate * d = c.daughter( i );
00156     assert( d != 0 );
00157     if ( accept( * d ) ) {
00158       cout << pre << "+-> ";
00159       string prepre( pre );
00160       if ( i == ndau - 1 ) prepre += "    ";
00161       else prepre += "|   ";
00162       printDecay( * d, prepre );
00163     }
00164   }
00165 }
00166 
00167 #include "FWCore/Framework/interface/MakerMacros.h"
00168 
00169 DEFINE_FWK_MODULE( ParticleTreeDrawer );
00170