CMS 3D CMS Logo

ParticleTreeDrawer.cc
Go to the documentation of this file.
1 /* class ParticleTreeDrawer
2  *
3  * \author Luca Lista, INFN
4  */
5 #include <sstream>
6 
13 
15 public:
17 private:
18  std::string getParticleName( int id ) const;
19  void analyze( const edm::Event &, const edm::EventSetup&) override;
21  void printDecay( const reco::Candidate &, const std::string & pre ) const;
26  typedef std::vector<int> vint;
27  vint status_;
29  void printInfo( const reco::Candidate & ) const;
31  bool accept( const reco::Candidate & ) const;
33  bool hasValidDaughters( const reco::Candidate & ) const;
35  std::vector<const reco::Candidate *> cands_;
36 };
37 
43 #include <iostream>
44 #include <algorithm>
45 using namespace std;
46 using namespace edm;
47 using namespace reco;
48 
50  srcToken_( consumes<edm::View<reco::Candidate> >(cfg.getParameter<InputTag>( "src" ) ) ),
51  printP4_( cfg.getUntrackedParameter<bool>( "printP4", false ) ),
52  printPtEtaPhi_( cfg.getUntrackedParameter<bool>( "printPtEtaPhi", false ) ),
53  printVertex_( cfg.getUntrackedParameter<bool>( "printVertex", false ) ),
54  printStatus_( cfg.getUntrackedParameter<bool>( "printStatus", false ) ),
55  printIndex_( cfg.getUntrackedParameter<bool>( "printIndex", false ) ),
56  status_( cfg.getUntrackedParameter<vint>( "status", vint() ) ) {
57 }
58 
60  if ( status_.empty() ) return true;
61  return find( status_.begin(), status_.end(), c.status() ) != status_.end();
62 }
63 
65  size_t ndau = c.numberOfDaughters();
66  for( size_t i = 0; i < ndau; ++ i )
67  if ( accept( * c.daughter( i ) ) )
68  return true;
69  return false;
70 }
71 
73 {
74  const ParticleData * pd = pdt_->particle( id );
75  if (!pd) {
76  std::ostringstream ss;
77  ss << "P" << id;
78  return ss.str();
79  } else
80  return pd->name();
81 }
82 
83 void ParticleTreeDrawer::analyze( const Event & event, const EventSetup & es ) {
84  es.getData( pdt_ );
86  event.getByToken( srcToken_, particles );
87  cands_.clear();
88  for( View<Candidate>::const_iterator p = particles->begin();
89  p != particles->end(); ++ p ) {
90  cands_.push_back( & * p );
91  }
92  for( View<Candidate>::const_iterator p = particles->begin();
93  p != particles->end(); ++ p ) {
94  if ( accept( * p ) ) {
95  if ( p->mother() == nullptr ) {
96  cout << "-- decay tree: --" << endl;
97  printDecay( * p, "" );
98  }
99  }
100  }
101 }
102 
104  if ( printP4_ ) cout << " (" << c.px() << ", " << c.py() << ", " << c.pz() << "; " << c.energy() << ")";
105  if ( printPtEtaPhi_ ) cout << " [" << c.pt() << ": " << c.eta() << ", " << c.phi() << "]";
106  if ( printVertex_ ) cout << " {" << c.vx() << ", " << c.vy() << ", " << c.vz() << "}";
107  if ( printStatus_ ) cout << "{status: " << c.status() << "}";
108  if ( printIndex_ ) {
109  int idx = -1;
110  vector<const Candidate *>::const_iterator found = find( cands_.begin(), cands_.end(), & c );
111  if ( found != cands_.end() ) {
112  idx = found - cands_.begin();
113  cout << " <idx: " << idx << ">";
114  }
115  }
116 }
117 
118 void ParticleTreeDrawer::printDecay( const Candidate & c, const string & pre ) const {
119  cout << getParticleName( c.pdgId() );
120  printInfo( c );
121  cout << endl;
122 
123  size_t ndau = c.numberOfDaughters(), validDau = 0;
124  for( size_t i = 0; i < ndau; ++ i )
125  if ( accept( * c.daughter( i ) ) )
126  ++ validDau;
127  if ( validDau == 0 ) return;
128 
129  bool lastLevel = true;
130  for( size_t i = 0; i < ndau; ++ i ) {
131  if ( hasValidDaughters( * c.daughter( i ) ) ) {
132  lastLevel = false;
133  break;
134  }
135  }
136 
137  if ( lastLevel ) {
138  cout << pre << "+-> ";
139  size_t vd = 0;
140  for( size_t i = 0; i < ndau; ++ i ) {
141  const Candidate * d = c.daughter( i );
142  if ( accept( * d ) ) {
143  cout << getParticleName( d->pdgId() );
144  printInfo( * d );
145  if ( vd != validDau - 1 )
146  cout << " ";
147  vd ++;
148  }
149  }
150  cout << endl;
151  return;
152  }
153 
154  for( size_t i = 0; i < ndau; ++i ) {
155  const Candidate * d = c.daughter( i );
156  assert( d != nullptr );
157  if ( accept( * d ) ) {
158  cout << pre << "+-> ";
159  string prepre( pre );
160  if ( i == ndau - 1 ) prepre += " ";
161  else prepre += "| ";
162  printDecay( * d, prepre );
163  }
164  }
165 }
166 
168 
170 
171 
virtual double pz() const =0
z coordinate of momentum vector
virtual double vx() const =0
x coordinate of vertex position
void printInfo(const reco::Candidate &) const
print 4 momenta
std::vector< const reco::Candidate * > cands_
pointer to collection
edm::ESHandle< ParticleDataTable > pdt_
std::string getParticleName(int id) const
edm::EDGetTokenT< edm::View< reco::Candidate > > srcToken_
virtual const Candidate * daughter(size_type i) const =0
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
#define nullptr
virtual double vy() const =0
y coordinate of vertex position
std::vector< int > vint
accepted status codes
void printDecay(const reco::Candidate &, const std::string &pre) const
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
bool getData(T &iHolder) const
Definition: EventSetup.h:111
void analyze(const edm::Event &, const edm::EventSetup &) override
virtual int status() const =0
status word
virtual double energy() const =0
energy
virtual double py() const =0
y coordinate of momentum vector
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::vector< DeviationSensor2D * > vd
bool accept(const reco::Candidate &) const
accept candidate
virtual int pdgId() const =0
PDG identifier.
ParticleTreeDrawer(const edm::ParameterSet &)
bool hasValidDaughters(const reco::Candidate &) const
has valid daughters in the chain
bool printP4_
print parameters
HepPDT::ParticleData ParticleData
virtual double eta() const =0
momentum pseudorapidity
virtual double pt() const =0
transverse momentum
fixed size matrix
HLT enums.
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
virtual double vz() const =0
z coordinate of vertex position
virtual double px() const =0
x coordinate of momentum vector
virtual size_type numberOfDaughters() const =0
number of daughters
virtual double phi() const =0
momentum azimuthal angle
Definition: event.py:1