CMS 3D CMS Logo

All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 
18 private:
19  std::string getParticleName(int id) const;
20  void analyze(const edm::Event &, const edm::EventSetup &) override;
22  void printDecay(const reco::Candidate &, const std::string &pre) const;
27  typedef std::vector<int> vint;
28  vint status_;
30  void printInfo(const reco::Candidate &) const;
32  bool accept(const reco::Candidate &) const;
34  bool hasValidDaughters(const reco::Candidate &) const;
36  std::vector<const reco::Candidate *> cands_;
37 };
38 
44 #include <iostream>
45 #include <algorithm>
46 using namespace std;
47 using namespace edm;
48 using namespace reco;
49 
51  : srcToken_(consumes<edm::View<reco::Candidate> >(cfg.getParameter<InputTag>("src"))),
52  printP4_(cfg.getUntrackedParameter<bool>("printP4", false)),
53  printPtEtaPhi_(cfg.getUntrackedParameter<bool>("printPtEtaPhi", false)),
54  printVertex_(cfg.getUntrackedParameter<bool>("printVertex", false)),
55  printStatus_(cfg.getUntrackedParameter<bool>("printStatus", false)),
56  printIndex_(cfg.getUntrackedParameter<bool>("printIndex", false)),
57  status_(cfg.getUntrackedParameter<vint>("status", vint())) {}
58 
60  if (status_.empty())
61  return true;
62  return find(status_.begin(), status_.end(), c.status()) != status_.end();
63 }
64 
66  size_t ndau = c.numberOfDaughters();
67  for (size_t i = 0; i < ndau; ++i)
68  if (accept(*c.daughter(i)))
69  return true;
70  return false;
71 }
72 
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 
84  es.getData(pdt_);
86  event.getByToken(srcToken_, particles);
87  cands_.clear();
88  for (View<Candidate>::const_iterator p = particles->begin(); p != particles->end(); ++p) {
89  cands_.push_back(&*p);
90  }
91  for (View<Candidate>::const_iterator p = particles->begin(); p != particles->end(); ++p) {
92  if (accept(*p)) {
93  if (p->mother() == nullptr) {
94  cout << "-- decay tree: --" << endl;
95  printDecay(*p, "");
96  }
97  }
98  }
99 }
100 
102  if (printP4_)
103  cout << " (" << c.px() << ", " << c.py() << ", " << c.pz() << "; " << c.energy() << ")";
104  if (printPtEtaPhi_)
105  cout << " [" << c.pt() << ": " << c.eta() << ", " << c.phi() << "]";
106  if (printVertex_)
107  cout << " {" << c.vx() << ", " << c.vy() << ", " << c.vz() << "}";
108  if (printStatus_)
109  cout << "{status: " << c.status() << "}";
110  if (printIndex_) {
111  int idx = -1;
112  vector<const Candidate *>::const_iterator found = find(cands_.begin(), cands_.end(), &c);
113  if (found != cands_.end()) {
114  idx = found - cands_.begin();
115  cout << " <idx: " << idx << ">";
116  }
117  }
118 }
119 
120 void ParticleTreeDrawer::printDecay(const Candidate &c, const string &pre) const {
121  cout << getParticleName(c.pdgId());
122  printInfo(c);
123  cout << endl;
124 
125  size_t ndau = c.numberOfDaughters(), validDau = 0;
126  for (size_t i = 0; i < ndau; ++i)
127  if (accept(*c.daughter(i)))
128  ++validDau;
129  if (validDau == 0)
130  return;
131 
132  bool lastLevel = true;
133  for (size_t i = 0; i < ndau; ++i) {
134  if (hasValidDaughters(*c.daughter(i))) {
135  lastLevel = false;
136  break;
137  }
138  }
139 
140  if (lastLevel) {
141  cout << pre << "+-> ";
142  size_t vd = 0;
143  for (size_t i = 0; i < ndau; ++i) {
144  const Candidate *d = c.daughter(i);
145  if (accept(*d)) {
146  cout << getParticleName(d->pdgId());
147  printInfo(*d);
148  if (vd != validDau - 1)
149  cout << " ";
150  vd++;
151  }
152  }
153  cout << endl;
154  return;
155  }
156 
157  for (size_t i = 0; i < ndau; ++i) {
158  const Candidate *d = c.daughter(i);
159  assert(d != nullptr);
160  if (accept(*d)) {
161  cout << pre << "+-> ";
162  string prepre(pre);
163  if (i == ndau - 1)
164  prepre += " ";
165  else
166  prepre += "| ";
167  printDecay(*d, prepre);
168  }
169  }
170 }
171 
173 
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) ...
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:19
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
bool getData(T &iHolder) const
Definition: EventSetup.h:113
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
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
d
Definition: ztail.py:151
std::vector< DeviationSensor2D * > vd
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