CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
BPhysicsSpectrum Class Reference

#include <BPhysicsSpectrum.h>

Inheritance diagram for BPhysicsSpectrum:
one::DQMEDAnalyzer< T > one::dqmimplementation::DQMBaseClass< T... >

Public Member Functions

void analyze (edm::Event const &, edm::EventSetup const &) override
 
void bookHistograms (DQMStore::IBooker &i, edm::Run const &, edm::EventSetup const &) override
 
 BPhysicsSpectrum (const edm::ParameterSet &)
 
void dqmBeginRun (const edm::Run &r, const edm::EventSetup &c) override
 
 ~BPhysicsSpectrum () override
 
- Public Member Functions inherited from one::DQMEDAnalyzer< T >
 DQMEDAnalyzer ()=default
 
 DQMEDAnalyzer (DQMEDAnalyzer< T... > const &)=delete
 
 DQMEDAnalyzer (DQMEDAnalyzer< T... > &&)=delete
 
 ~DQMEDAnalyzer () override=default
 

Private Attributes

edm::InputTag genparticleCollection_
 
edm::EDGetTokenT< reco::GenParticleCollectiongenparticleCollectionToken_
 
MonitorElementmass
 
double mass_max
 
double mass_min
 
std::string name
 
MonitorElementNobj
 
std::vector< int > Particles
 

Detailed Description

Definition at line 37 of file BPhysicsSpectrum.h.

Constructor & Destructor Documentation

BPhysicsSpectrum::BPhysicsSpectrum ( const edm::ParameterSet iPSet)
explicit

Definition at line 13 of file BPhysicsSpectrum.cc.

References genparticleCollection_, genparticleCollectionToken_, edm::ParameterSet::getParameter(), and Particles.

13  :
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 {
20  genparticleCollectionToken_=consumes<reco::GenParticleCollection>(genparticleCollection_);
21  Particles = iPSet.getParameter<std::vector<int> >("pdgids");
22 }
T getParameter(std::string const &) const
edm::EDGetTokenT< reco::GenParticleCollection > genparticleCollectionToken_
edm::InputTag genparticleCollection_
std::vector< int > Particles
BPhysicsSpectrum::~BPhysicsSpectrum ( )
override

Definition at line 24 of file BPhysicsSpectrum.cc.

24 {}

Member Function Documentation

void BPhysicsSpectrum::analyze ( edm::Event const &  ,
edm::EventSetup const &   
)
override

Definition at line 34 of file BPhysicsSpectrum.cc.

References funct::abs(), MonitorElement::Fill(), genparticleCollectionToken_, GenHFHadronMatcher_cfi::genParticles, edm::Event::getByToken(), mps_fire::i, mass, Nobj, and Particles.

34  {
36  iEvent.getByToken(genparticleCollectionToken_, genParticles );
37  for (reco::GenParticleCollection::const_iterator iter = genParticles->begin(); iter != genParticles->end(); ++iter) {
38  for(unsigned int i=0;i<Particles.size();i++){
39  if(abs(iter->pdgId())==abs(Particles[i])){
40  Nobj->Fill(0.5,1.0);
41  mass->Fill(iter->mass(),1.0);
42  }
43  }
44  }
45 }
MonitorElement * Nobj
MonitorElement * mass
void Fill(long long x)
int iEvent
Definition: GenABIO.cc:230
edm::EDGetTokenT< reco::GenParticleCollection > genparticleCollectionToken_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< int > Particles
void BPhysicsSpectrum::bookHistograms ( DQMStore::IBooker i,
edm::Run const &  ,
edm::EventSetup const &   
)
override

Definition at line 28 of file BPhysicsSpectrum.cc.

References DQMHelper::book1dHisto(), mass, mass_max, mass_min, name, Nobj, and DQMStore::IBooker::setCurrentFolder().

28  {
29  DQMHelper dqm(&i); i.setCurrentFolder("Generator/BPhysics");
30  Nobj = dqm.book1dHisto("NSpectrum"+name, "NSpectrum"+name, 1, 0., 1,"bin","Number of "+name);
31  mass = dqm.book1dHisto(name+"Mass","Mass Spectrum", 100, mass_min, mass_max,"Mass (GeV)","Number of Events");
32 }
MonitorElement * Nobj
MonitorElement * mass
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:268
void BPhysicsSpectrum::dqmBeginRun ( const edm::Run r,
const edm::EventSetup c 
)
override

Definition at line 26 of file BPhysicsSpectrum.cc.

26 {}

Member Data Documentation

edm::InputTag BPhysicsSpectrum::genparticleCollection_
private

Definition at line 48 of file BPhysicsSpectrum.h.

Referenced by BPhysicsSpectrum().

edm::EDGetTokenT<reco::GenParticleCollection> BPhysicsSpectrum::genparticleCollectionToken_
private

Definition at line 49 of file BPhysicsSpectrum.h.

Referenced by analyze(), and BPhysicsSpectrum().

MonitorElement* BPhysicsSpectrum::mass
private
double BPhysicsSpectrum::mass_max
private

Definition at line 51 of file BPhysicsSpectrum.h.

Referenced by bookHistograms().

double BPhysicsSpectrum::mass_min
private

Definition at line 51 of file BPhysicsSpectrum.h.

Referenced by bookHistograms().

std::string BPhysicsSpectrum::name
private
MonitorElement * BPhysicsSpectrum::Nobj
private

Definition at line 47 of file BPhysicsSpectrum.h.

Referenced by analyze(), and bookHistograms().

std::vector<int> BPhysicsSpectrum::Particles
private

Definition at line 52 of file BPhysicsSpectrum.h.

Referenced by analyze(), and BPhysicsSpectrum().