8 #ifndef Pythia8_LHAHelaconia_H
9 #define Pythia8_LHAHelaconia_H
11 #include "Pythia8/Pythia.h"
35 class LHAupHelaconia :
public LHAup {
40 LHAupHelaconia(Pythia *pythiaIn,
string dirIn =
"helaconiarun",
41 string exeIn =
"ho_cluster");
48 bool readString(
string line);
51 void setEvents(
int eventsIn);
54 bool setSeed(
int seedIn,
int runsIn = 30081);
60 bool setEvent(
int = 0);
67 bool execute(
string line);
70 bool run(
int eventsIn,
int seedIn = -1);
73 bool reader(
bool init);
76 int convert(
int idIn);
79 void errorMsg(
string messageIn) {
81 int times = messages[messageIn];
82 ++messages[messageIn];
84 if (times < TIMESTOPRINT) cout <<
" PYTHIA " << messageIn << endl;
94 int events, seed, runs, nRuns, nId, nQ, nR, nL, nJ;
95 string dir, exe, lhegz;
102 map<string, int> messages;
104 static const int TIMESTOPRINT = 1;
112 LHAupHelaconia::LHAupHelaconia(Pythia *pythiaIn,
string dirIn,
string exeIn) :
113 pythia(pythiaIn), lhef(0), events(10000), seed(-1), runs(30081),
114 nRuns(0), nId(443), nQ(4), nR(0), nL(0), nJ(3),
115 dir(dirIn), exe(exeIn), lhegz(dirIn +
"/events.lhe"), mQ(-2) {
116 mkdir(dir.c_str(), 0777);
117 if (pythia) pythia->readString(
"Beams:frameType = 5");
118 pythia->settings.addMode(
"Onia:state", -1,
false,
false, 0, 0);
125 LHAupHelaconia::~LHAupHelaconia() {
127 if (lhef)
delete lhef;
130 cout <<
"\n *------- LHAupHelaconia Error and Warning Messages Statistics"
131 <<
" --------------------------------------------------* \n"
134 <<
" | times message "
140 map<string, int>::iterator messageEntry = messages.begin();
141 if (messageEntry == messages.end())
142 cout <<
" | 0 no errors or warnings to report "
144 while (messageEntry != messages.end()) {
146 string temp = messageEntry->first;
147 int len = temp.length();
148 temp.insert( len, max(0, 102 - len),
' ');
149 cout <<
" | " << setw(6) << messageEntry->second <<
" "
157 <<
" *------- End LHAupHelaconia Error and Warning Messages "
158 <<
"Statistics ----------------------------------------------* "
173 bool LHAupHelaconia::readString(
string line) {
175 size_t n = line.find(
"state");
176 if (line.find(
"8)") != string::npos) mQ = -1;
177 if (n != string::npos && pythia) {
178 pythia->settings.readString(
"Onia:" + line.substr(n));
179 nId = abs(pythia->settings.mode(
"Onia:state"));
180 nQ = int(nId/1e2) % 10;
181 nR = int(nId/1e5) % 10;
182 nL = int(nId/1e4) % 10;
183 nJ = int(nId/1e0) % 10;
184 }
else lines.push_back(line);
209 bool LHAupHelaconia::setSeed(
int seedIn,
int runsIn) {
211 if (!pythia)
return false;
214 seed = pythia->settings.mode(
"Random:seed");
216 errorMsg(
"Error from LHAupHelaconia::setSeed: the given "
217 "Pythia seed is less than 1.");
return false;}
220 if (seed * runs > 30081 * 30081) {
221 errorMsg(
"Error from LHAupHelaconia::setSeed: the given seed "
222 "exceeds the HelacOnia limit.");
return false;}
232 void LHAupHelaconia::setEvents(
int eventsIn) {events = eventsIn;}
238 bool LHAupHelaconia::execute(
string line) {
return system(line.c_str()) != -1;}
244 bool LHAupHelaconia::run(
int eventsIn,
int seedIn) {
247 if (!pythia)
return false;
249 errorMsg(
"Error from LHAupHelaconia::run: maximum number "
250 "of allowed runs exceeded.");
return false;}
251 if (seed < 0 && !setSeed(seed, runs))
return false;
252 if (seedIn < 0) seedIn = (seed - 1) * runs + nRuns + 1;
256 mQ = (pythia->particleData.m0(nId)
257 + pythia->settings.parm(
"Onia:massSplit"))/2.0;
260 if (!pythia)
return false;
261 fstream config((dir +
"/generate.py").c_str(), ios::out);
262 for (
int iLine = 0; iLine < (int)lines.size(); ++iLine)
263 config << lines[iLine] <<
"\n";
264 config <<
"set seed = " << seedIn <<
"\n"
266 <<
"set unwevt = " << eventsIn <<
"\n"
267 <<
"set preunw = " << 3.0/2.0*eventsIn <<
"\n";
268 if (mQ > 0) config <<
"set " << (nQ == 4 ?
"c" :
"b") <<
"mass = " << mQ
270 config <<
"launch\n";
274 fstream shuffle((dir +
"/shuffle.py").c_str(), ios::out);
276 "import random, os\n"
277 "random.seed(" << seedIn <<
"); tag, pre, post, events = '', [], [], []\n"
278 "for line in open('events.lhe').readlines():\n"
279 " if line.strip().startswith('<'):\n"
280 " tag = line.strip()\n"
281 " if tag == '<event>': events += ['<event>\\n']; continue\n"
282 " if tag == '</event>': events[-1] += '</event>\\n'; continue\n"
283 " if tag == '<event>': events[-1] += line\n"
284 " elif len(events) == 0: pre += [line]\n"
285 " else: post += [line]\n"
286 "random.shuffle(events); os.unlink('events.lhe')\n"
287 "open('events.lhe', 'w').writelines(pre + events + post)\n";
291 if (!execute(
"rm -rf " + dir +
"/PROC* " + lhegz))
return false;
292 if (!execute(
"cd " + dir +
"; cat generate.py | " + exe))
return false;
293 if (!execute(
"cd " + dir +
"; ln -s PROC_HO_0/P0_calc_0/output/*.lhe "
294 "events.lhe;# python shuffle.py"))
return false;
295 if (access(lhegz.c_str(), F_OK) == -1)
return false;
305 bool LHAupHelaconia::reader(
bool init) {
308 if (!pythia)
return false;
309 if (lhef)
delete lhef;
310 bool setScales(pythia->settings.flag(
"Beams:setProductionScalesFromLHEF"));
311 lhef =
new LHAupLHEF(infoPtr, lhegz.c_str(), NULL,
false, setScales);
312 if (!lhef->setInit()) {
313 errorMsg(
"Error from LHAupHelaconia::reader: failed to "
314 "initialize the LHEF reader");
return false;}
315 if (lhef->sizeProc() != 1) {
316 errorMsg(
"Error from LHAupHelaconia::reader: number of "
317 "processes is not 1");
return false;}
322 double sig(lhef->xSec(0)), err(lhef->xErr(0));
325 setBeamA(lhef->idBeamA(), lhef->eBeamA(), lhef->pdfGroupBeamA(),
326 lhef->pdfSetBeamA());
327 setBeamB(lhef->idBeamB(), lhef->eBeamB(), lhef->pdfGroupBeamB(),
328 lhef->pdfSetBeamB());
329 setStrategy(lhef->strategy());
330 addProcess(lhef->idProcess(0), sig, err, lhef->xMax(0));
331 xSecSumSave = sig; xErrSumSave = err;
341 int LHAupHelaconia::convert(
int idIn) {
343 if (abs(idIn) < 9900000)
return idIn;
344 else idIn = abs(idIn) - 9900000;
346 if (idIn == nQ*110 + 3) nS = 0;
347 else if (idIn == nQ*110 + 1) nS = 1;
348 return 9900000 + 10000*nQ + 1000*nS + 100*nR + 10*nL + nJ;
356 bool LHAupHelaconia::setInit() {
359 if (!pythia)
return false;
360 if (!run(events))
return false;
361 if (!reader(
true))
return false;
371 bool LHAupHelaconia::setEvent(
int) {
374 if (!pythia)
return false;
376 errorMsg(
"Error from LHAupHelaconia::setEvent: LHAupLHEF "
377 "object not correctly initialized");
return false;}
378 if (!lhef->fileFound()) {
379 errorMsg(
"Error from LHAupHelaconia::setEvent: LHEF "
380 "event file was not found");
return false;}
381 if (!lhef->setEvent()) {
382 if (!run(events))
return false;
383 if (!reader(
false))
return false;
388 particlesSave.clear();
390 for (
int ip = 1; ip < lhef->sizePart(); ++ip) {
391 mom1 = lhef->mother1(ip);
392 mom2 = lhef->mother2(ip);
393 particlesSave.push_back
394 (LHAParticle(convert(lhef->id(ip)),
395 lhef->status(ip), mom1, mom2, lhef->col1(ip),
396 lhef->col2(ip), lhef->px(ip), lhef->py(ip), lhef->pz(ip),
397 lhef->e(ip), lhef->m(ip), lhef->tau(ip), lhef->spin(ip),
399 if (mom1 > 0 && mom1 < (
int)particlesSave.size() && mom2 == 0)
400 particlesSave[mom1 - 1].statusPart = 2;
404 setProcess(lhef->idProcess(), lhef->weight(), lhef->scale(),
405 lhef->alphaQED(), lhef->alphaQCD());
406 for (
int ip = 0; ip < (int)particlesSave.size(); ++ip)
407 addParticle(particlesSave[ip]);
408 setIdX(lhef->id1(), lhef->id2(), lhef->x1(), lhef->x2());
409 setPdf(lhef->id1pdf(), lhef->id2pdf(), lhef->x1pdf(), lhef->x2pdf(),
410 lhef->scalePDF(), lhef->pdf1(), lhef->pdf2(), lhef->pdfIsSet());
419 #endif // Pythia8_LHAHelaconia_H