src
Validation
EventGenerator
plugins
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
7
#include "
Validation/EventGenerator/interface/BPhysicsValidation.h
"
8
9
#include "
FWCore/Framework/interface/MakerMacros.h
"
10
11
using namespace
edm
;
12
13
BPhysicsValidation::BPhysicsValidation
(
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
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
26
BPhysicsValidation::~BPhysicsValidation
() {}
27
28
void
BPhysicsValidation::bookHistograms
(
DQMStore::IBooker
&
i
,
edm::Run
const
&,
edm::EventSetup
const
&) {
29
DQMHelper
dqm
(&
i
);
30
i
.setCurrentFolder(
"Generator/BPhysics"
);
31
Nobj
=
dqm
.book1dHisto(
"N"
+
name
,
"N"
+
name
, 1, 0., 1,
"bin"
,
"Number of "
+
name
);
32
particle
.
Configure
(
i
);
33
for
(
unsigned
int
j
= 0;
j
<
daughters
.size();
j
++) {
34
daughters
[
j
].Configure(
i
);
35
}
36
}
37
38
void
BPhysicsValidation::analyze
(
const
edm::Event
&
iEvent
,
const
edm::EventSetup
& iSetup) {
39
edm::Handle<reco::GenParticleCollection>
genParticles
;
40
iEvent
.getByToken(
genparticleCollectionToken_
,
genParticles
);
41
for
(reco::GenParticleCollection::const_iterator iter =
genParticles
->begin(); iter !=
genParticles
->end(); ++iter) {
42
if
(
abs
(iter->pdgId()) ==
abs
(
particle
.
PDGID
())) {
43
Nobj
->
Fill
(0.5, 1.0);
44
particle
.
Fill
(&(*iter), 1.0);
45
FillDaughters
(&(*iter));
46
}
47
}
48
}
49
50
void
BPhysicsValidation::FillDaughters
(
const
reco::GenParticle
*
p
) {
51
int
mpdgid =
p
->pdgId();
52
for
(
unsigned
int
i
= 0;
i
<
p
->numberOfDaughters();
i
++) {
53
const
reco::GenParticle
* dau =
static_cast<
const
reco::GenParticle
*
>
(
p
->daughter(
i
));
54
int
pdgid
= dau->
pdgId
();
55
for
(
unsigned
int
i
= 0;
i
<
daughters
.size();
i
++) {
56
if
(
abs
(mpdgid) !=
abs
(
daughters
[
i
].PDGID()) &&
daughters
[
i
].PDGID() ==
pdgid
)
57
daughters
[
i
].
Fill
(dau, 1.0);
58
// note: use abs when comparing to mother to avoid mixing
59
}
60
FillDaughters
(dau);
61
}
62
}
BPhysicsValidation::ParticleMonitor::Configure
void Configure(DQMStore::IBooker &i)
Definition:
BPhysicsValidation.h:49
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition:
ParameterSet.h:303
mps_fire.i
i
Definition:
mps_fire.py:428
dqmiolumiharvest.j
j
Definition:
dqmiolumiharvest.py:66
BPhysicsValidation::Nobj
MonitorElement * Nobj
Definition:
BPhysicsValidation.h:87
BPhysicsValidation::BPhysicsValidation
BPhysicsValidation(const edm::ParameterSet &)
Definition:
BPhysicsValidation.cc:13
MakerMacros.h
edm::Handle< reco::GenParticleCollection >
BPhysicsValidation::~BPhysicsValidation
~BPhysicsValidation() override
Definition:
BPhysicsValidation.cc:26
BPhysicsValidation::genparticleCollectionToken_
edm::EDGetTokenT< reco::GenParticleCollection > genparticleCollectionToken_
Definition:
BPhysicsValidation.h:83
std
Definition:
JetResolutionObject.h:76
EgammaValidation_cff.pdgid
pdgid
Definition:
EgammaValidation_cff.py:29
AlCaHLTBitMon_QueryRunRegistry.string
string string
Definition:
AlCaHLTBitMon_QueryRunRegistry.py:256
edm::ParameterSet::getUntrackedParameter
T getUntrackedParameter(std::string const &, T const &) const
reco::LeafCandidate::pdgId
int pdgId() const final
PDG identifier.
Definition:
LeafCandidate.h:176
dqm::impl::MonitorElement::Fill
void Fill(long long x)
Definition:
MonitorElement.h:292
dqm::implementation::IBooker
Definition:
DQMStore.h:43
iEvent
int iEvent
Definition:
GenABIO.cc:224
BPhysicsValidation::particle
ParticleMonitor particle
Definition:
BPhysicsValidation.h:85
AJJGenJetFilter_cfi.genParticles
genParticles
Definition:
AJJGenJetFilter_cfi.py:5
HcalObjRepresent::Fill
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
Definition:
HcalObjRepresent.h:1053
BPhysicsValidation::bookHistograms
void bookHistograms(DQMStore::IBooker &i, edm::Run const &, edm::EventSetup const &) override
Definition:
BPhysicsValidation.cc:28
BPhysicsValidation::genparticleCollection_
edm::InputTag genparticleCollection_
Definition:
BPhysicsValidation.h:82
funct::abs
Abs< T >::type abs(const T &t)
Definition:
Abs.h:22
DQMHelper
Definition:
DQMHelper.h:15
edm::EventSetup
Definition:
EventSetup.h:59
BPhysicsValidation.h
BPhysicsValidation::FillDaughters
void FillDaughters(const reco::GenParticle *p)
Definition:
BPhysicsValidation.cc:50
BPhysicsValidation::ParticleMonitor
Definition:
BPhysicsValidation.h:43
edm
HLT enums.
Definition:
AlignableModifier.h:19
edm::InputTag
Definition:
InputTag.h:15
reco::GenParticle
Definition:
GenParticle.h:21
BPhysicsValidation::ParticleMonitor::Fill
void Fill(const reco::GenParticle *p, double weight)
Definition:
BPhysicsValidation.h:64
edm::ParameterSet
Definition:
ParameterSet.h:47
BPhysicsValidation::daughters
std::vector< ParticleMonitor > daughters
Definition:
BPhysicsValidation.h:86
BPhysicsValidation::name
std::string name
Definition:
BPhysicsValidation.h:84
edm::Event
Definition:
Event.h:73
dqm
Definition:
DQMStore.h:18
BPhysicsValidation::ParticleMonitor::PDGID
int PDGID()
Definition:
BPhysicsValidation.h:72
AlCaHLTBitMon_ParallelJobs.p
def p
Definition:
AlCaHLTBitMon_ParallelJobs.py:153
edm::Run
Definition:
Run.h:45
Skims_PA_cff.name
name
Definition:
Skims_PA_cff.py:17
BPhysicsValidation::analyze
void analyze(edm::Event const &, edm::EventSetup const &) override
Definition:
BPhysicsValidation.cc:38
Generated for CMSSW Reference Manual by
1.8.14