00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014 #include "FWCore/Framework/interface/Frameworkfwd.h"
00015 #include "FWCore/Framework/interface/EDAnalyzer.h"
00016 #include "FWCore/Framework/interface/ESHandle.h"
00017 #include "FWCore/Framework/interface/Event.h"
00018 #include "FWCore/Framework/interface/EventSetup.h"
00019 #include "FWCore/Framework/interface/MakerMacros.h"
00020 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00021 #include "FWCore/ParameterSet/interface/InputTag.h"
00022 #include "SimGeneral/HepPDTRecord/interface/ParticleDataTable.h"
00023 #include "DataFormats/Candidate/interface/Candidate.h"
00024 #include "DataFormats/Candidate/interface/CandidateFwd.h"
00025 #include "DataFormats/Common/interface/Ref.h"
00026
00027
00028 #include <memory>
00029 #include <string>
00030 #include <iostream>
00031 #include <sstream>
00032
00033
00034 using namespace std;
00035 using namespace reco;
00036 using namespace edm;
00037
00038 class ParticleListDrawer : public edm::EDAnalyzer {
00039 public:
00040 explicit ParticleListDrawer(const edm::ParameterSet & );
00041 ~ParticleListDrawer() {};
00042 void analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup);
00043
00044 private:
00045 std::string getParticleName( int id ) const;
00046
00047 edm::InputTag src_;
00048 edm::ESHandle<ParticleDataTable> pdt_;
00049 unsigned int maxEventsToPrint_;
00050 unsigned int nEventAnalyzed_;
00051 bool printOnlyHardInteraction_;
00052 bool useMessageLogger_;
00053 };
00054
00055 ParticleListDrawer::ParticleListDrawer(const edm::ParameterSet & pset) :
00056 src_(pset.getParameter<InputTag>("src")),
00057 maxEventsToPrint_ (pset.getUntrackedParameter<int>("maxEventsToPrint",1)),
00058 nEventAnalyzed_(0),
00059 printOnlyHardInteraction_(pset.getUntrackedParameter<bool>("printOnlyHardInteraction", false)),
00060 useMessageLogger_(pset.getUntrackedParameter<bool>("useMessageLogger", false)) {
00061 }
00062
00063 std::string ParticleListDrawer::getParticleName(int id) const
00064 {
00065 const ParticleData * pd = pdt_->particle( id );
00066 if (!pd) {
00067 std::ostringstream ss;
00068 ss << "P" << id;
00069 return ss.str();
00070 } else
00071 return pd->name();
00072 }
00073
00074 void ParticleListDrawer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
00075 Handle<reco::CandidateView> particles;
00076 iEvent.getByLabel (src_, particles );
00077 if (!particles.isValid()) {
00078 cerr << "[ParticleListDrawer] caught std::exception " << endl;
00079 return;
00080 } else {
00081 iSetup.getData( pdt_ );
00082 }
00083
00084 if(maxEventsToPrint_ < 0 || nEventAnalyzed_ < maxEventsToPrint_) {
00085 ostringstream out;
00086 char buf[256];
00087
00088 out << "[ParticleListDrawer] analysing event " << iEvent.id() << endl;
00089
00090 out << endl;
00091 out << "**********************" << endl;
00092 out << "* GenEvent *" << endl;
00093 out << "**********************" << endl;
00094
00095 snprintf(buf, 256, " idx | ID - Name |Stat| Mo1 Mo2 Da1 Da2 |nMo nDa| pt eta phi | px py pz m |\n");
00096 out << buf;
00097
00098 int idx = -1;
00099 int iMo1 = -1;
00100 int iMo2 = -1;
00101 int iDa1 = -1;
00102 int iDa2 = -1;
00103 vector<const reco::Candidate *> cands;
00104 vector<const Candidate *>::const_iterator found = cands.begin();
00105 for(CandidateView::const_iterator p = particles->begin();
00106 p != particles->end(); ++ p) {
00107 cands.push_back(&*p);
00108 }
00109
00110 for(CandidateView::const_iterator p = particles->begin();
00111 p != particles->end();
00112 p ++) {
00113 if (printOnlyHardInteraction_ && p->status() != 3) continue;
00114
00115
00116 int id = p->pdgId();
00117 string particleName = getParticleName(id);
00118
00119
00120 idx = p - particles->begin();
00121
00122
00123 iMo1 = -1;
00124 iMo2 = -1;
00125 iDa1 = -1;
00126 iDa2 = -1;
00127 int nMo = p->numberOfMothers();
00128 int nDa = p->numberOfDaughters();
00129
00130 found = find(cands.begin(), cands.end(), p->mother(0));
00131 if(found != cands.end()) iMo1 = found - cands.begin() ;
00132
00133 found = find(cands.begin(), cands.end(), p->mother(nMo-1));
00134 if(found != cands.end()) iMo2 = found - cands.begin() ;
00135
00136 found = find(cands.begin(), cands.end(), p->daughter(0));
00137 if(found != cands.end()) iDa1 = found - cands.begin() ;
00138
00139 found = find(cands.begin(), cands.end(), p->daughter(nDa-1));
00140 if(found != cands.end()) iDa2 = found - cands.begin() ;
00141
00142 char buf[256];
00143 snprintf(buf, 256,
00144 " %4d | %5d - %10s | %2d | %4d %4d %4d %4d | %2d %2d | %7.3f %10.3f %6.3f | %10.3f %10.3f %10.3f %8.3f |\n",
00145 idx,
00146 p->pdgId(),
00147 particleName.c_str(),
00148 p->status(),
00149 iMo1,iMo2,iDa1,iDa2,nMo,nDa,
00150 p->pt(),
00151 p->eta(),
00152 p->phi(),
00153 p->px(),
00154 p->py(),
00155 p->pz(),
00156 p->mass()
00157 );
00158 out << buf;
00159 }
00160 nEventAnalyzed_++;
00161
00162 if (useMessageLogger_)
00163 LogVerbatim("ParticleListDrawer") << out.str();
00164 else
00165 cout << out.str();
00166 }
00167 }
00168
00169 DEFINE_FWK_MODULE(ParticleListDrawer);
00170