CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_4/src/PhysicsTools/HepMCCandAlgos/plugins/ParticleDecayDrawer.cc

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