8 #ifndef Pythia8_LHAMadgraph_H
9 #define Pythia8_LHAMadgraph_H
11 #include "Pythia8/Pythia.h"
12 #include "Pythia8Plugins/JetMatching.h"
13 #include "Pythia8Plugins/GeneratorInput.h"
70 enum Stage{Auto, Configure, Generate, Launch};
74 string dirIn =
"madgraphrun",
string exeIn =
"mg5_aMC");
80 bool readString(
string line, Stage stage = Auto);
83 void addCard(
string src,
string dst);
86 void setEvents(
int eventsIn);
89 bool setSeed(
int seedIn,
int runsIn = 30081);
92 void setJets(
int jetsIn);
98 bool setEvent(
int = 0);
103 bool execute(
string line);
115 bool run(
int eventsIn,
int seedIn = -1);
118 bool reader(
bool init);
126 int events, seed, runs, nRuns, jets;
127 bool match, amcatnlo;
128 string dir, exe, lhegz;
130 vector< pair<string, string> > cards;
133 vector<string> configureLines, generateLines, launchLines;
136 vector<bool>
override;
144 LHAupMadgraph::LHAupMadgraph(
Pythia *pythiaIn,
bool matchIn,
string dirIn,
146 pythia(pythiaIn), lhef(0), hook(0), events(10000), seed(-1), runs(30081),
147 nRuns(0), jets(-1), match(matchIn), dir(dirIn), exe(exeIn),
148 lhegz(dirIn +
"/events.lhe.gz"), override(vector<bool>(3, false)) {
149 mkdir(dir.c_str(), 0777);
150 if (pythia) pythia->readString(
"Beams:frameType = 5");
157 LHAupMadgraph::~LHAupMadgraph() {
if (lhef)
delete lhef;
if (hook)
delete hook;}
180 bool LHAupMadgraph::readString(
string line, Stage stage) {
182 if (line.substr(0, 4) ==
" set") launchLines.push_back(line);
183 else if (line.substr(0, 10) ==
"configure ")
184 configureLines.push_back(line.substr(10));
185 else if (line.substr(0, 6) !=
"output" && line.substr(0, 6) !=
"launch")
186 generateLines.push_back(line);
188 }
else if (stage == Configure) {
189 override[Configure] =
true;
if (line !=
"") configureLines.push_back(line);
190 }
else if (stage == Generate) {
191 override[Generate] =
true; generateLines.push_back(line);
192 }
else if (stage == Launch) {
193 override[Launch] =
true; launchLines.push_back(line);
209 void LHAupMadgraph::addCard(
string src,
string dst) {
210 cards.push_back(make_pair(src, dst));
233 bool LHAupMadgraph::setSeed(
int seedIn,
int runsIn) {
235 if (!pythia)
return false;
238 seed = pythia->settings.mode(
"Random:seed");
240 pythia->info.errorMsg(
"Error from LHAupMadgraph::setSeed: the given "
241 "Pythia seed is less than 1.");
return false;}
244 if (seed * runs > 30081 * 30081) {
245 pythia->info.errorMsg(
"Error from LHAupMadgraph::setSeed: the given seed "
246 "exceeds the MadGraph limit.");
return false;}
256 void LHAupMadgraph::setEvents(
int eventsIn) {events = eventsIn;}
266 void LHAupMadgraph::setJets(
int jetsIn) {jets = jetsIn;}
272 bool LHAupMadgraph::execute(
string line) {
return system(line.c_str()) != -1;}
282 bool LHAupMadgraph::configure() {
284 if (
override[Configure] && configureLines.size() == 0)
return true;
285 mkdir((dir +
"/.mg5").c_str(), 0777);
286 fstream config((dir +
"/.mg5/mg5_configuration.txt").c_str(), ios::out);
287 for (
int iLine = 0; iLine < (int)configureLines.size(); ++iLine)
288 config << configureLines[iLine] <<
"\n";
289 if (!
override[Configure])
290 config <<
"automatic_html_opening = False\n"
291 <<
"auto_update = 0\n";
306 bool LHAupMadgraph::generate() {
309 if (!pythia)
return false;
310 fstream config((dir +
"/generate.py").c_str(), ios::out);
311 for (
int iLine = 0; iLine < (int)generateLines.size(); ++iLine)
312 config << generateLines[iLine] <<
"\n";
313 if (!
override[Generate]) config <<
"output " << dir <<
"/tmp -f -nojpeg\n";
317 fstream orig((dir +
"/.mg5/mg5_configuration.txt").c_str(), ios::in);
318 char *home = getenv(
"HOME");
319 setenv(
"HOME", dir.c_str(), 1);
320 bool success = execute(exe +
" " + dir +
"/generate.py");
321 setenv(
"HOME", home, 1);
322 if (!success) {orig.close();
return false;}
323 else if (access((dir +
"/tmp/Cards/run_card.dat").c_str(), F_OK) == -1) {
324 pythia->info.errorMsg(
"Error from LHAupMadgraph::generate: MadGraph "
325 "failed to produce run_card.dat");
326 orig.close();
return false;
327 }
else execute(
"mv " + dir +
"/tmp/* " + dir +
"; rmdir " + dir +
"/tmp");
331 access((dir +
"/Cards/amcatnlo_configuration.txt").c_str(), F_OK) != -1;
333 fstream copy((dir +
"/Cards/" + (amcatnlo ?
"amcatnlo" :
"me5") +
334 "_configuration.txt").c_str(), ios::out);
335 copy << orig.rdbuf(); copy.close();
340 for (
int iCard = 0; iCard < (int)cards.size(); ++iCard) {
341 ifstream src((cards[iCard].first).c_str(), ios::binary);
342 ofstream dst((dir +
"/Cards/" + cards[iCard].second).c_str(), ios::binary);
356 bool LHAupMadgraph::launch() {
359 if (!pythia)
return false;
360 fstream config((dir +
"/launch.py").c_str(), ios::out);
361 if (!
override[Launch]) {
362 config <<
"launch " << dir <<
" -n run";
363 if (amcatnlo) config <<
" -p\n" <<
"set parton_shower PYTHIA8\n"
364 <<
"set ickkw 3\n" <<
"set nevents 0\n"
365 <<
"set req_acc 0.001\n";
366 else config <<
" -s parton\n" <<
"set ickkw 1\n" <<
"set gridpack True\n";
370 for (
int iLine = 0; iLine < (int)launchLines.size(); ++iLine)
371 config << launchLines[iLine] <<
"\n";
372 if (!
override[Launch]) config <<
"done\n";
377 string line =
"cd " + dir +
"/MCatNLO/lib; LINKS=`ls`; for LINK in $LINKS;"
378 " do TARG=`readlink $LINK`; if [[ $TARG = ../* ]]; then "
379 "rm $LINK; ln -s ${TARG:3} $LINK; fi; done";
380 if (!execute(line)) {
381 pythia->info.errorMsg(
"Error from LHAupMadgraph::launch: failed to "
382 "link aMC@NLO libraries");
return false;}
386 if (!execute(exe +
" " + dir +
"/launch.py"))
return false;
388 if (access((dir +
"/SubProcesses/results.dat").c_str(), F_OK) == -1) {
389 pythia->info.errorMsg(
"Error from LHAupMadgraph::launch: aMC@NLO failed "
390 "to produce results.dat");
return false;}
391 fstream script((dir +
"/run.sh").c_str(), ios::out);
392 script <<
"#!/usr/bin/env bash\n"
393 <<
"sed -i \"s/.*= *nevents/$1 = nevents/g\" ./Cards/run_card.dat\n"
394 <<
"sed -i \"s/.*= *iseed/$2 = iseed/g\" ./Cards/run_card.dat\n"
395 <<
"./bin/generate_events --parton --nocompile --only_generation "
396 "--force --name run\n" <<
"mv Events/run/events.lhe.gz ./\n";
397 script.close(); execute(
"chmod 755 " + dir +
"/run.sh");
399 string gpk =
"run_gridpack.tar.gz";
400 if (access((dir +
"/" + gpk).c_str(), F_OK) == -1) {
401 pythia->info.errorMsg(
"Error from LHAupMadgraph::launch: MadEvent failed"
402 " to produce " + gpk);
return false;}
403 string line =
"cd " + dir +
"; tar -xzf " + gpk +
"; cd madevent/lib; "
404 "LINK=`readlink libLHAPDF.a`; if [[ $LINK = ../* ]]; then "
405 "rm libLHAPDF.a; ln -s ../$LINK libLHAPDF.a; fi; cd ../; "
406 "./bin/compile dynamic; ./bin/clean4grid";
407 if (!execute(line)) {
408 pythia->info.errorMsg(
"Error from LHAupMadgraph::launch: failed to "
409 "compile MadEvent code");
return false;}
419 bool LHAupMadgraph::run(
int eventsIn,
int seedIn) {
421 if (!pythia)
return false;
423 pythia->info.errorMsg(
"Error from LHAupMadgraph::run: maximum number "
424 "of allowed runs exceeded.");
return false;}
425 if (access((dir +
"/run.sh").c_str(), F_OK) == -1)
return false;
426 if (seed < 0 && !setSeed(seed, runs))
return false;
427 if (seedIn < 0) seedIn = (seed - 1) * runs + nRuns + 1;
429 line <<
"cd " + dir +
"; ./run.sh " << eventsIn <<
" " << seedIn;
430 if (!amcatnlo) line <<
"; rm -rf ./madevent/Events/*";
431 if (!execute(line.str()))
return false;
432 if (access(lhegz.c_str(), F_OK) == -1)
return false;
442 bool LHAupMadgraph::reader(
bool init) {
445 if (!pythia)
return false;
446 if (lhef)
delete lhef;
447 bool setScales(pythia->settings.flag(
"Beams:setProductionScalesFromLHEF"));
448 lhef =
new LHAupLHEF(&pythia->info, lhegz.c_str(), NULL,
false, setScales);
449 if (!lhef->setInit()) {
450 pythia->info.errorMsg(
"Error from LHAupMadgraph::reader: failed to "
451 "initialize the LHEF reader");
return false;}
452 if (lhef->sizeProc() != 1) {
453 pythia->info.errorMsg(
"Error from LHAupMadgraph::reader: number of "
454 "processes is not 1");
return false;}
459 double sig(lhef->xSec(0)), err(lhef->xErr(0));
461 fstream results((dir +
"/madevent/SubProcesses/"
462 "run_results.dat").c_str(), ios::in);
463 string v; vector<double> vs;
464 while (std::getline(results, v,
' ')) vs.push_back(atof(v.c_str()));
466 pythia->info.errorMsg(
"Error from LHAupMadgraph::reader: could not "
467 "extract cross-section");
return false;}
468 sig = vs[0]; err = vs[1];
472 setBeamA(lhef->idBeamA(), lhef->eBeamA(), lhef->pdfGroupBeamA(),
473 lhef->pdfSetBeamA());
474 setBeamB(lhef->idBeamB(), lhef->eBeamB(), lhef->pdfGroupBeamB(),
475 lhef->pdfSetBeamB());
476 setStrategy(lhef->strategy());
477 addProcess(lhef->idProcess(0), sig, err, lhef->xMax(0));
478 xSecSumSave = sig; xErrSumSave = err;
491 bool LHAupMadgraph::setInit() {
494 if (!pythia)
return false;
495 if (access((dir +
"/run.sh").c_str(), F_OK) == -1) {
497 pythia->info.errorMsg(
"Error from LHAupMadgraph::setInit: failed to "
498 "create the MadGraph configuration");
return false;}
500 pythia->info.errorMsg(
"Error from LHAupMadgraph::setInit: failed to "
501 "generate the MadGraph process");
return false;}
503 pythia->info.errorMsg(
"Error from LHAupMadgraph::setInit: failed to "
504 "launch the MadGraph process");
return false;}
507 access((dir +
"/Cards/amcatnlo_configuration.txt").c_str(), F_OK) != -1;
513 ifstream card((dir +
"/Cards/run_card.dat").c_str());
514 string str((istreambuf_iterator<char>(card)), istreambuf_iterator<char>());
520 int iLine = jets < 0 ? 0 : generateLines.size();
521 for (; iLine < (int)generateLines.size(); ++iLine) {
522 string line = generateLines[iLine];
523 size_t found = line.find(
">");
524 if (found == string::npos)
continue;
525 else line = line.substr(found);
526 stringstream sline(line);
string p;
int n(0);
527 while (std::getline(sline, p,
' ') && p !=
",")
528 if (p ==
"j" || p ==
"g" || p ==
"u" || p ==
"d" || p ==
"c" ||
529 p ==
"s" || p ==
"b" || p ==
"u~" || p ==
"d~" || p ==
"c~" ||
530 p ==
"s~" || p ==
"b~") ++n;
531 if (n > jets) jets = n;
535 double etaj = mad.getParam(
"etaj");
536 Settings &set = pythia->settings;
537 set.flag(
"JetMatching:merge",
true);
538 set.mode(
"JetMatching:scheme", 1);
539 set.flag(
"JetMatching:setMad",
false);
540 set.mode(
"JetMatching:nQmatch", mad.getParamAsInt(
"maxjetflavor"));
541 set.parm(
"JetMatching:qCut", mad.getParam(
"ptj"));
542 set.parm(
"JetMatching:etaJetMax", etaj > 0 ? etaj : 100);
543 set.mode(
"JetMatching:nJetMax", jets);
544 set.parm(
"Check:epTolErr", 1e-2);
548 set.parm(
"JetMatching:coneRadius", mad.getParam(
"jetradius"));
549 set.mode(
"JetMatching:slowJetPower", mad.getParam(
"jetalgo"));
550 set.parm(
"JetMatching:qCutME", mad.getParam(
"ptj"));
551 set.mode(
"JetMatching:jetAlgorithm", 2);
552 set.flag(
"JetMatching:doFxFx",
true);
553 set.flag(
"SpaceShower:MEcorrections",
false);
554 set.parm(
"TimeShower:pTmaxMatch", 1);
555 set.parm(
"TimeShower:pTmaxFudge", 1);
556 set.flag(
"TimeShower:MEcorrections",
false);
557 set.flag(
"TimeShower:globalRecoil",
true);
558 set.flag(
"TimeShower:limitPTmaxGlobal",
true);
559 set.mode(
"TimeShower:nMaxGlobalRecoil", 1);
560 set.mode(
"TimeShower:globalRecoilMode", 2);
563 }
else set.parm(
"JetMatching:clFact", mad.getParam(
"alpsfact"));
567 pythia->setUserHooksPtr(hook);
571 if (!run(events))
return false;
572 if (!reader(
true))
return false;
582 bool LHAupMadgraph::setEvent(
int) {
585 if (!pythia)
return false;
587 pythia->info.errorMsg(
"Error from LHAupMadgraph::setEvent: LHAupLHEF "
588 "object not correctly initialized");
return false;}
589 if (!lhef->fileFound()) {
590 pythia->info.errorMsg(
"Error from LHAupMadgraph::setEvent: LHEF "
591 "event file was not found");
return false;}
592 if (!lhef->setEvent()) {
593 if (!run(events))
return false;
594 if (!reader(
false))
return false;
599 setProcess(lhef->idProcess(), lhef->weight(), lhef->scale(),
600 lhef->alphaQED(), lhef->alphaQCD());
601 for (
int ip = 1; ip < lhef->sizePart(); ++ip)
602 addParticle(lhef->id(ip), lhef->status(ip), lhef->mother1(ip),
603 lhef->mother2(ip), lhef->col1(ip), lhef->col2(ip),
604 lhef->px(ip), lhef->py(ip), lhef->pz(ip), lhef->e(ip),
605 lhef->m(ip), lhef->tau(ip), lhef->spin(ip), lhef->scale(ip));
606 setIdX(lhef->id1(), lhef->id2(), lhef->x1(), lhef->x2());
607 setPdf(lhef->id1pdf(), lhef->id2pdf(), lhef->x1pdf(), lhef->x2pdf(),
608 lhef->scalePDF(), lhef->pdf1(), lhef->pdf2(), lhef->pdfIsSet());
617 #endif // Pythia8_LHAMadgraph_H