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