Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 11:27:17

0001 ///////////////////////////////////////
0002 //
0003 // class Validation: Class to fill dqm monitor elements from existing EDM file
0004 //
0005 ///////////////////////////////////////
0006 
0007 #include "Validation/EventGenerator/interface/BPhysicsValidation.h"
0008 
0009 #include "FWCore/Framework/interface/MakerMacros.h"
0010 
0011 using namespace edm;
0012 
0013 BPhysicsValidation::BPhysicsValidation(const edm::ParameterSet& iPSet)
0014     : genparticleCollection_(iPSet.getParameter<edm::InputTag>("genparticleCollection")),
0015       // do not include weights right now to allow for running on aod
0016       name(iPSet.getParameter<std::string>("name")),
0017       particle(name, iPSet) {
0018   genparticleCollectionToken_ = consumes<reco::GenParticleCollection>(genparticleCollection_);
0019   std::vector<std::string> daughterNames = iPSet.getParameter<std::vector<std::string> >("daughters");
0020   for (unsigned int i = 0; i < daughterNames.size(); i++) {
0021     std::string curSet = daughterNames[i];
0022     daughters.push_back(ParticleMonitor(name + curSet, iPSet.getUntrackedParameter<ParameterSet>(curSet)));
0023   }
0024 }
0025 
0026 BPhysicsValidation::~BPhysicsValidation() {}
0027 
0028 void BPhysicsValidation::bookHistograms(DQMStore::IBooker& i, edm::Run const&, edm::EventSetup const&) {
0029   DQMHelper dqm(&i);
0030   i.setCurrentFolder("Generator/BPhysics");
0031   Nobj = dqm.book1dHisto("N" + name, "N" + name, 1, 0., 1, "bin", "Number of " + name);
0032   particle.Configure(i);
0033   for (unsigned int j = 0; j < daughters.size(); j++) {
0034     daughters[j].Configure(i);
0035   }
0036 }
0037 
0038 void BPhysicsValidation::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0039   edm::Handle<reco::GenParticleCollection> genParticles;
0040   iEvent.getByToken(genparticleCollectionToken_, genParticles);
0041   for (reco::GenParticleCollection::const_iterator iter = genParticles->begin(); iter != genParticles->end(); ++iter) {
0042     if (abs(iter->pdgId()) == abs(particle.PDGID())) {
0043       Nobj->Fill(0.5, 1.0);
0044       particle.Fill(&(*iter), 1.0);
0045       FillDaughters(&(*iter));
0046     }
0047   }
0048 }
0049 
0050 void BPhysicsValidation::FillDaughters(const reco::GenParticle* p) {
0051   int mpdgid = p->pdgId();
0052   for (unsigned int i = 0; i < p->numberOfDaughters(); i++) {
0053     const reco::GenParticle* dau = static_cast<const reco::GenParticle*>(p->daughter(i));
0054     int pdgid = dau->pdgId();
0055     for (unsigned int i = 0; i < daughters.size(); i++) {
0056       if (abs(mpdgid) != abs(daughters[i].PDGID()) && daughters[i].PDGID() == pdgid)
0057         daughters[i].Fill(dau, 1.0);
0058       // note: use abs when comparing to mother to avoid mixing
0059     }
0060     FillDaughters(dau);
0061   }
0062 }