CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
BPhysicsValidation.h
Go to the documentation of this file.
1 #ifndef BPhysicsValidation_H
2 #define BPhysicsValidation_H
3 
4 /*class BPhysicsValidation
5  *
6  * Class to fill Event Generator dqm monitor elements; works on HepMCProduct
7  *
8  *
9  */
10 #include <iostream>
11 #include "TMath.h"
12 // framework & common header files
13 
17 
22 
23 //DQM services
27 
31 
33 
35 public:
36  explicit BPhysicsValidation(const edm::ParameterSet &);
37  ~BPhysicsValidation() override;
38 
39  void bookHistograms(DQMStore::IBooker &i, edm::Run const &, edm::EventSetup const &) override;
40  void analyze(edm::Event const &, edm::EventSetup const &) override;
41 
42 private:
44  public:
46  : p(p_), name(name_), pdgid(p.getParameter<int>("pdgid")){};
48 
51  double mass_min = p.getParameter<double>("massmin");
52  double mass_max = p.getParameter<double>("massmax");
53  DQMHelper dqm(&i);
54  i.setCurrentFolder("Generator/BPhysics");
55  // Number of analyzed events
56  pt = dqm.book1dHisto(name + "PT", "P_{t} of the " + pname + "s", 100, 0., 100, "P_{t} (GeV)", "Number of Events");
57  eta = dqm.book1dHisto(name + "ETA", "#eta of the " + pname + "s", 100, -5., 5., "#eta", "Number of Events");
58  phi = dqm.book1dHisto(
59  name + "PHI", "#phi of the " + pname + "s", 100, 0, 2 * TMath::Pi(), "#phi", "Number of Events");
60  mass = dqm.book1dHisto(
61  name + "MASS", "Mass of the " + pname + "s", 100, mass_min, mass_max, "Mass (GeV)", "Number of Events");
62  }
63 
64  void Fill(const reco::GenParticle *p, double weight) {
65  if (abs(p->pdgId()) == abs(pdgid)) {
66  pt->Fill(p->pt(), weight);
67  eta->Fill(p->eta(), weight);
68  phi->Fill(p->phi(), weight);
69  mass->Fill(p->mass(), weight);
70  }
71  }
72  int PDGID() { return pdgid; }
73 
74  private:
77  int pdgid;
79  };
80 
81  void FillDaughters(const reco::GenParticle *p);
86  std::vector<ParticleMonitor> daughters;
88 };
89 
90 #endif
const double Pi
void Configure(DQMStore::IBooker &i)
ParticleMonitor(std::string name_, const edm::ParameterSet &p_)
double pt() const final
transverse momentum
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
MonitorElement * Nobj
BPhysicsValidation(const edm::ParameterSet &)
~BPhysicsValidation() override
edm::EDGetTokenT< reco::GenParticleCollection > genparticleCollectionToken_
int pdgId() const final
PDG identifier.
void Fill(long long x)
ParticleMonitor particle
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
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
void FillDaughters(const reco::GenParticle *p)
double mass() const final
mass
void Fill(const reco::GenParticle *p, double weight)
std::vector< ParticleMonitor > daughters
int weight
Definition: histoStyle.py:51
double phi() const final
momentum azimuthal angle
Definition: Run.h:45
void analyze(edm::Event const &, edm::EventSetup const &) override
double eta() const final
momentum pseudorapidity