CMS 3D CMS Logo

Public Member Functions | Private Member Functions | Private Attributes

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.

                                                                   :
  src_( cfg.getParameter<InputTag>( "src" ) ),
  printP4_( cfg.getUntrackedParameter<bool>( "printP4", false ) ),
  printPtEtaPhi_( cfg.getUntrackedParameter<bool>( "printPtEtaPhi", false ) ),
  printVertex_( cfg.getUntrackedParameter<bool>( "printVertex", false ) ) {
}

Member Function Documentation

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

accept candidate

Referenced by decay().

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 gather_cfg::cout, decay(), edm::EventSetup::getData(), j, m, reco::Candidate::mother(), n, AlCaHLTBitMon_ParallelJobs::p, pdt_, select(), createPayload::skip, and src_.

                                                                              {  
  es.getData( pdt_ );
  Handle<View<Candidate> > particles;
  event.getByLabel( src_, particles );
  list<const Candidate *> skip;
  vector<const Candidate *> nodes, moms;
  for( View<Candidate>::const_iterator p = particles->begin();
       p != particles->end(); ++ p ) {
    if( p->numberOfMothers() > 1 ) {
      if ( select( * p ) ) {
        skip.push_back( & * p );
        nodes.push_back( & * p );
        for( size_t j = 0; j < p->numberOfMothers(); ++ j ) {
          const Candidate * mom = p->mother( j );
          const Candidate * grandMom;
          while ( ( grandMom = mom->mother() ) != 0 )
            mom = grandMom;
          if ( select( * mom ) ) {
            moms.push_back( mom );
          }
        }
      }
    }
  }
  cout << "-- decay: --" << endl;
  if( moms.size() > 0 ) {
    if ( moms.size() > 1 )
      for( size_t m = 0; m < moms.size(); ++ m ) {
        string dec = decay( * moms[ m ], skip );
        if ( ! dec.empty() )
          cout << "{ " << dec << " } ";
      }
    else 
      cout << decay( * moms[ 0 ], skip );
  }
  if ( nodes.size() > 0 ) {
    cout << "-> ";
    if ( nodes.size() > 1 ) {
      for( size_t n = 0; n < nodes.size(); ++ n ) {    
        skip.remove( nodes[ n ] );
        string dec = decay( * nodes[ n ], skip );
        if ( ! dec.empty() ) {
          if ( dec.find( "->", 0 ) != string::npos )
            cout << " ( " << dec << " )";
          else 
            cout << " " << dec;
        }
      }
    } else {
      skip.remove( nodes[ 0 ] );
      cout << decay( * nodes[ 0 ], skip );
    }
  }
  cout << endl;
}
string ParticleDecayDrawer::decay ( const reco::Candidate ,
std::list< const reco::Candidate * > &   
) const [private]

Definition at line 132 of file ParticleDecayDrawer.cc.

References accept(), trackerHits::c, reco::Candidate::daughter(), spr::find(), i, reco::Candidate::numberOfDaughters(), dbtoconf::out, reco::Candidate::pdgId(), pdt_, printP4(), and createPayload::skip.

Referenced by analyze().

                                                                          {
  string out;
  if ( find( skip.begin(), skip.end(), & c ) != skip.end() ) 
    return out;
  skip.push_back( & c );

  
  int id = c.pdgId();
  const ParticleData * pd = pdt_->particle( id );  
  assert( pd != 0 );
  out += ( pd->name() + printP4( c ) );
  
  size_t validDau = 0, ndau = c.numberOfDaughters();
  for( size_t i = 0; i < ndau; ++ i )
    if ( accept( * c.daughter( i ), skip ) )
      ++ validDau;
  if ( validDau == 0 ) return out;

  out += " ->";
  
  for( size_t i = 0; i < ndau; ++ i ) {
    const Candidate * d = c.daughter( i );
    if ( accept( * d, skip ) ) {
      string dec = decay( * d, skip );
      if ( dec.find( "->", 0 ) != string::npos )
        out += ( " ( " + dec + " )" );
      else
        out += ( " " + dec );
    }
  }
  return out;
}
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().

                                                                           {
  size_t ndau = c.numberOfDaughters();
  for( size_t i = 0; i < ndau; ++ i )
    if ( select( * c.daughter( i ) ) )
      return true;
  return false;
}
string ParticleDecayDrawer::printP4 ( const reco::Candidate c) const [private]

print 4 momenta

Definition at line 124 of file ParticleDecayDrawer.cc.

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

Referenced by decay().

                                                               {
  ostringstream cout;
  if ( printP4_ ) cout << " (" << c.px() << ", " << c.py() << ", " << c.pz() << "; " << c.energy() << ")"; 
  if ( printPtEtaPhi_ ) cout << " [" << c.pt() << ": " << c.eta() << ", " << c.phi() << "]";
  if ( printVertex_ ) cout << " {" << c.vx() << ", " << c.vy() << ", " << c.vz() << "}";
  return cout.str();
}
bool ParticleDecayDrawer::select ( const reco::Candidate c) const [private]

select candidate

Definition at line 56 of file ParticleDecayDrawer.cc.

References reco::Candidate::status().

Referenced by EcalZeroSuppressor< EBDataFrame >::accept(), analyze(), and hasValidDaughters().

                                                                {
  return c.status() == 3;
}

Member Data Documentation

Definition at line 18 of file ParticleDecayDrawer.cc.

Referenced by analyze(), and decay().

print parameters

Definition at line 20 of file ParticleDecayDrawer.cc.

Referenced by printP4().

Definition at line 20 of file ParticleDecayDrawer.cc.

Referenced by printP4().

Definition at line 20 of file ParticleDecayDrawer.cc.

Referenced by printP4().

Definition at line 16 of file ParticleDecayDrawer.cc.

Referenced by analyze().