Validation
EventGenerator
plugins
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
7
#include "
Validation/EventGenerator/interface/BPhysicsSpectrum.h
"
8
9
#include "
FWCore/Framework/interface/MakerMacros.h
"
10
11
using namespace
edm
;
12
13
BPhysicsSpectrum::BPhysicsSpectrum
(
const
edm::ParameterSet
& iPSet)
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
23
BPhysicsSpectrum::~BPhysicsSpectrum
() {}
24
25
void
BPhysicsSpectrum::bookHistograms
(
DQMStore::IBooker
&
i
,
edm::Run
const
&,
edm::EventSetup
const
&) {
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
32
void
BPhysicsSpectrum::analyze
(
const
edm::Event
&
iEvent
,
const
edm::EventSetup
& iSetup) {
33
edm::Handle<reco::GenParticleCollection>
genParticles
;
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
}
mps_fire.i
i
Definition:
mps_fire.py:355
genParticles2HepMC_cfi.genParticles
genParticles
Definition:
genParticles2HepMC_cfi.py:4
edm::Run
Definition:
Run.h:45
edm
HLT enums.
Definition:
AlignableModifier.h:19
BPhysicsSpectrum::genparticleCollection_
edm::InputTag genparticleCollection_
Definition:
BPhysicsSpectrum.h:44
BPhysicsSpectrum::genparticleCollectionToken_
edm::EDGetTokenT< reco::GenParticleCollection > genparticleCollectionToken_
Definition:
BPhysicsSpectrum.h:45
BPhysicsSpectrum::mass_max
double mass_max
Definition:
BPhysicsSpectrum.h:47
BPhysicsSpectrum::BPhysicsSpectrum
BPhysicsSpectrum(const edm::ParameterSet &)
Definition:
BPhysicsSpectrum.cc:13
edm::Handle< reco::GenParticleCollection >
MakerMacros.h
dqm::impl::MonitorElement::Fill
void Fill(long long x)
Definition:
MonitorElement.h:290
BPhysicsSpectrum::name
std::string name
Definition:
BPhysicsSpectrum.h:46
BPhysicsSpectrum.h
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition:
AlCaHLTBitMon_QueryRunRegistry.py:256
edm::ParameterSet
Definition:
ParameterSet.h:36
BPhysicsSpectrum::analyze
void analyze(edm::Event const &, edm::EventSetup const &) override
Definition:
BPhysicsSpectrum.cc:32
iEvent
int iEvent
Definition:
GenABIO.cc:224
edm::EventSetup
Definition:
EventSetup.h:57
BPhysicsSpectrum::Particles
std::vector< int > Particles
Definition:
BPhysicsSpectrum.h:48
DQMHelper
Definition:
DQMHelper.h:15
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
BPhysicsSpectrum::~BPhysicsSpectrum
~BPhysicsSpectrum() override
Definition:
BPhysicsSpectrum.cc:23
std
Definition:
JetResolutionObject.h:76
BPhysicsSpectrum::mass
MonitorElement * mass
Definition:
BPhysicsSpectrum.h:43
Skims_PA_cff.name
name
Definition:
Skims_PA_cff.py:17
dqm::implementation::IBooker
Definition:
DQMStore.h:43
dqm
Definition:
DQMStore.h:18
funct::abs
Abs< T >::type abs(const T &t)
Definition:
Abs.h:22
BPhysicsSpectrum::mass_min
double mass_min
Definition:
BPhysicsSpectrum.h:47
edm::Event
Definition:
Event.h:73
BPhysicsSpectrum::Nobj
MonitorElement * Nobj
Definition:
BPhysicsSpectrum.h:43
edm::InputTag
Definition:
InputTag.h:15
BPhysicsSpectrum::bookHistograms
void bookHistograms(DQMStore::IBooker &i, edm::Run const &, edm::EventSetup const &) override
Definition:
BPhysicsSpectrum.cc:25
Generated for CMSSW Reference Manual by
1.8.16