CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ParticleListDrawer.cc
Go to the documentation of this file.
1 #include <memory>
2 #include <string>
3 #include <iostream>
4 #include <sstream>
5 
18 
37 using namespace std;
38 using namespace reco;
39 using namespace edm;
40 
42  public:
43  explicit ParticleListDrawer(const edm::ParameterSet & );
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_;
58 };
59 
61  src_(pset.getParameter<InputTag>("src")),
62  srcToken_(consumes<reco::CandidateView>(src_)),
63  maxEventsToPrint_ (pset.getUntrackedParameter<int>("maxEventsToPrint",1)),
64  nEventAnalyzed_(0),
65  printOnlyHardInteraction_(pset.getUntrackedParameter<bool>("printOnlyHardInteraction", false)),
66  printVertex_(pset.getUntrackedParameter<bool>("printVertex", false)),
67  useMessageLogger_(pset.getUntrackedParameter<bool>("useMessageLogger", false)) {
68 }
69 
71 {
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
91  << "[ParticleListDrawer] analysing particle collection " << src_.label() << endl;
92 
93  snprintf(buf, 256, " idx | ID - Name |Stat| Mo1 Mo2 Da1 Da2 |nMo nDa| pt eta phi | px py pz m |");
94  out << buf;
95  if (printVertex_) {
96  snprintf(buf, 256, " vx vy vz |");
97  out << buf;
98  }
99  out << endl;
100 
101  int idx = -1;
102  int iMo1 = -1;
103  int iMo2 = -1;
104  int iDa1 = -1;
105  int iDa2 = -1;
106  vector<const reco::Candidate *> cands;
107  vector<const Candidate *>::const_iterator found = cands.begin();
108  for(CandidateView::const_iterator p = particles->begin();
109  p != particles->end(); ++ p) {
110  cands.push_back(&*p);
111  }
112 
113  for(CandidateView::const_iterator p = particles->begin();
114  p != particles->end();
115  p ++) {
116  if (printOnlyHardInteraction_ && p->status() != 3) 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()) iMo1 = found - cands.begin() ;
135 
136  found = find(cands.begin(), cands.end(), p->mother(nMo-1));
137  if(found != cands.end()) iMo2 = found - cands.begin() ;
138 
139  found = find(cands.begin(), cands.end(), p->daughter(0));
140  if(found != cands.end()) iDa1 = found - cands.begin() ;
141 
142  found = find(cands.begin(), cands.end(), p->daughter(nDa-1));
143  if(found != cands.end()) iDa2 = found - cands.begin() ;
144 
145  char buf[256];
146  snprintf(buf, 256,
147  " %4d | %5d - %10s | %2d | %4d %4d %4d %4d | %2d %2d | %7.3f %10.3f %6.3f | %10.3f %10.3f %10.3f %8.3f |",
148  idx,
149  p->pdgId(),
150  particleName.c_str(),
151  p->status(),
152  iMo1,iMo2,iDa1,iDa2,nMo,nDa,
153  p->pt(),
154  p->eta(),
155  p->phi(),
156  p->px(),
157  p->py(),
158  p->pz(),
159  p->mass()
160  );
161  out << buf;
162 
163  if (printVertex_) {
164  snprintf(buf, 256, " %10.3f %10.3f %10.3f |",
165  p->vertex().x(),
166  p->vertex().y(),
167  p->vertex().z());
168  out << buf;
169  }
170 
171  out << endl;
172  }
173  nEventAnalyzed_++;
174 
175  if (useMessageLogger_)
176  LogVerbatim("ParticleListDrawer") << out.str();
177  else
178  cout << out.str();
179  }
180 }
181 
183 
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:81
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
void getData(T &iHolder) const
Definition: EventSetup.h:67
edm::EDGetTokenT< reco::CandidateView > srcToken_
int iEvent
Definition: GenABIO.cc:243
ParticleListDrawer(const edm::ParameterSet &)
unsigned int nEventAnalyzed_
HepPDT::ParticleData ParticleData
void analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) override
tuple out
Definition: dbtoconf.py:99
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
std::string const & label() const
Definition: InputTag.h:42
edm::ESHandle< ParticleDataTable > pdt_
std::string getParticleName(int id) const
tuple cout
Definition: gather_cfg.py:121
volatile std::atomic< bool > shutdown_flag false
Module to analyze the particle listing as provided by common event generators.