Validation
EventGenerator
plugins
WeightManager.cc
Go to the documentation of this file.
1
#include "
Validation/EventGenerator/interface/WeightManager.h
"
2
#include "
FWCore/Framework/interface/Event.h
"
3
#include "
FWCore/ParameterSet/interface/ParameterSet.h
"
4
5
#include "
DataFormats/Common/interface/Handle.h
"
6
#include "
SimDataFormats/GeneratorProducts/interface/HepMCProduct.h
"
7
#include "
SimDataFormats/GeneratorProducts/interface/GenEventInfoProduct.h
"
8
#include "
Validation/EventGenerator/interface/DQMHelper.h
"
9
10
using namespace
edm
;
11
12
WeightManager::WeightManager
(
const
ParameterSet
& iConfig,
edm::ConsumesCollector
iC)
13
: _useHepMC(iConfig.getParameter<
bool
>(
"UseWeightFromHepMC"
)) {
14
if
(
_useHepMC
) {
15
_hepmcCollection
= iConfig.
getParameter
<
InputTag
>(
"hepmcCollection"
);
16
hepmcCollectionToken_
= iC.
consumes
<
HepMCProduct
>(
_hepmcCollection
);
17
}
else
{
18
_genEventInfos
= iConfig.
getParameter
<std::vector<InputTag> >(
"genEventInfos"
);
19
for
(
unsigned
int
i
= 0;
i
<
_genEventInfos
.size();
i
++)
20
genEventInfosTokens_
.push_back(iC.
consumes
<std::vector<InputTag> >(
_genEventInfos
[
i
]));
21
}
22
}
23
24
double
WeightManager::weight
(
const
Event
&
iEvent
) {
25
if
(
_useHepMC
) {
26
edm::Handle<HepMCProduct>
evt;
27
iEvent
.getByToken(
hepmcCollectionToken_
, evt);
28
const
HepMC::GenEvent
* myGenEvent = evt->
GetEvent
();
29
30
double
weight
= 1.;
31
if
(!myGenEvent->weights().empty())
32
weight
= myGenEvent->weights()[0];
33
return
weight
;
34
}
else
{
35
double
weight
= 1.;
36
for
(
unsigned
int
i
= 0;
i
<
genEventInfosTokens_
.size(); ++
i
) {
37
edm::Handle<GenEventInfoProduct>
info
;
38
iEvent
.getByToken(
genEventInfosTokens_
[
i
],
info
);
39
weight
*=
info
->weight();
40
}
41
return
weight
;
42
}
43
}
Handle.h
electrons_cff.bool
bool
Definition:
electrons_cff.py:372
mps_fire.i
i
Definition:
mps_fire.py:355
edm
HLT enums.
Definition:
AlignableModifier.h:19
info
static const TGPicture * info(bool iBackgroundIsBlack)
Definition:
FWCollectionSummaryWidget.cc:152
WeightManager::WeightManager
WeightManager(const edm::ParameterSet &, edm::ConsumesCollector iC)
Definition:
WeightManager.cc:12
edm::Handle
Definition:
AssociativeIterator.h:50
HepMC::GenEvent
Definition:
hepmc_rootio.cc:9
DQMHelper.h
edm::ConsumesCollector::consumes
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
Definition:
ConsumesCollector.h:49
WeightManager.h
edm::ParameterSet
Definition:
ParameterSet.h:36
GenEventInfoProduct.h
Event.h
iEvent
int iEvent
Definition:
GenABIO.cc:224
edm::HepMCProduct::GetEvent
const HepMC::GenEvent * GetEvent() const
Definition:
HepMCProduct.h:34
WeightManager::_genEventInfos
std::vector< edm::InputTag > _genEventInfos
Definition:
WeightManager.h:29
WeightManager::weight
double weight(const edm::Event &)
Definition:
WeightManager.cc:24
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
WeightManager::_hepmcCollection
edm::InputTag _hepmcCollection
Definition:
WeightManager.h:30
WeightManager::_useHepMC
bool _useHepMC
Definition:
WeightManager.h:28
WeightManager::genEventInfosTokens_
std::vector< edm::EDGetTokenT< std::vector< edm::InputTag > > > genEventInfosTokens_
Definition:
WeightManager.h:33
edm::HepMCProduct
Definition:
HepMCProduct.h:18
ParameterSet.h
HepMCProduct.h
edm::Event
Definition:
Event.h:73
edm::InputTag
Definition:
InputTag.h:15
edm::ConsumesCollector
Definition:
ConsumesCollector.h:39
weight
Definition:
weight.py:1
WeightManager::hepmcCollectionToken_
edm::EDGetTokenT< edm::HepMCProduct > hepmcCollectionToken_
Definition:
WeightManager.h:32
Generated for CMSSW Reference Manual by
1.8.16