CMS 3D CMS Logo

ParticleDecayDrawer Class Reference

Inheritance diagram for ParticleDecayDrawer:

edm::EDAnalyzer

List of all members.

Public Member Functions

 ParticleDecayDrawer (const edm::ParameterSet &)

Private Member Functions

bool accept (const reco::Candidate &, const std::list< const reco::Candidate * > &) const
 accept candidate
void analyze (const edm::Event &, const edm::EventSetup &)
std::string decay (const reco::Candidate &, std::list< const reco::Candidate * > &) const
bool hasValidDaughters (const reco::Candidate &) const
 has valid daughters in the chain
std::string printP4 (const reco::Candidate &) const
 print 4 momenta
bool select (const reco::Candidate &) const
 select candidate

Private Attributes

edm::ESHandle< ParticleDataTablepdt_
bool printP4_
 print parameters
bool printPtEtaPhi_
bool printVertex_
edm::InputTag src_


Detailed Description

Definition at line 11 of file ParticleDecayDrawer.cc.


Constructor & Destructor Documentation

ParticleDecayDrawer::ParticleDecayDrawer ( const edm::ParameterSet cfg  ) 

Definition at line 44 of file ParticleDecayDrawer.cc.

00044                                                                    :
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 }


Member Function Documentation

bool ParticleDecayDrawer::accept ( const reco::Candidate ,
const std::list< const reco::Candidate * > &   
) const [private]

accept candidate

void ParticleDecayDrawer::analyze ( const edm::Event event,
const edm::EventSetup es 
) [private, virtual]

Implements edm::EDAnalyzer.

Definition at line 68 of file ParticleDecayDrawer.cc.

References GenMuonPlsPt100GeV_cfg::cout, decay(), lat::endl(), edm::EventSetup::getData(), j, m, reco::Candidate::mother(), n, p, pdt_, select(), and src_.

00068                                                                               {  
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     } else {
00116       skip.remove( nodes[ 0 ] );
00117       cout << decay( * nodes[ 0 ], skip );
00118     }
00119   }
00120   cout << endl;
00121 }

std::string ParticleDecayDrawer::decay ( const reco::Candidate ,
std::list< const reco::Candidate * > &   
) const [private]

Referenced by analyze().

bool ParticleDecayDrawer::hasValidDaughters ( const reco::Candidate c  )  const [private]

has valid daughters in the chain

Definition at line 60 of file ParticleDecayDrawer.cc.

References reco::Candidate::daughter(), i, reco::Candidate::numberOfDaughters(), and select().

00060                                                                            {
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 }

string ParticleDecayDrawer::printP4 ( const reco::Candidate c  )  const [private]

print 4 momenta

Definition at line 123 of file ParticleDecayDrawer.cc.

References GenMuonPlsPt100GeV_cfg::cout, reco::Particle::energy(), reco::Particle::eta(), reco::Particle::phi(), printP4_, printPtEtaPhi_, printVertex_, reco::Particle::pt(), reco::Particle::px(), reco::Particle::py(), reco::Particle::pz(), reco::Particle::vx(), reco::Particle::vy(), and reco::Particle::vz().

00123                                                                {
00124   ostringstream cout;
00125   if ( printP4_ ) cout << " (" << c.px() << ", " << c.py() << ", " << c.pz() << "; " << c.energy() << ")"; 
00126   if ( printPtEtaPhi_ ) cout << " [" << c.pt() << ": " << c.eta() << ", " << c.phi() << "]";
00127   if ( printVertex_ ) cout << " {" << c.vx() << ", " << c.vy() << ", " << c.vz() << "}";
00128   return cout.str();
00129 }

bool ParticleDecayDrawer::select ( const reco::Candidate c  )  const [private]

select candidate

Definition at line 56 of file ParticleDecayDrawer.cc.

References reco::Particle::status().

Referenced by analyze(), and hasValidDaughters().

00056                                                                 {
00057   return c.status() == 3;
00058 }


Member Data Documentation

edm::ESHandle<ParticleDataTable> ParticleDecayDrawer::pdt_ [private]

Definition at line 18 of file ParticleDecayDrawer.cc.

Referenced by analyze().

bool ParticleDecayDrawer::printP4_ [private]

print parameters

Definition at line 20 of file ParticleDecayDrawer.cc.

Referenced by printP4().

bool ParticleDecayDrawer::printPtEtaPhi_ [private]

Definition at line 20 of file ParticleDecayDrawer.cc.

Referenced by printP4().

bool ParticleDecayDrawer::printVertex_ [private]

Definition at line 20 of file ParticleDecayDrawer.cc.

Referenced by printP4().

edm::InputTag ParticleDecayDrawer::src_ [private]

Definition at line 16 of file ParticleDecayDrawer.cc.

Referenced by analyze().


The documentation for this class was generated from the following file:
Generated on Tue Jun 9 18:29:27 2009 for CMSSW by  doxygen 1.5.4