IOMC
ParticleGuns
src
FlatRandomMultiParticlePGunProducer.cc
Go to the documentation of this file.
1
#include <ostream>
2
3
#include "
IOMC/ParticleGuns/interface/FlatRandomMultiParticlePGunProducer.h
"
4
5
#include "
SimDataFormats/GeneratorProducts/interface/HepMCProduct.h
"
6
#include "
SimDataFormats/GeneratorProducts/interface/GenEventInfoProduct.h
"
7
8
#include "
FWCore/Framework/interface/Event.h
"
9
#include "
FWCore/MessageLogger/interface/MessageLogger.h
"
10
#include "
FWCore/ParameterSet/interface/ParameterSet.h
"
11
#include "
FWCore/ServiceRegistry/interface/Service.h
"
12
#include "
FWCore/Utilities/interface/RandomNumberGenerator.h
"
13
#include "
FWCore/Utilities/interface/EDMException.h
"
14
#include "CLHEP/Random/RandFlat.h"
15
16
using namespace
edm
;
17
18
FlatRandomMultiParticlePGunProducer::FlatRandomMultiParticlePGunProducer
(
const
ParameterSet
&
pset
)
19
:
BaseFlatGunProducer
(
pset
) {
20
ParameterSet
pgunParams =
pset
.getParameter<
ParameterSet
>(
"PGunParameters"
);
21
fProbParticle_
= pgunParams.
getParameter
<std::vector<double> >(
"ProbParts"
);
22
fMinP_
= pgunParams.
getParameter
<
double
>(
"MinP"
);
23
fMaxP_
= pgunParams.
getParameter
<
double
>(
"MaxP"
);
24
if
(
fProbParticle_
.size() !=
fPartIDs
.size())
25
throw
cms::Exception
(
"Configuration"
) <<
"Not all probabilities given for all particle types "
26
<<
fProbParticle_
.size() <<
":"
<<
fPartIDs
.size() <<
" need them to match\n"
;
27
28
produces<HepMCProduct>(
"unsmeared"
);
29
produces<GenEventInfoProduct>();
30
31
edm::LogVerbatim
(
"ParticleGun"
) <<
"FlatRandomMultiParticlePGun is initialzed for "
<<
fPartIDs
.size()
32
<<
" particles in momentum range "
<<
fMinP_
<<
":"
<<
fMaxP_
;
33
for
(
unsigned
int
k
= 0;
k
<
fPartIDs
.size(); ++
k
)
34
edm::LogVerbatim
(
"ParticleGun"
) <<
" ["
<<
k
<<
"] "
<<
fPartIDs
[
k
] <<
":"
<<
fProbParticle_
[
k
];
35
36
for
(
unsigned
int
k
= 1;
k
<
fProbParticle_
.size(); ++
k
)
37
fProbParticle_
[
k
] +=
fProbParticle_
[
k
- 1];
38
for
(
unsigned
int
k
= 0;
k
<
fProbParticle_
.size(); ++
k
) {
39
fProbParticle_
[
k
] /=
fProbParticle_
[
fProbParticle_
.size() - 1];
40
edm::LogVerbatim
(
"ParticleGun"
) <<
"Corrected probaility for "
<<
fPartIDs
[
k
] <<
":"
<<
fProbParticle_
[
k
];
41
}
42
}
43
44
FlatRandomMultiParticlePGunProducer::~FlatRandomMultiParticlePGunProducer
() {}
45
46
void
FlatRandomMultiParticlePGunProducer::produce
(
edm::Event
&
e
,
const
edm::EventSetup
& es) {
47
edm::Service<edm::RandomNumberGenerator>
rng;
48
CLHEP::HepRandomEngine* engine = &rng->
getEngine
(
e
.streamID());
49
50
if
(
fVerbosity
> 0)
51
edm::LogVerbatim
(
"ParticleGun"
) <<
"FlatRandomMultiParticlePGunProducer: Begin New Event Generation"
;
52
53
// event loop (well, another step in it...)
54
// no need to clean up GenEvent memory - done in HepMCProduct
55
// here re-create fEvt (memory)
56
//
57
fEvt
=
new
HepMC::GenEvent
();
58
59
// now actualy, cook up the event from PDGTable and gun parameters
60
//
61
62
// 1st, primary vertex
63
//
64
HepMC::GenVertex* Vtx =
new
HepMC::GenVertex(HepMC::FourVector(0., 0., 0.));
65
66
// loop over particles
67
//
68
int
barcode(0),
PartID
(
fPartIDs
[0]);
69
double
r1
= CLHEP::RandFlat::shoot(engine, 0., 1.);
70
for
(
unsigned
int
ip = 0; ip <
fPartIDs
.size(); ip++) {
71
if
(
r1
<=
fProbParticle_
[ip])
72
break
;
73
PartID
=
fPartIDs
[ip];
74
}
75
if
(
fVerbosity
> 0)
76
edm::LogVerbatim
(
"ParticleGun"
) <<
"Random "
<<
r1
<<
" PartID "
<<
PartID
;
77
78
double
mom = CLHEP::RandFlat::shoot(engine,
fMinP_
,
fMaxP_
);
79
double
eta
= CLHEP::RandFlat::shoot(engine,
fMinEta
,
fMaxEta
);
80
double
phi
= CLHEP::RandFlat::shoot(engine,
fMinPhi
,
fMaxPhi
);
81
const
HepPDT::ParticleData
* PData =
fPDGTable
->particle(
HepPDT::ParticleID
(
abs
(
PartID
)));
82
double
mass
= PData->mass().value();
83
double
energy
=
sqrt
(mom * mom +
mass
*
mass
);
84
double
theta
= 2. * atan(
exp
(-
eta
));
85
double
px
= mom *
sin
(
theta
) *
cos
(
phi
);
86
double
py
= mom *
sin
(
theta
) *
sin
(
phi
);
87
double
pz = mom *
cos
(
theta
);
88
89
HepMC::FourVector
p
(
px
,
py
, pz,
energy
);
90
HepMC::GenParticle
* Part =
new
HepMC::GenParticle
(
p
,
PartID
, 1);
91
barcode++;
92
Part->suggest_barcode(barcode);
93
Vtx->add_particle_out(Part);
94
95
if
(
fAddAntiParticle
) {
96
HepMC::FourVector ap(-
px
, -
py
, -pz,
energy
);
97
int
APartID = (
PartID
== 22 ||
PartID
== 23) ?
PartID
: -
PartID
;
98
HepMC::GenParticle
* APart =
new
HepMC::GenParticle
(ap, APartID, 1);
99
barcode++;
100
APart->suggest_barcode(barcode);
101
Vtx->add_particle_out(APart);
102
}
103
104
fEvt
->add_vertex(Vtx);
105
fEvt
->set_event_number(
e
.id().event());
106
fEvt
->set_signal_process_id(20);
107
108
if
(
fVerbosity
> 1)
109
fEvt
->print();
110
111
std::unique_ptr<HepMCProduct> BProduct(
new
HepMCProduct
());
112
BProduct->addHepMCData(
fEvt
);
113
e
.put(
std::move
(BProduct),
"unsmeared"
);
114
115
std::unique_ptr<GenEventInfoProduct>
genEventInfo
(
new
GenEventInfoProduct
(
fEvt
));
116
e
.put(
std::move
(
genEventInfo
));
117
118
if
(
fVerbosity
> 0)
119
edm::LogVerbatim
(
"ParticleGun"
) <<
" FlatRandomMultiParticlePGunProducer : Event Generation Done"
;
120
}
edm::BaseFlatGunProducer::fMaxPhi
double fMaxPhi
Definition:
BaseFlatGunProducer.h:45
GenEventInfoProduct
Definition:
GenEventInfoProduct.h:17
edm::FlatRandomMultiParticlePGunProducer::FlatRandomMultiParticlePGunProducer
FlatRandomMultiParticlePGunProducer(const ParameterSet &pset)
Definition:
FlatRandomMultiParticlePGunProducer.cc:18
edm::RandomNumberGenerator::getEngine
virtual CLHEP::HepRandomEngine & getEngine(StreamID const &)=0
Use this engine in event methods.
edm::FlatRandomMultiParticlePGunProducer::fProbParticle_
std::vector< double > fProbParticle_
Definition:
FlatRandomMultiParticlePGunProducer.h:18
MessageLogger.h
edm::BaseFlatGunProducer::fAddAntiParticle
bool fAddAntiParticle
Definition:
BaseFlatGunProducer.h:61
edm::FlatRandomMultiParticlePGunProducer::~FlatRandomMultiParticlePGunProducer
~FlatRandomMultiParticlePGunProducer() override
Definition:
FlatRandomMultiParticlePGunProducer.cc:44
edm::BaseFlatGunProducer::fMaxEta
double fMaxEta
Definition:
BaseFlatGunProducer.h:43
multPhiCorr_741_25nsDY_cfi.py
py
Definition:
multPhiCorr_741_25nsDY_cfi.py:12
edm
HLT enums.
Definition:
AlignableModifier.h:19
RandomNumberGenerator.h
AlCaHLTBitMon_ParallelJobs.p
p
Definition:
AlCaHLTBitMon_ParallelJobs.py:153
edm::BaseFlatGunProducer::fEvt
HepMC::GenEvent * fEvt
Definition:
BaseFlatGunProducer.h:48
ZgammaFilter_cfi.HepMCProduct
HepMCProduct
Definition:
ZgammaFilter_cfi.py:9
edm::BaseFlatGunProducer::fMinPhi
double fMinPhi
Definition:
BaseFlatGunProducer.h:44
edm::FlatRandomMultiParticlePGunProducer::fMaxP_
double fMaxP_
Definition:
FlatRandomMultiParticlePGunProducer.h:20
HepMC::GenEvent
Definition:
hepmc_rootio.cc:9
ParticleData
HepPDT::ParticleData ParticleData
Definition:
ParticleDataTable.h:9
funct::sin
Sin< T >::type sin(const T &t)
Definition:
Sin.h:22
EDMException.h
edm::BaseFlatGunProducer::fPartIDs
std::vector< int > fPartIDs
Definition:
BaseFlatGunProducer.h:41
funct::cos
Cos< T >::type cos(const T &t)
Definition:
Cos.h:22
Service.h
PVValHelper::eta
Definition:
PVValidationHelpers.h:70
mathSSE::sqrt
T sqrt(T t)
Definition:
SSEVec.h:19
theta
Geom::Theta< T > theta() const
Definition:
Basic3DVectorLD.h:150
HCALHighEnergyHPDFilter_cfi.energy
energy
Definition:
HCALHighEnergyHPDFilter_cfi.py:5
dqmdumpme.k
k
Definition:
dqmdumpme.py:60
edm::ParameterSet
Definition:
ParameterSet.h:47
GenEventInfoProduct.h
Event.h
FlatRandomMultiParticlePGunProducer.h
edm::Service< edm::RandomNumberGenerator >
edm::EventSetup
Definition:
EventSetup.h:58
genfragment_ptgun_cfg.PartID
PartID
Definition:
genfragment_ptgun_cfg.py:6
edm::BaseFlatGunProducer
Definition:
BaseFlatGunProducer.h:26
edm::BaseFlatGunProducer::fPDGTable
ESHandle< HepPDT::ParticleDataTable > fPDGTable
Definition:
BaseFlatGunProducer.h:57
DDAxes::phi
multPhiCorr_741_25nsDY_cfi.px
px
Definition:
multPhiCorr_741_25nsDY_cfi.py:10
GenParticle.GenParticle
GenParticle
Definition:
GenParticle.py:18
eostools.move
def move(src, dest)
Definition:
eostools.py:511
genParticles2HepMC_cfi.genEventInfo
genEventInfo
Definition:
genParticles2HepMC_cfi.py:6
diffTwoXMLs.r1
r1
Definition:
diffTwoXMLs.py:53
edm::LogVerbatim
Log< level::Info, true > LogVerbatim
Definition:
MessageLogger.h:128
edm::FlatRandomMultiParticlePGunProducer::fMinP_
double fMinP_
Definition:
FlatRandomMultiParticlePGunProducer.h:19
EgHLTOffHistBins_cfi.mass
mass
Definition:
EgHLTOffHistBins_cfi.py:34
edm::FlatRandomMultiParticlePGunProducer::produce
void produce(Event &e, const EventSetup &es) override
Definition:
FlatRandomMultiParticlePGunProducer.cc:46
edm::BaseFlatGunProducer::fVerbosity
int fVerbosity
Definition:
BaseFlatGunProducer.h:59
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition:
ParameterSet.h:303
cms::Exception
Definition:
Exception.h:70
funct::abs
Abs< T >::type abs(const T &t)
Definition:
Abs.h:22
ParameterSet.h
HepMCProduct.h
JetChargeProducer_cfi.exp
exp
Definition:
JetChargeProducer_cfi.py:6
edm::Event
Definition:
Event.h:73
edm::Log
Definition:
MessageLogger.h:70
LHEGenericFilter_cfi.ParticleID
ParticleID
Definition:
LHEGenericFilter_cfi.py:6
edm::BaseFlatGunProducer::fMinEta
double fMinEta
Definition:
BaseFlatGunProducer.h:42
muonDTDigis_cfi.pset
pset
Definition:
muonDTDigis_cfi.py:27
MillePedeFileConverter_cfg.e
e
Definition:
MillePedeFileConverter_cfg.py:37
Generated for CMSSW Reference Manual by
1.8.16