IOMC
ParticleGuns
src
FlatRandomEGunProducer.cc
Go to the documentation of this file.
1
/*
2
* \author Julia Yarba
3
*/
4
5
#include <ostream>
6
7
#include "
IOMC/ParticleGuns/interface/FlatRandomEGunProducer.h
"
8
9
#include "
SimDataFormats/GeneratorProducts/interface/HepMCProduct.h
"
10
#include "
SimDataFormats/GeneratorProducts/interface/GenEventInfoProduct.h
"
11
12
#include "
FWCore/Framework/interface/Event.h
"
13
#include "
FWCore/ParameterSet/interface/ParameterSet.h
"
14
#include "
FWCore/ServiceRegistry/interface/Service.h
"
15
#include "
FWCore/Utilities/interface/RandomNumberGenerator.h
"
16
17
#include "CLHEP/Random/RandFlat.h"
18
19
using namespace
edm
;
20
using namespace
std
;
21
22
FlatRandomEGunProducer::FlatRandomEGunProducer
(
const
ParameterSet
&
pset
) :
BaseFlatGunProducer
(
pset
) {
23
ParameterSet
defpset;
24
// ParameterSet pgun_params = pset.getParameter<ParameterSet>("PGunParameters") ;
25
ParameterSet
pgun_params =
pset
.getParameter<
ParameterSet
>(
"PGunParameters"
);
26
27
// doesn't seem necessary to check if pset is empty - if this
28
// is the case, default values will be taken for params
29
fMinE
= pgun_params.
getParameter
<
double
>(
"MinE"
);
30
fMaxE
= pgun_params.
getParameter
<
double
>(
"MaxE"
);
31
32
produces<HepMCProduct>(
"unsmeared"
);
33
produces<GenEventInfoProduct>();
34
35
cout
<<
"Internal FlatRandomEGun is initialzed"
<< endl;
36
// cout << "It is going to generate " << remainingEvents() << "events" << endl ;
37
}
38
39
FlatRandomEGunProducer::~FlatRandomEGunProducer
() {
40
// no need to cleanup fEvt since it's done in HepMCProduct
41
}
42
43
void
FlatRandomEGunProducer::produce
(
Event
&
e
,
const
EventSetup
& es) {
44
edm::Service<edm::RandomNumberGenerator>
rng;
45
CLHEP::HepRandomEngine* engine = &rng->
getEngine
(
e
.streamID());
46
47
if
(
fVerbosity
> 0) {
48
cout
<<
" FlatRandomEGunProducer : Begin New Event Generation"
<< endl;
49
}
50
51
// event loop (well, another step in it...)
52
53
// no need to clean up GenEvent memory - done in HepMCProduct
54
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 = 1;
69
for
(
unsigned
int
ip = 0; ip <
fPartIDs
.size(); ip++) {
70
double
energy
= CLHEP::RandFlat::shoot(engine,
fMinE
,
fMaxE
);
71
double
eta
= CLHEP::RandFlat::shoot(engine,
fMinEta
,
fMaxEta
);
72
double
phi = CLHEP::RandFlat::shoot(engine,
fMinPhi
,
fMaxPhi
);
73
int
PartID
=
fPartIDs
[ip];
74
const
HepPDT::ParticleData
* PData =
fPDGTable
->particle(
HepPDT::ParticleID
(
abs
(
PartID
)));
75
double
mass
= PData->mass().value();
76
double
mom2 =
energy
*
energy
-
mass
*
mass
;
77
double
mom = 0.;
78
if
(mom2 > 0.) {
79
mom =
sqrt
(mom2);
80
}
else
{
81
mom = 0.;
82
}
83
double
theta
= 2. * atan(
exp
(-
eta
));
84
double
px
= mom *
sin
(
theta
) *
cos
(phi);
85
double
py
= mom *
sin
(
theta
) *
sin
(phi);
86
double
pz = mom *
cos
(
theta
);
87
88
HepMC::FourVector
p
(
px
,
py
, pz,
energy
);
89
HepMC::GenParticle
* Part =
new
HepMC::GenParticle
(
p
,
PartID
, 1);
90
Part->suggest_barcode(barcode);
91
barcode++;
92
Vtx->add_particle_out(Part);
93
94
if
(
fAddAntiParticle
) {
95
HepMC::FourVector ap(-
px
, -
py
, -pz,
energy
);
96
int
APartID = -
PartID
;
97
if
(
PartID
== 22 ||
PartID
== 23) {
98
APartID =
PartID
;
99
}
100
HepMC::GenParticle
* APart =
new
HepMC::GenParticle
(ap, APartID, 1);
101
APart->suggest_barcode(barcode);
102
barcode++;
103
Vtx->add_particle_out(APart);
104
}
105
}
106
fEvt
->add_vertex(Vtx);
107
fEvt
->set_event_number(
e
.id().event());
108
fEvt
->set_signal_process_id(20);
109
110
if
(
fVerbosity
> 0) {
111
fEvt
->print();
112
}
113
114
unique_ptr<HepMCProduct> BProduct(
new
HepMCProduct
());
115
BProduct->addHepMCData(
fEvt
);
116
e
.put(
std::move
(BProduct),
"unsmeared"
);
117
118
unique_ptr<GenEventInfoProduct>
genEventInfo
(
new
GenEventInfoProduct
(
fEvt
));
119
e
.put(
std::move
(
genEventInfo
));
120
121
if
(
fVerbosity
> 0) {
122
// for testing purpose only
123
//fEvt->print() ; // for some strange reason, it prints NO event info
124
// after it's been put into edm::Event...
125
cout
<<
" FlatRandomEGunProducer : Event Generation Done "
<< endl;
126
}
127
}
edm::BaseFlatGunProducer::fMaxPhi
double fMaxPhi
Definition:
BaseFlatGunProducer.h:45
GenEventInfoProduct
Definition:
GenEventInfoProduct.h:17
edm::RandomNumberGenerator::getEngine
virtual CLHEP::HepRandomEngine & getEngine(StreamID const &)=0
Use this engine in event methods.
edm::FlatRandomEGunProducer::FlatRandomEGunProducer
FlatRandomEGunProducer(const ParameterSet &pset)
Definition:
FlatRandomEGunProducer.cc:22
edm::BaseFlatGunProducer::fAddAntiParticle
bool fAddAntiParticle
Definition:
BaseFlatGunProducer.h:61
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
gather_cfg.cout
cout
Definition:
gather_cfg.py:144
edm::BaseFlatGunProducer::fMinPhi
double fMinPhi
Definition:
BaseFlatGunProducer.h:44
edm::FlatRandomEGunProducer::fMaxE
double fMaxE
Definition:
FlatRandomEGunProducer.h:25
edm::FlatRandomEGunProducer::fMinE
double fMinE
Definition:
FlatRandomEGunProducer.h:24
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
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
edm::ParameterSet
Definition:
ParameterSet.h:47
GenEventInfoProduct.h
Event.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
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
std
Definition:
JetResolutionObject.h:76
genParticles2HepMC_cfi.genEventInfo
genEventInfo
Definition:
genParticles2HepMC_cfi.py:6
EgHLTOffHistBins_cfi.mass
mass
Definition:
EgHLTOffHistBins_cfi.py:34
edm::BaseFlatGunProducer::fVerbosity
int fVerbosity
Definition:
BaseFlatGunProducer.h:59
edm::FlatRandomEGunProducer::produce
void produce(Event &e, const EventSetup &es) override
Definition:
FlatRandomEGunProducer.cc:43
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition:
ParameterSet.h:303
FlatRandomEGunProducer.h
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::FlatRandomEGunProducer::~FlatRandomEGunProducer
~FlatRandomEGunProducer() override
Definition:
FlatRandomEGunProducer.cc:39
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