CMS 3D CMS Logo

ParticleDecayDrawer.cc
Go to the documentation of this file.
1 /* class ParticleDecayDrawer
2  *
3  * \author Luca Lista, INFN
4  */
11 
13 public:
15 
16 private:
17  void analyze(const edm::Event &, const edm::EventSetup &) override;
19  std::string decay(const reco::Candidate &, std::list<const reco::Candidate *> &) const;
24  std::string printP4(const reco::Candidate &) const;
26  bool accept(const reco::Candidate &, const std::list<const reco::Candidate *> &) const;
28  bool select(const reco::Candidate &) const;
30  bool hasValidDaughters(const reco::Candidate &) const;
31 };
32 
38 #include <iostream>
39 #include <sstream>
40 #include <algorithm>
41 using namespace std;
42 using namespace edm;
43 using namespace reco;
44 
46  : srcToken_(consumes<edm::View<reco::Candidate> >(cfg.getParameter<InputTag>("src"))),
47  printP4_(cfg.getUntrackedParameter<bool>("printP4", false)),
48  printPtEtaPhi_(cfg.getUntrackedParameter<bool>("printPtEtaPhi", false)),
49  printVertex_(cfg.getUntrackedParameter<bool>("printVertex", false)) {}
50 
51 bool ParticleDecayDrawer::accept(const reco::Candidate &c, const list<const Candidate *> &skip) const {
52  if (find(skip.begin(), skip.end(), &c) != skip.end())
53  return false;
54  return select(c);
55 }
56 
57 bool ParticleDecayDrawer::select(const reco::Candidate &c) const { return c.status() == 3; }
58 
60  size_t ndau = c.numberOfDaughters();
61  for (size_t i = 0; i < ndau; ++i)
62  if (select(*c.daughter(i)))
63  return true;
64  return false;
65 }
66 
68  es.getData(pdt_);
70  event.getByToken(srcToken_, particles);
71  list<const Candidate *> skip;
72  vector<const Candidate *> nodes, moms;
73  for (View<Candidate>::const_iterator p = particles->begin(); p != particles->end(); ++p) {
74  if (p->numberOfMothers() > 1) {
75  if (select(*p)) {
76  skip.push_back(&*p);
77  nodes.push_back(&*p);
78  for (size_t j = 0; j < p->numberOfMothers(); ++j) {
79  const Candidate *mom = p->mother(j);
80  const Candidate *grandMom;
81  while ((grandMom = mom->mother()) != nullptr)
82  mom = grandMom;
83  if (select(*mom)) {
84  moms.push_back(mom);
85  }
86  }
87  }
88  }
89  }
90  cout << "-- decay: --" << endl;
91  if (!moms.empty()) {
92  if (moms.size() > 1)
93  for (size_t m = 0; m < moms.size(); ++m) {
94  string dec = decay(*moms[m], skip);
95  if (!dec.empty())
96  cout << "{ " << dec << " } ";
97  }
98  else
99  cout << decay(*moms[0], skip);
100  }
101  if (!nodes.empty()) {
102  cout << "-> ";
103  if (nodes.size() > 1) {
104  for (size_t n = 0; n < nodes.size(); ++n) {
105  skip.remove(nodes[n]);
106  string dec = decay(*nodes[n], skip);
107  if (!dec.empty()) {
108  if (dec.find("->", 0) != string::npos)
109  cout << " ( " << dec << " )";
110  else
111  cout << " " << dec;
112  }
113  }
114  } else {
115  skip.remove(nodes[0]);
116  cout << decay(*nodes[0], skip);
117  }
118  }
119  cout << endl;
120 }
121 
123  ostringstream cout;
124  if (printP4_)
125  cout << " (" << c.px() << ", " << c.py() << ", " << c.pz() << "; " << c.energy() << ")";
126  if (printPtEtaPhi_)
127  cout << " [" << c.pt() << ": " << c.eta() << ", " << c.phi() << "]";
128  if (printVertex_)
129  cout << " {" << c.vx() << ", " << c.vy() << ", " << c.vz() << "}";
130  return cout.str();
131 }
132 
133 string ParticleDecayDrawer::decay(const Candidate &c, list<const Candidate *> &skip) const {
134  string out;
135  if (find(skip.begin(), skip.end(), &c) != skip.end())
136  return out;
137  skip.push_back(&c);
138 
139  int id = c.pdgId();
140  const ParticleData *pd = pdt_->particle(id);
141  assert(pd != nullptr);
142  out += (pd->name() + printP4(c));
143 
144  size_t validDau = 0, ndau = c.numberOfDaughters();
145  for (size_t i = 0; i < ndau; ++i)
146  if (accept(*c.daughter(i), skip))
147  ++validDau;
148  if (validDau == 0)
149  return out;
150 
151  out += " ->";
152 
153  for (size_t i = 0; i < ndau; ++i) {
154  const Candidate *d = c.daughter(i);
155  if (accept(*d, skip)) {
156  string dec = decay(*d, skip);
157  if (dec.find("->", 0) != string::npos)
158  out += (" ( " + dec + " )");
159  else
160  out += (" " + dec);
161  }
162  }
163  return out;
164 }
165 
167 
bool hasValidDaughters(const reco::Candidate &) const
has valid daughters in the chain
virtual double pz() const =0
z coordinate of momentum vector
virtual double vx() const =0
x coordinate of vertex position
bool accept(const reco::Candidate &, const std::list< const reco::Candidate * > &) const
accept candidate
std::string printP4(const reco::Candidate &) const
print 4 momenta
virtual const Candidate * daughter(size_type i) const =0
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
virtual const Candidate * mother(size_type i=0) const =0
return pointer to mother
virtual double vy() const =0
y coordinate of vertex position
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
edm::ESHandle< ParticleDataTable > pdt_
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
bool getData(T &iHolder) const
Definition: EventSetup.h:113
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::string decay(const reco::Candidate &, std::list< const reco::Candidate * > &) const
bool select(const reco::Candidate &) const
select candidate
virtual int pdgId() const =0
PDG identifier.
ParticleDecayDrawer(const edm::ParameterSet &)
bool printP4_
print parameters
void analyze(const edm::Event &, const edm::EventSetup &) override
HepPDT::ParticleData ParticleData
d
Definition: ztail.py:151
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
edm::EDGetTokenT< edm::View< reco::Candidate > > srcToken_
Definition: event.py:1