1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
|
#include <iostream>
#include <cmath>
#include <string>
#include "GeneratorInterface/ReggeGribovPartonMCInterface/interface/ReggeGribovPartonMCHadronizer.h"
#include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
#include "SimDataFormats/GeneratorProducts/interface/GenEventInfoProduct.h"
#include "SimDataFormats/GeneratorProducts/interface/GenRunInfoProduct.h"
#include "CLHEP/Random/RandomEngine.h"
#include <HepMC/GenCrossSection.h>
#include <HepMC/GenEvent.h>
#include <HepMC/GenVertex.h>
#include <HepMC/GenParticle.h>
#include <HepMC/HeavyIon.h>
#include <HepMC/PdfInfo.h>
#include <HepMC/Units.h>
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/Run.h"
#include "FWCore/MessageLogger/interface/MessageLogger.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "FWCore/Framework/interface/Frameworkfwd.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/MakerMacros.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
using namespace edm;
using namespace std;
using namespace gen;
#include "GeneratorInterface/ReggeGribovPartonMCInterface/interface/IO_EPOS.h"
#include "GeneratorInterface/ReggeGribovPartonMCInterface/interface/EPOS_Wrapper.h"
EPOS::IO_EPOS conv;
static CLHEP::HepRandomEngine* reggeGribovRandomEngine;
extern "C" {
float gen::rangen_() {
float a = float(reggeGribovRandomEngine->flat());
return a;
}
double gen::drangen_(int* idummy) {
double a = reggeGribovRandomEngine->flat();
return a;
}
}
ReggeGribovPartonMCHadronizer::ReggeGribovPartonMCHadronizer(const ParameterSet& pset)
: BaseHadronizer(pset),
pset_(pset),
m_BeamMomentum(pset.getParameter<double>("beammomentum")),
m_TargetMomentum(pset.getParameter<double>("targetmomentum")),
m_BeamID(pset.getParameter<int>("beamid")),
m_TargetID(pset.getParameter<int>("targetid")),
m_HEModel(pset.getParameter<int>("model")),
m_bMin(pset.getParameter<double>("bmin")),
m_bMax(pset.getParameter<double>("bmax")),
m_ParamFileName(pset.getUntrackedParameter<string>("paramFileName")),
m_SkipNuclFrag(pset.getParameter<bool>("skipNuclFrag")),
m_NEvent(0),
m_NParticles(0),
m_ImpactParameter(0.),
m_IsInitialized(false) {
int nevet = 1; //needed for CS
int noTables = 0; //don't calculate tables
int LHEoutput = 0; //no lhe
int dummySeed = 123;
char dummyName[] = "dummy";
crmc_set_f_(nevet,
dummySeed,
m_BeamMomentum,
m_TargetMomentum,
m_BeamID,
m_TargetID,
m_HEModel,
noTables,
LHEoutput,
dummyName,
m_ParamFileName.fullPath().c_str());
//additionally initialise tables
initializeTablePaths();
//change impact paramter
nucl2_.bminim = float(m_bMin);
nucl2_.bmaxim = float(m_bMax);
}
//_____________________________________________________________________
ReggeGribovPartonMCHadronizer::~ReggeGribovPartonMCHadronizer() {
// destructor
}
//_____________________________________________________________________
void ReggeGribovPartonMCHadronizer::doSetRandomEngine(CLHEP::HepRandomEngine* v) { reggeGribovRandomEngine = v; }
//_____________________________________________________________________
bool ReggeGribovPartonMCHadronizer::generatePartonsAndHadronize() {
int iout = 0, ievent = 0;
crmc_f_(iout,
ievent,
m_NParticles,
m_ImpactParameter,
m_PartID[0],
m_PartPx[0],
m_PartPy[0],
m_PartPz[0],
m_PartEnergy[0],
m_PartMass[0],
m_PartStatus[0]);
LogDebug("ReggeGribovPartonMCInterface") << "event generated" << endl;
const bool isHeavyIon = (m_TargetID + m_BeamID > 2);
if (isHeavyIon)
conv.set_trust_beam_particles(false);
conv.set_skip_nuclear_fragments(m_SkipNuclFrag);
HepMC::GenEvent* evt = conv.read_next_event();
evt->set_event_number(m_NEvent++);
int sig_id = -1;
switch (int(c2evt_.typevt)) // if negative typevt mini plasma was created by event (except -4)
{
case 0:
break; //unknown for qgsjetII
case 1:
sig_id = 101;
break;
case -1:
sig_id = 101;
break;
case 2:
sig_id = 105;
break;
case -2:
sig_id = 105;
break;
case 3:
sig_id = 102;
break;
case -3:
sig_id = 102;
break;
case 4:
sig_id = 103;
break;
case -4:
sig_id = 104;
break;
default:
LogDebug("ReggeGribovPartonMCInterface") << "Signal ID not recognised for setting HEPEVT" << endl;
}
evt->set_signal_process_id(sig_id); //an integer ID uniquely specifying the signal process (i.e. MSUB in Pythia)
#ifdef HEPMC_HAS_CROSS_SECTION
// set cross section information for this event
HepMC::GenCrossSection theCrossSection;
theCrossSection.set_cross_section(double(hadr5_.sigineaa) * 1e9); //required in pB
evt->set_cross_section(theCrossSection);
#endif
if (isHeavyIon) //other than pp
{
HepMC::HeavyIon ion(cevt_.kohevt, // Number of hard scatterings
cevt_.npjevt, // Number of projectile participants
cevt_.ntgevt, // Number of target participants
cevt_.kolevt, // Number of NN (nucleon-nucleon) collisions
cevt_.npnevt + cevt_.ntnevt, // Number of spectator neutrons
cevt_.nppevt + cevt_.ntpevt, // Number of spectator protons
-1, // Number of N-Nwounded collisions
-1, // Number of Nwounded-N collisons
-1, // Number of Nwounded-Nwounded collisions
cevt_.bimevt, // Impact Parameter(fm) of collision
cevt_.phievt, // Azimuthal angle of event plane
c2evt_.fglevt, // eccentricity of participating nucleons
hadr5_.sigine * 1e9); // nucleon-nucleon inelastic (in pB)
evt->set_heavy_ion(ion);
}
event().reset(evt);
//evt->print();
//EPOS::EPOS_Wrapper::print_hepcom();
return true;
}
//_____________________________________________________________________
bool ReggeGribovPartonMCHadronizer::hadronize() { return false; }
bool ReggeGribovPartonMCHadronizer::decay() { return true; }
bool ReggeGribovPartonMCHadronizer::residualDecay() { return true; }
void ReggeGribovPartonMCHadronizer::finalizeEvent() { return; }
void ReggeGribovPartonMCHadronizer::statistics() { return; }
const char* ReggeGribovPartonMCHadronizer::classname() const { return "gen::ReggeGribovPartonMCHadronizer"; }
bool ReggeGribovPartonMCHadronizer::declareStableParticles(const std::vector<int>&) { return true; }
bool ReggeGribovPartonMCHadronizer::initializeForInternalPartons() {
if (!m_IsInitialized) {
//use set parameters to init models
crmc_init_f_();
m_IsInitialized = true;
}
return true;
}
bool ReggeGribovPartonMCHadronizer::initializeTablePaths() {
//epos
string path_fnii(FileInPath("GeneratorInterface/ReggeGribovPartonMCInterface/data/epos.initl").fullPath());
string path_fnie(FileInPath("GeneratorInterface/ReggeGribovPartonMCInterface/data/epos.iniev").fullPath());
string path_fnrj(FileInPath("GeneratorInterface/ReggeGribovPartonMCInterface/data/epos.inirj").fullPath());
string path_fncs(FileInPath("GeneratorInterface/ReggeGribovPartonMCInterface/data/epos.inics").fullPath());
if (path_fnii.length() >= 500)
LogError("ReggeGribovPartonMCInterface") << "table path too long" << endl;
else
nfname_.nfnii = path_fnii.length();
if (path_fnie.length() >= 500)
LogError("ReggeGribovPartonMCInterface") << "table path too long" << endl;
else
nfname_.nfnie = path_fnie.length();
if (path_fnrj.length() >= 500)
LogError("ReggeGribovPartonMCInterface") << "table path too long" << endl;
else
nfname_.nfnrj = path_fnrj.length();
if (path_fncs.length() >= 500)
LogError("ReggeGribovPartonMCInterface") << "table path too long" << endl;
else
nfname_.nfncs = path_fncs.length();
strcpy(fname_.fnii, path_fnii.c_str());
strcpy(fname_.fnie, path_fnie.c_str());
strcpy(fname_.fnrj, path_fnrj.c_str());
strcpy(fname_.fncs, path_fncs.c_str());
//qgsjet
string path_fndat(FileInPath("GeneratorInterface/ReggeGribovPartonMCInterface/data/qgsjet.dat").fullPath());
string path_fnncs(FileInPath("GeneratorInterface/ReggeGribovPartonMCInterface/data/qgsjet.ncs").fullPath());
if (path_fndat.length() >= 500)
LogError("ReggeGribovPartonMCInterface") << "table path too long" << endl;
else
qgsnfname_.nfndat = path_fndat.length();
if (path_fnncs.length() >= 500)
LogError("ReggeGribovPartonMCInterface") << "table path too long" << endl;
else
qgsnfname_.nfnncs = path_fnncs.length();
strcpy(qgsfname_.fndat, path_fndat.c_str());
strcpy(qgsfname_.fnncs, path_fnncs.c_str());
qgsfname_.ifdat = 1; //option to overwrite the normal path
qgsfname_.ifncs = 2;
//qgsjetII
string path_fniidat(FileInPath("GeneratorInterface/ReggeGribovPartonMCInterface/data/qgsdat-II-04.lzma").fullPath());
string path_fniincs(FileInPath("GeneratorInterface/ReggeGribovPartonMCInterface/data/sectnu-II-04").fullPath());
if (path_fniidat.length() >= 500)
LogError("ReggeGribovPartonMCInterface") << "table path too long" << endl;
else
qgsiinfname_.nfniidat = path_fniidat.length();
if (path_fniincs.length() >= 500)
LogError("ReggeGribovPartonMCInterface") << "table path too long" << endl;
else
qgsiinfname_.nfniincs = path_fniincs.length();
strcpy(qgsiifname_.fniidat, path_fniidat.c_str());
strcpy(qgsiifname_.fniincs, path_fniincs.c_str());
qgsiifname_.ifiidat = 1; //option to overwrite the normal path
qgsiifname_.ifiincs = 2;
return true;
}
|