CMS 3D CMS Logo

ParticleListDrawer.cc
Go to the documentation of this file.
1 #include <memory>
2 #include <string>
3 #include <iostream>
4 #include <sstream>
5 
19 
37 using namespace std;
38 using namespace reco;
39 using namespace edm;
40 
42 public:
43  explicit ParticleListDrawer(const edm::ParameterSet&);
44  ~ParticleListDrawer() override{};
45  void analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) override;
46 
47 private:
48  std::string getParticleName(int id) const;
49 
53  int maxEventsToPrint_; // Must be signed, because -1 is used for no limit
54  unsigned int nEventAnalyzed_;
59 };
60 
62  : src_(pset.getParameter<InputTag>("src")),
63  srcToken_(consumes<reco::CandidateView>(src_)),
64  maxEventsToPrint_(pset.getUntrackedParameter<int>("maxEventsToPrint", 1)),
65  nEventAnalyzed_(0),
66  printOnlyHardInteraction_(pset.getUntrackedParameter<bool>("printOnlyHardInteraction", false)),
67  printVertex_(pset.getUntrackedParameter<bool>("printVertex", false)),
68  printFlags_(pset.getUntrackedParameter<bool>("printFlags", false)),
69  useMessageLogger_(pset.getUntrackedParameter<bool>("useMessageLogger", false)) {}
70 
72  const ParticleData* pd = pdt_->particle(id);
73  if (!pd) {
74  std::ostringstream ss;
75  ss << "P" << id;
76  return ss.str();
77  } else
78  return pd->name();
79 }
80 
83  iEvent.getByToken(srcToken_, particles);
84  iSetup.getData(pdt_);
85 
86  if (maxEventsToPrint_ < 0 || nEventAnalyzed_ < static_cast<unsigned int>(maxEventsToPrint_)) {
87  ostringstream out;
88  char buf[256];
89 
90  out << endl << "[ParticleListDrawer] analysing particle collection " << src_.label() << endl;
91 
92  snprintf(buf,
93  sizeof(buf),
94  " idx | ID - Name |Stat| Mo1 Mo2 Da1 Da2 |nMo nDa| pt eta phi | px "
95  " py pz m |");
96  out << buf;
97  if (printVertex_) {
98  snprintf(buf, sizeof(buf), " vx vy vz |");
99  out << buf;
100  }
101  out << endl;
102 
103  int idx = -1;
104  int iMo1 = -1;
105  int iMo2 = -1;
106  int iDa1 = -1;
107  int iDa2 = -1;
108  vector<const reco::Candidate*> cands;
109  vector<const Candidate*>::const_iterator found = cands.begin();
110  for (CandidateView::const_iterator p = particles->begin(); p != particles->end(); ++p) {
111  cands.push_back(&*p);
112  }
113 
114  for (CandidateView::const_iterator p = particles->begin(); p != particles->end(); p++) {
115  if (printOnlyHardInteraction_ && p->status() != 3)
116  continue;
117 
118  // Particle Name
119  int id = p->pdgId();
120  string particleName = getParticleName(id);
121 
122  // Particle Index
123  idx = p - particles->begin();
124 
125  // Particles Mothers and Daighters
126  iMo1 = -1;
127  iMo2 = -1;
128  iDa1 = -1;
129  iDa2 = -1;
130  int nMo = p->numberOfMothers();
131  int nDa = p->numberOfDaughters();
132 
133  found = find(cands.begin(), cands.end(), p->mother(0));
134  if (found != cands.end())
135  iMo1 = found - cands.begin();
136 
137  found = find(cands.begin(), cands.end(), p->mother(nMo - 1));
138  if (found != cands.end())
139  iMo2 = found - cands.begin();
140 
141  found = find(cands.begin(), cands.end(), p->daughter(0));
142  if (found != cands.end())
143  iDa1 = found - cands.begin();
144 
145  found = find(cands.begin(), cands.end(), p->daughter(nDa - 1));
146  if (found != cands.end())
147  iDa2 = found - cands.begin();
148 
149  char buf[2400];
150  snprintf(
151  buf,
152  sizeof(buf),
153  " %4d | %5d - %10s | %2d | %4d %4d %4d %4d | %2d %2d | %7.3f %10.3f %6.3f | %10.3f %10.3f %10.3f %8.3f |",
154  idx,
155  p->pdgId(),
156  particleName.c_str(),
157  p->status(),
158  iMo1,
159  iMo2,
160  iDa1,
161  iDa2,
162  nMo,
163  nDa,
164  p->pt(),
165  p->eta(),
166  p->phi(),
167  p->px(),
168  p->py(),
169  p->pz(),
170  p->mass());
171  out << buf;
172 
173  if (printVertex_) {
174  snprintf(buf, sizeof(buf), " %10.3f %10.3f %10.3f |", p->vertex().x(), p->vertex().y(), p->vertex().z());
175  out << buf;
176  }
177 
178  if (printFlags_) {
179  const reco::GenParticle* gp = dynamic_cast<const reco::GenParticle*>(&*p);
180  if (!gp)
181  throw cms::Exception("Unsupported", "Status flags can be printed only for reco::GenParticle objects\n");
182  if (gp->isPromptFinalState())
183  out << " PromptFinalState";
185  out << " DirectPromptTauDecayProductFinalState";
186  if (gp->isHardProcess())
187  out << " HardProcess";
188  if (gp->fromHardProcessFinalState())
189  out << " HardProcessFinalState";
190  if (gp->fromHardProcessBeforeFSR())
191  out << " HardProcessBeforeFSR";
192  if (gp->statusFlags().isFirstCopy())
193  out << " FirstCopy";
194  if (gp->isLastCopy())
195  out << " LastCopy";
196  if (gp->isLastCopyBeforeFSR())
197  out << " LastCopyBeforeFSR";
198  }
199 
200  out << endl;
201  }
202  nEventAnalyzed_++;
203 
204  if (useMessageLogger_)
205  LogVerbatim("ParticleListDrawer") << out.str();
206  else
207  cout << out.str();
208  }
209 }
210 
bool isPromptFinalState() const
Definition: GenParticle.h:51
bool isLastCopyBeforeFSR() const
Definition: GenParticle.h:106
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
bool isLastCopy() const
Definition: GenParticle.h:101
bool isHardProcess() const
Definition: GenParticle.h:72
bool fromHardProcessBeforeFSR() const
Definition: GenParticle.h:94
bool isFirstCopy() const
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
example_stream void analyze(const edm::Event &, const edm::EventSetup &) override
bool getData(T &iHolder) const
Definition: EventSetup.h:113
edm::EDGetTokenT< reco::CandidateView > srcToken_
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
const_iterator begin() const
ParticleListDrawer(const edm::ParameterSet &)
unsigned int nEventAnalyzed_
HepPDT::ParticleData ParticleData
void analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) override
const GenStatusFlags & statusFlags() const
Definition: GenParticle.h:38
bool fromHardProcessFinalState() const
Definition: GenParticle.h:75
std::string const & label() const
Definition: InputTag.h:36
edm::ESHandle< ParticleDataTable > pdt_
fixed size matrix
HLT enums.
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
std::string getParticleName(int id) const
bool isDirectPromptTauDecayProductFinalState() const
Definition: GenParticle.h:59
~ParticleListDrawer() override
const_iterator end() const
Module to analyze the particle listing as provided by common event generators.