7 #ifndef Pythia8_PowhegProcs_H
8 #define Pythia8_PowhegProcs_H
10 #include "Pythia8/Pythia.h"
11 #include "Pythia8Plugins/PowhegHooks.h"
34 string pdfIn =
"",
bool random =
true);
40 bool readString(
string line);
43 bool readFile(
string name);
55 typedef void DeleteLHAupPowheg(
LHAup*);
58 string proc, dir, pdf;
61 map<string, string> settings;
100 PowhegProcs::PowhegProcs(
Pythia *pythiaPtrIn,
string procIn,
string dirIn,
101 string pdfIn,
bool random) : lhaup(0), proc(procIn), dir(dirIn), pdf(pdfIn),
102 pythia(pythiaPtrIn), lib(0) {
105 NewLHAupPowheg *sym(0);
106 const char* error(0);
109 lib = dlopen((
"libpythia8powheg" + proc +
".so").c_str(), RTLD_LAZY);
112 pythia->info.errorMsg(
"Error from PowhegProcs::PowhegProcs: "
119 sym = (NewLHAupPowheg*)dlsym(lib,
"newLHAupPowheg");
122 pythia->info.errorMsg(
"Error from PowhegProcs::PowhegProcs: "
126 pythia->settings.addWord(
"POWHEG:dir", dir);
127 pythia->settings.addFlag(
"POWHEG:pythiaRandom", random);
128 if (sym) lhaup = sym(pythia);
131 pythia->setLHAupPtr(lhaup);
132 pythia->setUserHooksPtr(&hooks);
140 PowhegProcs::~PowhegProcs() {
144 DeleteLHAupPowheg *sym(0);
145 sym = (DeleteLHAupPowheg*)dlsym(lib,
"deleteLHAupPowheg");
150 if (lib) {dlclose(lib); dlerror();}
159 bool PowhegProcs::readString(
string line) {
162 if (!pythia)
return false;
163 if (line.find_first_not_of(
" \n\t\v\b\r\f\a") == string::npos)
return true;
164 int firstChar = line.find_first_not_of(
" \n\t\v\b\r\f\a");
165 int lastChar = line.find_last_not_of(
" \n\t\v\b\r\f\a");
166 line = line.substr(firstChar, lastChar + 1 - firstChar);
169 firstChar = line.find_first_of(
" \t\f\v\n\r");
170 string key = toLower( line.substr(0, firstChar),
false);
174 && key.find_first_of(
"abcdedfghijklmnopqrtsuvwxyz") == 0) {
175 map<string, string>::iterator setting = settings.find(key);
176 if (setting != settings.end()) {
177 pythia->info.errorMsg(
"Warning from PowhegProcs::readString: replacing "
178 "previous POWHEG setting for " + key +
".");
179 setting->second = line;
180 }
else settings[key] = line;
190 bool PowhegProcs::readFile(
string name) {
192 fstream config(name.c_str(), ios::in);
string line;
193 while (getline(config, line,
'\n')) readString(line);
203 bool PowhegProcs::init() {
207 fstream pdfin(pdf.c_str(), ios::in | ios::binary);
208 fstream pdfout((dir +
"/" + pdf.substr(0, pdf.find_last_of(
"/"))).c_str(),
209 ios::out | ios::binary);
210 pdfout << pdfin.rdbuf();
216 fstream config((dir +
"/" +
"powheg.input").c_str(), ios::out);
217 for (map<string, string>::iterator setting = settings.begin();
218 setting != settings.end(); ++setting) config << setting->second <<
"\n";
228 #endif // Pythia8_PowhegProcs_H