Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:23:25

0001 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0002 #include "FWCore/Framework/interface/Event.h"
0003 #include "FWCore/Framework/interface/EventSetup.h"
0004 #include "DataFormats/Common/interface/ValidHandle.h"
0005 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0006 #include "FWCore/Utilities/interface/EDGetToken.h"
0007 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0008 #include "FWCore/Utilities/interface/InputTag.h"
0009 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
0010 #include "SimGeneral/HepPDTRecord/interface/ParticleDataTable.h"
0011 #include "FWCore/Utilities/interface/EDMException.h"
0012 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0013 
0014 class DummyHepMCAnalyzer : public edm::one::EDAnalyzer<> {
0015 private:
0016   const bool dumpHepMC_;
0017   const bool dumpPDF_;
0018   const bool checkPDG_;
0019   const edm::EDGetTokenT<edm::HepMCProduct> srcToken_;
0020   const edm::ESGetToken<HepPDT::ParticleDataTable, PDTRecord> pdtToken_;
0021 
0022 public:
0023   explicit DummyHepMCAnalyzer(const edm::ParameterSet& cfg)
0024       : dumpHepMC_(cfg.getUntrackedParameter<bool>("dumpHepMC")),
0025         dumpPDF_(cfg.getUntrackedParameter<bool>("dumpPDF")),
0026         checkPDG_(cfg.getUntrackedParameter<bool>("checkPDG")),
0027         srcToken_(consumes<edm::HepMCProduct>(cfg.getParameter<edm::InputTag>("src"))),
0028         pdtToken_(esConsumes<HepPDT::ParticleDataTable, PDTRecord>()) {}
0029 
0030   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0031     edm::ParameterSetDescription desc;
0032     desc.addUntracked<bool>("dumpHepMC", true);
0033     desc.addUntracked<bool>("dumpPDF", false);
0034     desc.addUntracked<bool>("checkPDG", false);
0035     desc.add<edm::InputTag>("src", edm::InputTag("generatorSmeared"))
0036         ->setComment("Input generated HepMC event after vtx smearing");
0037     descriptions.add("dummyHepMCAnalyzer", desc);
0038   }
0039 
0040   void analyze(const edm::Event& evt, const edm::EventSetup& es) override {
0041     auto hepMC = makeValid(evt.getHandle(srcToken_));
0042 
0043     const HepMC::GenEvent* mc = hepMC->GetEvent();
0044     edm::LogPrint("HepMCAnalyzer") << "\n particles #: " << mc->particles_size();
0045 
0046     if (dumpPDF_) {
0047       edm::LogPrint("HepMCAnalyzer") << "\n PDF info: " << mc->pdf_info();
0048     }
0049 
0050     if (dumpHepMC_) {
0051       mc->print(std::cout);
0052     }
0053 
0054     if (checkPDG_) {
0055       const auto& fPDGTable = es.getHandle(pdtToken_);
0056       for (HepMC::GenEvent::particle_const_iterator part = mc->particles_begin(); part != mc->particles_end(); ++part) {
0057         const HepPDT::ParticleData* PData = fPDGTable->particle(HepPDT::ParticleID((*part)->pdg_id()));
0058         if (!PData) {
0059           edm::LogWarning("HepMCAnalyzer") << "Missing entry in particle table for PDG code = " << (*part)->pdg_id();
0060         }
0061       }
0062     }
0063   }
0064 };
0065 
0066 #include "FWCore/Framework/interface/MakerMacros.h"
0067 
0068 DEFINE_FWK_MODULE(DummyHepMCAnalyzer);