8 #ifndef Pythia8_LHAHelaconia_H 9 #define Pythia8_LHAHelaconia_H 11 #include "Pythia8/Pythia.h" 41 string exeIn =
"ho_cluster");
54 bool setSeed(
int seedIn,
int runsIn = 30081);
70 bool run(
int eventsIn,
int seedIn = -1);
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;
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);
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 ----------------------------------------------* " 175 size_t n = line.find(
"state");
176 if (line.find(
"8)") != string::npos) mQ = -1;
177 if (n != string::npos && pythia) {
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);
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;}
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;
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;
308 if (!pythia)
return false;
309 if (lhef)
delete lhef;
310 bool setScales(pythia->
settings.
flag(
"Beams:setProductionScalesFromLHEF"));
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());
330 addProcess(lhef->idProcess(0), sig, err, lhef->xMax(0));
331 xSecSumSave = sig; xErrSumSave = err;
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;
359 if (!pythia)
return false;
360 if (!
run(events))
return false;
361 if (!
reader(
true))
return false;
374 if (!pythia)
return false;
376 errorMsg(
"Error from LHAupHelaconia::setEvent: LHAupLHEF " 377 "object not correctly initialized");
return false;}
379 errorMsg(
"Error from LHAupHelaconia::setEvent: LHEF " 380 "event file was not found");
return false;}
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
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)
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());
void setEvents(int eventsIn)
Set the number of events to generate per run.
Definition: LHAHelaconia.h:232
void setBeamA(int idIn, double eIn, int pdfGroupIn=0, int pdfSetIn=0)
Input beam info.
Definition: LesHouches.h:222
void addProcess(int idProcIn, double xSecIn=1., double xErrIn=0., double xMaxIn=1.)
Input process info.
Definition: LesHouches.h:233
bool flag(string keyIn)
Give back current value, with check that key exists.
Definition: Settings.cc:1478
bool fileFound()
Confirm that file was found and opened as expected.
Definition: LesHouches.h:408
void setStrategy(int strategyIn)
Input process weight strategy.
Definition: LesHouches.h:230
ParticleData particleData
ParticleData: the particle data table/database.
Definition: Pythia.h:335
bool readString(string line, bool warn=true)
Read in one update from a single line.
Definition: Settings.cc:374
void setIdX(int id1In, int id2In, double x1In, double x2In)
Input info on flavour and x values of hard-process initiators.
Definition: LesHouches.h:263
bool setSeed(int seedIn, int runsIn=30081)
Set the random seed and maximum runs.
Definition: LHAHelaconia.h:209
A class for the particles stored in LHAup.
Definition: LesHouches.h:48
Settings settings
Settings: databases of flags/modes/parms/words to control run.
Definition: Pythia.h:332
int id1() const
Give back info on flavour and x values of hard-process initiators.
Definition: LesHouches.h:158
Definition: LesHouches.h:78
void addParticle(LHAParticle particleIn)
Input particle info, one particle at the time.
Definition: LesHouches.h:252
void setProcess(int idProcIn=0, double weightIn=1., double scaleIn=0., double alphaQEDIn=0.0073, double alphaQCDIn=0.12)
Input info on the selected process.
Definition: LesHouches.h:243
A derived class with information read from a Les Houches Event File.
Definition: LesHouches.h:346
void setPdf(int id1pdfIn, int id2pdfIn, double x1pdfIn, double x2pdfIn, double scalePDFIn, double pdf1In, double pdf2In, bool pdfIsSetIn)
Optionally input info on parton density values of event.
Definition: LesHouches.h:267
int convert(int idIn)
Convert the color octet HelacOnia ID to a Pythia 8 ID.
Definition: LHAHelaconia.h:341
A derived class from LHAup which generates events with HelacOnia.
Definition: LHAHelaconia.h:35
int sizePart() const
Give back info on separate particle.
Definition: LesHouches.h:141
void errorMsg(string messageIn)
Print a message the first few times. Insert in database.
Definition: LHAHelaconia.h:79
bool readString(string line)
Read a HelacOnia command string.
Definition: LHAHelaconia.h:173
LHAupHelaconia(Pythia *pythiaIn, string dirIn="helaconiarun", string exeIn="ho_cluster")
Constructor.
Definition: LHAHelaconia.h:112
bool setEvent(int=0)
Set the event information.
Definition: LHAHelaconia.h:371
bool execute(string line)
Execute a system command.
Definition: LHAHelaconia.h:238
~LHAupHelaconia()
Destructor.
Definition: LHAHelaconia.h:125
bool pdfIsSet() const
Optional: give back info on parton density values of event.
Definition: LesHouches.h:164
Header for classes to set beam momentum and interaction vertex spread.
Definition: Analysis.h:20
The Pythia class contains the top-level routines to generate an event.
Definition: Pythia.h:71
bool setEvent(int=0)
Routine for doing the job of reading and setting info on next event.
Definition: LesHouches.h:420
Info * infoPtr
Pointer to various information on the generation.
Definition: LesHouches.h:218
bool reader(bool init)
Create the LHEF reader.
Definition: LHAHelaconia.h:305
bool readString(string, bool warn=true, int subrun=SUBRUNDEFAULT)
Read in one update for a setting or particle data from a single line.
Definition: Pythia.cc:365
bool run(int eventsIn, int seedIn=-1)
Run HelacOnia.
Definition: LHAHelaconia.h:244
void listInit()
Print the initialization info; useful to check that setting it worked.
Definition: LesHouches.cc:33
bool setInit()
Set the initialization information.
Definition: LHAHelaconia.h:356