CMS 3D CMS Logo

BPhysicsValidation.cc
Go to the documentation of this file.
1 //
3 // class Validation: Class to fill dqm monitor elements from existing EDM file
4 //
6 
8 
10 
11 using namespace edm;
12 
14  : genparticleCollection_(iPSet.getParameter<edm::InputTag>("genparticleCollection")),
15  // do not include weights right now to allow for running on aod
16  name(iPSet.getParameter<std::string>("name")),
17  particle(name, iPSet) {
18  genparticleCollectionToken_ = consumes<reco::GenParticleCollection>(genparticleCollection_);
19  std::vector<std::string> daughterNames = iPSet.getParameter<std::vector<std::string> >("daughters");
20  for (unsigned int i = 0; i < daughterNames.size(); i++) {
21  std::string curSet = daughterNames[i];
22  daughters.push_back(ParticleMonitor(name + curSet, iPSet.getUntrackedParameter<ParameterSet>(curSet)));
23  }
24 }
25 
27 
29 
31  DQMHelper dqm(&i);
32  i.setCurrentFolder("Generator/BPhysics");
33  Nobj = dqm.book1dHisto("N" + name, "N" + name, 1, 0., 1, "bin", "Number of " + name);
35  for (unsigned int j = 0; j < daughters.size(); j++) {
36  daughters[j].Configure(i);
37  }
38 }
39 
42  iEvent.getByToken(genparticleCollectionToken_, genParticles);
43  for (reco::GenParticleCollection::const_iterator iter = genParticles->begin(); iter != genParticles->end(); ++iter) {
44  if (abs(iter->pdgId()) == abs(particle.PDGID())) {
45  Nobj->Fill(0.5, 1.0);
46  particle.Fill(&(*iter), 1.0);
47  FillDaughters(&(*iter));
48  }
49  }
50 }
51 
53  int mpdgid = p->pdgId();
54  for (unsigned int i = 0; i < p->numberOfDaughters(); i++) {
55  const reco::GenParticle* dau = static_cast<const reco::GenParticle*>(p->daughter(i));
56  int pdgid = dau->pdgId();
57  for (unsigned int i = 0; i < daughters.size(); i++) {
58  if (abs(mpdgid) != abs(daughters[i].PDGID()) && daughters[i].PDGID() == pdgid)
59  daughters[i].Fill(dau, 1.0);
60  // note: use abs when comparing to mother to avoid mixing
61  }
62  FillDaughters(dau);
63  }
64 }
void Configure(DQMStore::IBooker &i)
T getParameter(std::string const &) const
int pdgId() const final
PDG identifier.
const Candidate * daughter(size_type) const override
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
T getUntrackedParameter(std::string const &, T const &) const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
MonitorElement * Nobj
BPhysicsValidation(const edm::ParameterSet &)
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:418
~BPhysicsValidation() override
edm::EDGetTokenT< reco::GenParticleCollection > genparticleCollectionToken_
void Fill(long long x)
int iEvent
Definition: GenABIO.cc:224
MonitorElement * book1dHisto(std::string name, std::string title, int n, double xmin, double xmax, std::string xaxis, std::string yaxis)
Definition: DQMHelper.cc:7
ParticleMonitor particle
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
void bookHistograms(DQMStore::IBooker &i, edm::Run const &, edm::EventSetup const &) override
edm::InputTag genparticleCollection_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void dqmBeginRun(const edm::Run &r, const edm::EventSetup &c) override
size_t numberOfDaughters() const override
number of daughters
void FillDaughters(const reco::GenParticle *p)
HLT enums.
void Fill(const reco::GenParticle *p, double weight)
std::vector< ParticleMonitor > daughters
Definition: Run.h:45
void analyze(edm::Event const &, edm::EventSetup const &) override