CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
BPhysicsSpectrum.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  mass_min(iPSet.getParameter<double>("massmin")),
18  mass_max(iPSet.getParameter<double>("massmax")) {
19  genparticleCollectionToken_ = consumes<reco::GenParticleCollection>(genparticleCollection_);
20  Particles = iPSet.getParameter<std::vector<int> >("pdgids");
21 }
22 
24 
26  DQMHelper dqm(&i);
27  i.setCurrentFolder("Generator/BPhysics");
28  Nobj = dqm.book1dHisto("NSpectrum" + name, "NSpectrum" + name, 1, 0., 1, "bin", "Number of " + name);
29  mass = dqm.book1dHisto(name + "Mass", "Mass Spectrum", 100, mass_min, mass_max, "Mass (GeV)", "Number of Events");
30 }
31 
34  iEvent.getByToken(genparticleCollectionToken_, genParticles);
35  for (reco::GenParticleCollection::const_iterator iter = genParticles->begin(); iter != genParticles->end(); ++iter) {
36  for (unsigned int i = 0; i < Particles.size(); i++) {
37  if (abs(iter->pdgId()) == abs(Particles[i])) {
38  Nobj->Fill(0.5, 1.0);
39  mass->Fill(iter->mass(), 1.0);
40  }
41  }
42  }
43 }
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
MonitorElement * Nobj
MonitorElement * mass
void bookHistograms(DQMStore::IBooker &i, edm::Run const &, edm::EventSetup const &) override
void Fill(long long x)
int iEvent
Definition: GenABIO.cc:224
edm::EDGetTokenT< reco::GenParticleCollection > genparticleCollectionToken_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
BPhysicsSpectrum(const edm::ParameterSet &)
edm::InputTag genparticleCollection_
~BPhysicsSpectrum() override
MonitorElement * book1dHisto(const std::string &name, const std::string &title, int n, double xmin, double xmax, const std::string &xaxis, const std::string &yaxis)
Definition: DQMHelper.cc:7
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
std::vector< int > Particles
void analyze(edm::Event const &, edm::EventSetup const &) override
Definition: Run.h:45