StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
LesHouches.cc
1 // LesHouches.cc is a part of the PYTHIA event generator.
2 // Copyright (C) 2020 Torbjorn Sjostrand.
3 // PYTHIA is licenced under the GNU GPL v2 or later, see COPYING for details.
4 // Please respect the MCnet Guidelines, see GUIDELINES for details.
5 
6 // Function definitions (not found in the header) for the LHAup and
7 // LHAupLHEF classes.
8 
9 #include "Pythia8/Pythia.h"
10 #include "Pythia8/LesHouches.h"
11 
12 // Access time information.
13 #include <ctime>
14 
15 namespace Pythia8 {
16 
17 //==========================================================================
18 
19 // LHAup class.
20 
21 //--------------------------------------------------------------------------
22 
23 // Constants: could be changed here if desired, but normally should not.
24 // These are of technical nature, as described for each.
25 
26 // LHA convention with cross section in pb may require conversion from mb.
27 const double LHAup::CONVERTMB2PB = 1e9;
28 
29 //--------------------------------------------------------------------------
30 
31 // Print the initialization info; to check it worked.
32 
33 void LHAup::listInit() {
34 
35  // Header.
36  cout << "\n -------- LHA initialization information ------------ \n";
37 
38  // Beam info.
39  cout << fixed << setprecision(3)
40  << "\n beam kind energy pdfgrp pdfset \n"
41  << " A " << setw(6) << idBeamASave
42  << setw(12) << eBeamASave
43  << setw(8) << pdfGroupBeamASave
44  << setw(8) << pdfSetBeamASave << "\n"
45  << " B " << setw(6) << idBeamBSave
46  << setw(12) << eBeamBSave
47  << setw(8) << pdfGroupBeamBSave
48  << setw(8) << pdfSetBeamBSave << "\n";
49 
50  // Event weighting strategy.
51  cout << "\n Event weighting strategy = " << setw(2)
52  << strategySave << "\n" ;
53 
54  // Process list.
55  cout << scientific << setprecision(4)
56  << "\n Processes, with strategy-dependent cross section info \n"
57  << " number xsec (pb) xerr (pb) xmax (pb) \n" ;
58  for (int ip = 0; ip < int(processes.size()); ++ip) {
59  cout << setw(8) << processes[ip].idProc
60  << setw(15) << processes[ip].xSecProc
61  << setw(15) << processes[ip].xErrProc
62  << setw(15) << processes[ip].xMaxProc << "\n";
63  }
64 
65  // Finished.
66  cout << "\n -------- End LHA initialization information -------- \n";
67 
68 }
69 
70 //--------------------------------------------------------------------------
71 
72 // Print the event info; to check it worked.
73 
74 void LHAup::listEvent() {
75 
76  // Header.
77  cout << "\n -------- LHA event information and listing -------------"
78  << "--------------------------------------------------------- \n";
79 
80  // Basic event info.
81  cout << scientific << setprecision(4)
82  << "\n process = " << setw(8) << idProc
83  << " weight = " << setw(12) << weightProc
84  << " scale = " << setw(12) << scaleProc << " (GeV) \n"
85  << " "
86  << " alpha_em = " << setw(12) << alphaQEDProc
87  << " alpha_strong = " << setw(12) << alphaQCDProc << "\n";
88 
89  // Particle list
90  cout << fixed << setprecision(3)
91  << "\n Participating Particles \n"
92  << " no id stat mothers colours p_x "
93  << "p_y p_z e m tau spin \n" ;
94  for (int ip = 1; ip < int(particles.size()); ++ip) {
95  cout << setw(6) << ip
96  << setw(10) << particles[ip].idPart
97  << setw(5) << particles[ip].statusPart
98  << setw(6) << particles[ip].mother1Part
99  << setw(6) << particles[ip].mother2Part
100  << setw(6) << particles[ip].col1Part
101  << setw(6) << particles[ip].col2Part
102  << setw(11) << particles[ip].pxPart
103  << setw(11) << particles[ip].pyPart
104  << setw(11) << particles[ip].pzPart
105  << setw(11) << particles[ip].ePart
106  << setw(11) << particles[ip].mPart
107  << setw(8) << particles[ip].tauPart
108  << setw(8) << particles[ip].spinPart << "\n";
109  }
110 
111  // PDF info - optional.
112  if (pdfIsSetSave) cout << "\n pdf: id1 =" << setw(5) << id1pdfSave
113  << " id2 =" << setw(5) << id2pdfSave
114  << " x1 =" << scientific << setw(10) << x1pdfSave
115  << " x2 =" << setw(10) << x2pdfSave
116  << " scalePDF =" << setw(10) << scalePDFSave
117  << " pdf1 =" << setw(10) << pdf1Save
118  << " pdf2 =" << setw(10) << pdf2Save << "\n";
119 
120  // Finished.
121  cout << "\n -------- End LHA event information and listing ---------"
122  << "--------------------------------------------------------- \n";
123 
124 }
125 
126 //--------------------------------------------------------------------------
127 
128 // Open and write header to a Les Houches Event File.
129 
130 bool LHAup::openLHEF(string fileNameIn) {
131 
132  // Open file for writing. Reset it to be empty.
133  fileName = fileNameIn;
134  const char* cstring = fileName.c_str();
135  osLHEF.open(cstring, ios::out | ios::trunc);
136  if (!osLHEF) {
137  infoPtr->errorMsg("Error in LHAup::openLHEF:"
138  " could not open file", fileName);
139  return false;
140  }
141 
142  // Read out current date and time.
143  time_t t = time(0);
144  strftime(dateNow,12,"%d %b %Y",localtime(&t));
145  strftime(timeNow,9,"%H:%M:%S",localtime(&t));
146 
147  // Write header.
148  osLHEF << "<LesHouchesEvents version=\"1.0\">\n"
149  << "<!--\n"
150  << " File written by Pythia8::LHAup on "
151  << dateNow << " at " << timeNow << "\n"
152  << "-->" << endl;
153 
154  // Done.
155  return true;
156 
157 }
158 
159 //--------------------------------------------------------------------------
160 
161 // Write initialization information to a Les Houches Event File.
162 
163 bool LHAup::initLHEF() {
164 
165  // Write information on beams.
166  osLHEF << "<init>\n" << scientific << setprecision(6)
167  << " " << idBeamASave << " " << idBeamBSave
168  << " " << eBeamASave << " " << eBeamBSave
169  << " " << pdfGroupBeamASave << " " << pdfGroupBeamBSave
170  << " " << pdfSetBeamASave << " " << pdfSetBeamBSave
171  << " " << strategySave << " " << processes.size() << "\n";
172 
173  // Write information on all the subprocesses.
174  for (int ip = 0; ip < int(processes.size()); ++ip)
175  osLHEF << " " << setw(13) << processes[ip].xSecProc
176  << " " << setw(13) << processes[ip].xErrProc
177  << " " << setw(13) << processes[ip].xMaxProc
178  << " " << setw(6) << processes[ip].idProc << "\n";
179 
180  // Done.
181  osLHEF << "</init>" << endl;
182  return true;
183 
184 }
185 
186 //--------------------------------------------------------------------------
187 
188 // Write event information to a Les Houches Event File.
189 // Normal mode is to line up event info in columns, but the non-verbose
190 // alternative saves space at the expense of human readability.
191 
192 bool LHAup::eventLHEF(bool verbose) {
193 
194  // Default verbose option.
195  if (verbose) {
196 
197  // Write information on process as such.
198  osLHEF << "<event>\n" << scientific << setprecision(6)
199  << " " << setw(5) << particles.size() - 1
200  << " " << setw(5) << idProc
201  << " " << setw(13) << weightProc
202  << " " << setw(13) << scaleProc
203  << " " << setw(13) << alphaQEDProc
204  << " " << setw(13) << alphaQCDProc << "\n";
205 
206  // Write information on the particles, excluding zeroth.
207  for (int ip = 1; ip < int(particles.size()); ++ip) {
208  LHAParticle& ptNow = particles[ip];
209  osLHEF << " " << setw(8) << ptNow.idPart
210  << " " << setw(5) << ptNow.statusPart
211  << " " << setw(5) << ptNow.mother1Part
212  << " " << setw(5) << ptNow.mother2Part
213  << " " << setw(5) << ptNow.col1Part
214  << " " << setw(5) << ptNow.col2Part << setprecision(10)
215  << " " << setw(17) << ptNow.pxPart
216  << " " << setw(17) << ptNow.pyPart
217  << " " << setw(17) << ptNow.pzPart
218  << " " << setw(17) << ptNow.ePart
219  << " " << setw(17) << ptNow.mPart << setprecision(6);
220  if (ptNow.tauPart == 0.) osLHEF << " 0.";
221  else osLHEF << " " << setw(13) << ptNow.tauPart;
222  if (ptNow.spinPart == 9.) osLHEF << " 9.";
223  else osLHEF << " " << setw(13) << ptNow.spinPart;
224  osLHEF << "\n";
225  }
226 
227  // Optionally write information on PDF values at hard interaction.
228  if (pdfIsSetSave) osLHEF << "#pdf"
229  << " " << setw(4) << id1pdfSave
230  << " " << setw(4) << id2pdfSave
231  << " " << setw(13) << x1pdfSave
232  << " " << setw(13) << x2pdfSave
233  << " " << setw(13) << scalePDFSave
234  << " " << setw(13) << pdf1Save
235  << " " << setw(13) << pdf2Save << "\n";
236 
237  // Optionally write information on shower scales, primarily in DPS events.
238  if (scaleShowersIsSetSave) osLHEF << "#scaleShowers"
239  << " " << setw(13) << scaleShowersSave[0]
240  << " " << setw(13) << scaleShowersSave[1] << "\n";
241 
242  // Alternative non-verbose option.
243  } else {
244 
245  // Write information on process as such.
246  osLHEF << "<event>\n" << scientific << setprecision(6)
247  << particles.size() - 1 << " " << idProc << " "
248  << weightProc << " " << scaleProc << " "
249  << alphaQEDProc << " " << alphaQCDProc << "\n";
250 
251  // Write information on the particles, excluding zeroth.
252  for (int ip = 1; ip < int(particles.size()); ++ip) {
253  LHAParticle& ptNow = particles[ip];
254  osLHEF << ptNow.idPart << " " << ptNow.statusPart
255  << " " << ptNow.mother1Part << " " << ptNow.mother2Part
256  << " " << ptNow.col1Part << " " << ptNow.col2Part
257  << setprecision(10) << " " << ptNow.pxPart
258  << " " << ptNow.pyPart << " " << ptNow.pzPart
259  << " " << ptNow.ePart << " " << ptNow.mPart
260  << setprecision(6);
261  if (ptNow.tauPart == 0.) osLHEF << " 0.";
262  else osLHEF << " " << setw(13) << ptNow.tauPart;
263  if (ptNow.spinPart == 9.) osLHEF << " 9.";
264  else osLHEF << " " << setw(13) << ptNow.spinPart;
265  osLHEF << "\n";
266  }
267 
268  // Optionally write information on PDF values at hard interaction.
269  if (pdfIsSetSave) osLHEF << "#pdf" << " " << id1pdfSave
270  << " " << id2pdfSave << " " << x1pdfSave << " " << x2pdfSave
271  << " " << scalePDFSave << " " << pdf1Save << " " << pdf2Save
272  << "\n";
273 
274  // Optionally write information on shower scales, primarily in DPS events.
275  if (scaleShowersIsSetSave) osLHEF << "#scaleShowers" << " "
276  << scaleShowersSave[0] << " " << scaleShowersSave[1] << "\n";
277  }
278 
279  // Done.
280  osLHEF << "</event>" << endl;
281  return true;
282 
283 }
284 
285 //--------------------------------------------------------------------------
286 
287 // Write end of a Les Houches Event File and close it.
288 
289 bool LHAup::closeLHEF(bool updateInit) {
290 
291  // Write an end to the file.
292  osLHEF << "</LesHouchesEvents>" << endl;
293  osLHEF.close();
294 
295  // Optionally update the cross section information.
296  if (updateInit) {
297  const char* cstring = fileName.c_str();
298  osLHEF.open(cstring, ios::in | ios::out);
299 
300  // Rewrite header; identically with what openLHEF did.
301  osLHEF << "<LesHouchesEvents version=\"1.0\">\n"
302  << "<!--\n"
303  << " File written by Pythia8::LHAup on "
304  << dateNow << " at " << timeNow << "\n"
305  << "-->" << endl;
306 
307  // Redo initialization information.
308  initLHEF();
309  osLHEF.close();
310  }
311 
312  // Done.
313  return true;
314 
315 }
316 
317 //--------------------------------------------------------------------------
318 
319 // Read in initialization information from a Les Houches Event File.
320 
321 bool LHAup::setInitLHEF(istream& is, bool readHeaders) {
322 
323  // Check that first line is consistent with proper LHEF file.
324  string line;
325  if (!getline(is, line)) return false;
326  if (line.find("<LesHouchesEvents") == string::npos) return false;
327  if (line.find("version=\"1.0\"" ) == string::npos ) return false;
328 
329  // What to search for if reading headers; if not reading
330  // headers then return to default behaviour
331  string headerTag = (readHeaders) ? "<header>" : "<init";
332 
333  // Loop over lines until an <init (or optionally <header>) tag
334  // is found first on a line.
335  string tag = " ";
336  do {
337  if (!getline(is, line)) return false;
338  if (line.find_first_not_of(" \n\t\v\b\r\f\a") != string::npos) {
339  istringstream getfirst(line);
340  getfirst >> tag;
341  if (!getfirst) return false;
342  }
343  } while (tag != "<init>" && tag != "<init" && tag != headerTag);
344 
345  // If header tag found, process if required
346  if (readHeaders == true && tag == headerTag) {
347  // Temporary local storage
348  map < string, string > headerMap;
349 
350  // Loop over lines until an <init> tag is found.
351  bool read = true, newKey = false;
352  string key = "base";
353  vector < string > keyVec;
354 
355  while (true) {
356  if (!getline(is, line)) return false;
357 
358  // Break lines containing multiple tags into two segments.
359  // (Could be generalized to multiple segments but this is
360  // sufficient to handle at least <tag>info</tag> on same line.
361  size_t firstTagEnd = line.find_first_of(">");
362  size_t secondTagBegin = line.find_first_of("<",firstTagEnd);
363  vector<string> lineVec;
364  if (firstTagEnd != string::npos && secondTagBegin != string::npos) {
365  lineVec.push_back(line.substr(0,secondTagBegin));
366  lineVec.push_back(line.substr(secondTagBegin,
367  line.size()-secondTagBegin));
368  }
369  else {
370  lineVec.push_back(line);
371  }
372 
373  // Loop over segments of current line
374  for (int iVec=0; iVec<int(lineVec.size()); ++iVec) {
375  line = lineVec[iVec];
376 
377  // Clean line to contain only valid characters
378  size_t posBeg = line.find_first_not_of(" \n\t\v\b\r\f\a");
379  size_t posEnd = line.find_last_not_of(" \n\t\v\b\r\f\a");
380  string lineClean = " ";
381  if (posBeg != string::npos && posEnd != string::npos && posBeg
382  < posEnd) {
383  lineClean = line.substr(posBeg, posEnd - posBeg + 1);
384  posBeg = 0;
385  posEnd = lineClean.size();
386  }
387 
388  // Check for empty line
389  if (lineClean == " " || posBeg >= posEnd) continue;
390 
391  // PZS Jan 2015: Allow multiple open/close tags on a single line.
392  size_t tagBeg = lineClean.find_first_of("<");
393  size_t tagEnd = lineClean.find_first_of(">");
394 
395  while (tagBeg != string::npos && tagBeg < tagEnd) {
396 
397  // Update remainder (non-tag) part of line, for later storage
398  posBeg = tagEnd+1;
399 
400  // Only take the first word of the tag,
401  tag = lineClean.substr(tagBeg + 1, tagEnd - tagBeg - 1);
402  istringstream getfirst(tag);
403  getfirst >> tag;
404 
405  // Prepare for next while iteration:
406  // Look for next tag on line and update posBeg and posEnd.
407  tagBeg = lineClean.find_first_of("<",tagEnd);
408  tagEnd = lineClean.find_first_of(">",tagBeg+1);
409 
410  // Tag present, so handle here
411  if (getfirst) {
412 
413  // Exit condition
414  if (tag == "init") break;
415 
416  // End of header block; keep reading until <init> tag,
417  // but do not store any further information
418  else if (tag == "/header") {
419  read = false;
420  continue;
421 
422  // Opening tag
423  } else if (tag[0] != '/') {
424  keyVec.push_back(tag);
425  newKey = true;
426  continue;
427 
428  // Closing tag that matches current key
429  } else if (tag == "/" + keyVec.back()) {
430  keyVec.pop_back();
431  newKey = true;
432  continue;
433 
434  // Also check for forgotten close tag: next-to-last element
435  } else if (keyVec.size() >= 2
436  && tag == "/" + keyVec[keyVec.size()-2]) {
437  infoPtr->errorMsg("Warning in LHAup::setInitLHEF:"
438  " corrupt LHEF end tag",keyVec.back());
439  keyVec.pop_back();
440  keyVec.pop_back();
441  newKey = true;
442  continue;
443  }
444 
445  } // if (getfirst)
446 
447  } // Loop over tags
448 
449  // Exit condition
450  if (tag == "init") break;
451 
452  // At this point the (rest of) the line is not a tag;
453  // If no longer reading anything, skip.
454  if (!read) continue;
455 
456  // Check for key change
457  if (newKey) {
458  if (keyVec.empty()) key = "base";
459  else key = keyVec[0];
460  for (size_t i = 1; i < keyVec.size(); i++)
461  key += "." + keyVec[i];
462  newKey = false;
463  }
464 
465  // Check if anything remains to store of this line
466  posBeg = line.find_first_not_of(" \n\t\v\b\r\f\a",posBeg);
467  if (posBeg == string::npos || posBeg > posEnd) continue;
468 
469  // Append information to local storage
470  headerMap[key] += line.substr(posBeg,posEnd - posBeg + 1) + "\n";
471 
472  } // Loop over line segments
473 
474  // Exit condition
475  if (tag == "init") break;
476 
477  } // while (true)
478 
479  // Copy information to info using LHAup::setInfoHeader
480  for (map < string, string >::iterator it = headerMap.begin();
481  it != headerMap.end(); it++)
482  setInfoHeader(it->first, it->second);
483 
484  } // if (readHeaders == true && tag == headerTag)
485 
486  // Read in first info line; done if empty.
487  if (!getline(is, line)) return false;
488  if (line.find("</init") != string::npos) return true;
489 
490  // Read in beam and strategy info, and store it.
491  int idbmupA, idbmupB;
492  double ebmupA, ebmupB;
493  int pdfgupA, pdfgupB, pdfsupA, pdfsupB, idwtup, nprup;
494  istringstream getbms(line);
495  getbms >> idbmupA >> idbmupB >> ebmupA >> ebmupB >> pdfgupA
496  >> pdfgupB >> pdfsupA >> pdfsupB >> idwtup >> nprup;
497  if (!getbms) return false;
498  setBeamA(idbmupA, ebmupA, pdfgupA, pdfsupA);
499  setBeamB(idbmupB, ebmupB, pdfgupB, pdfsupB);
500  setStrategy(idwtup);
501 
502  // Read in process info, one process at a time, and store it.
503  double xsecup, xerrup, xmaxup;
504  xSecSumSave = 0.;
505  xErrSumSave = 0.;
506  int lprup;
507  for (int ip = 0; ip < nprup; ++ip) {
508  if (!getline(is, line)) return false;
509  istringstream getpro(line);
510  getpro >> xsecup >> xerrup >> xmaxup >> lprup ;
511  if (!getpro) return false;
512  addProcess(lprup, xsecup, xerrup, xmaxup);
513  xSecSumSave += xsecup;
514  xErrSumSave += pow2(xerrup);
515  }
516  xErrSumSave = sqrt(xErrSumSave);
517 
518  // Reading worked.
519  return true;
520 
521 }
522 
523 //--------------------------------------------------------------------------
524 
525 // Read in event information from a Les Houches Event File,
526 // into a staging area where it can be reused by setOldEventLHEF.
527 
528 bool LHAup::setNewEventLHEF(istream& is) {
529 
530  // Loop over lines until an <event tag is found first on a line.
531  string line, tag;
532  do {
533  if (!getline(is, line)) return false;
534  if (line.find_first_not_of(" \n\t\v\b\r\f\a") != string::npos) {
535  istringstream getfirst(line);
536  getfirst >> tag;
537  if (!getfirst) return false;
538  }
539  } while (tag != "<event>" && tag != "<event");
540 
541  // Read in process info and store it.
542  if (!getline(is, line)) return false;
543  istringstream getpro(line);
544  getpro >> nupSave >> idprupSave >> xwgtupSave >> scalupSave
545  >> aqedupSave >> aqcdupSave;
546  if (!getpro) return false;
547 
548  // Reset particlesSave vector, add slot-0 empty particle.
549  particlesSave.clear();
550  particlesSave.push_back( LHAParticle() );
551 
552  // Read in particle info one by one, and store it.
553  // Note unusual C++ loop range, to better reflect LHA/Fortran standard.
554  // (Recall that process(...) above added empty particle at index 0.)
555  int idup, istup, mothup1, mothup2, icolup1, icolup2;
556  double pup1, pup2, pup3, pup4, pup5, vtimup, spinup;
557  for (int ip = 1; ip <= nupSave; ++ip) {
558  if (!getline(is, line)) return false;
559  istringstream getall(line);
560  getall >> idup >> istup >> mothup1 >> mothup2 >> icolup1 >> icolup2
561  >> pup1 >> pup2 >> pup3 >> pup4 >> pup5 >> vtimup >> spinup;
562  if (!getall) return false;
563  particlesSave.push_back( LHAParticle( idup, istup, mothup1, mothup2,
564  icolup1, icolup2, pup1, pup2, pup3, pup4, pup5, vtimup, spinup, -1.) );
565  }
566 
567  // Flavour and x values of hard-process initiators.
568  id1InSave = particlesSave[1].idPart;
569  id2InSave = particlesSave[2].idPart;
570  x1InSave = (eBeamASave > 0.) ? particlesSave[1].ePart / eBeamASave : 0.;
571  x2InSave = (eBeamBSave > 0.) ? particlesSave[2].ePart / eBeamBSave : 0.;
572 
573  // Continue parsing till </event>. Look for optional info on the way.
574  getPDFSave = false;
575  getScale = false;
576  getScaleShowers = false;
577  do {
578  if (!getline(is, line)) return false;
579  istringstream getinfo(line);
580  getinfo >> tag;
581  if (!getinfo) return false;
582 
583  // Extract PDF info if present.
584  if (tag == "#pdf" && !getPDFSave) {
585  getinfo >> id1pdfInSave >> id2pdfInSave >> x1pdfInSave >> x2pdfInSave
586  >> scalePDFInSave >> pdf1InSave >> pdf2InSave;
587  if (!getinfo) return false;
588  getPDFSave = true;
589 
590  // Extract shower scales info if present.
591  } else if (tag == "#scaleShowers") {
592  getinfo >> scaleShowersInSave[0] >> scaleShowersInSave[1];
593  if (!getinfo) return false;
594  getScaleShowers = true;
595 
596  // Extract scale info if present.
597  } else if (tag == "#" && !getScale) {
598  double scaleIn = 0;
599  for (int i = 3; i < int(particlesSave.size()); ++i)
600  if (particlesSave[i].statusPart == 1) {
601  if ( !(getinfo >> scaleIn) ) return false;
602  particlesSave[i].scalePart = scaleIn;
603  }
604  if (!getinfo) return false;
605  getScale = true;
606  }
607  } while (tag != "</event>" && tag != "</event");
608 
609  // Need id and x values even when no PDF info. Rest empty.
610  if (!getPDFSave) {
611  id1pdfInSave = id1InSave;
612  id2pdfInSave = id2InSave;
613  x1pdfInSave = x1InSave;
614  x2pdfInSave = x2InSave;
615  scalePDFInSave = 0.;
616  pdf1InSave = 0.;
617  pdf2InSave = 0.;
618  }
619 
620  // Reading worked.
621  return true;
622 
623 }
624 
625 //--------------------------------------------------------------------------
626 
627 // Make current event information read in by setNewEventLHEF.
628 
629 bool LHAup::setOldEventLHEF() {
630 
631  // Store saved event, optionally also parton density and scale information.
632  setProcess( idprupSave, xwgtupSave, scalupSave, aqedupSave, aqcdupSave);
633  for (int ip = 1; ip <= nupSave; ++ip) addParticle( particlesSave[ip] );
634  setIdX( id1InSave, id2InSave, x1InSave, x2InSave);
635  setPdf( id1pdfInSave, id2pdfInSave, x1pdfInSave, x2pdfInSave,
636  scalePDFInSave, pdf1InSave, pdf2InSave, getPDFSave);
637  if (getScaleShowers)
638  setScaleShowers( scaleShowersInSave[0], scaleShowersInSave[1]);
639 
640  // Done;
641  return true;
642 
643 }
644 
645 //--------------------------------------------------------------------------
646 
647 // Open a file using provided ifstream and return a pointer to an istream
648 // that can be used to process the file.
649 
650 istream* LHAup::openFile(const char *fn, ifstream &ifs) {
651  // Open the file
652  ifs.open(fn);
653  return (istream *) &ifs;
654 }
655 
656 //--------------------------------------------------------------------------
657 
658 // Companion method to 'openFile', above.
659 // Correctly deallocates memory if required before closing the file.
660 
661 void LHAup::closeFile(istream *&is, ifstream &ifs) {
662  // If the istream pointer is not NULL and is not the
663  // same as the ifstream, then delete pointer.
664  if (is && is != &ifs) delete is;
665  is = NULL;
666 
667  // Close the file
668  if (ifs.is_open()) ifs.close();
669 }
670 
671 //==========================================================================
672 
673 // LHAupLHEF class.
674 
675 //--------------------------------------------------------------------------
676 
677 // Routine for doing the job of reading and setting initialization info.
678 
679 bool LHAupLHEF::setInitLHEF( istream & isIn, bool readHead ) {
680 
681  // Done if there was a problem with initialising the reader
682  if (!reader.isGood) return false;
683 
684  // Construct header information (stored in comments strings or optional
685  // header file), so that reading of headers is possible.
686  string comments;
687  comments+="<LesHouchesEvents version =\"3.0\">\n";
688  comments+="<header>\n";
689  comments+=reader.headerComments;
690  comments+="</header>\n";
691  comments+="<init>\n";
692  comments+=reader.initComments;
693  comments+="</init>\n";
694  istringstream is1(comments);
695  bool useComments = (headerfile == NULL);
696  istream & iss((useComments ? is1 : isIn));
697 
698  // Check that first line is consistent with proper LHEF file.
699  string line;
700  if ( useComments && !getline(iss,line)) return false;
701  if (!useComments && !getLine(line)) return false;
702 
703  // What to search for if reading headers; if not reading
704  // headers then return to default behaviour
705  string headerTag = (readHead) ? "<header>" : "<init";
706 
707  // Loop over lines until an <init (or optionally <header>) tag
708  // is found first on a line.
709  string tag = " ";
710  do {
711  if ( useComments && !getline(iss,line)) return false;
712  if (!useComments && !getLine(line)) return false;
713  if (line.find_first_not_of(" \n\t\v\b\r\f\a") != string::npos) {
714  istringstream getfirst(line);
715  getfirst >> tag;
716  if (!getfirst) return false;
717  }
718  } while (tag != "<init>" && tag != "<init" && tag != headerTag);
719 
720  // If header tag found, process if required
721  if (readHead == true && tag == headerTag) {
722  // Temporary local storage
723  map < string, string > headerMap;
724 
725  // Loop over lines until an <init> tag is found.
726  bool read = true, newKey = false;
727  int commentDepth = 0;
728  string key = "base";
729  vector < string > keyVec;
730  while (true) {
731  if ( useComments && !getline(iss,line)) return false;
732  if (!useComments && !getLine(line)) return false;
733 
734  // Tell XML parser to ignore comment and CDATA blocks
735  // If we are currently inside a comment block, check for block end
736  if (commentDepth >= 1 && line.find("-->") != string::npos) {
737  commentDepth--;
738  size_t comBeg = line.find("-->")+2;
739  size_t comEnd = line.find_last_not_of("\n\t\v\b\r\f\a ");
740  if( comEnd == comBeg ) continue;
741  line = line.substr(comBeg,comEnd-comBeg+1);
742  }
743  if (commentDepth >= 1 && line.find("]]>") != string::npos)
744  commentDepth--;
745  // If the comment block did not end on this line, skip to next line
746  if (commentDepth >= 1) continue;
747  // Check for beginning of comment blocks (parse until comment begins)
748  if (line.find("<!--") != string::npos) {
749  if (line.find("-->") == string::npos) commentDepth++;
750  int comBeg = line.find("<!--");
751  line = line.substr(0,comBeg);
752  }
753  // Check for beginning of CDATA statement (parse until CDATA begins)
754  if (line.find("<![cdata[") != string::npos
755  || line.find("<![CDATA[") != string::npos) {
756  if (line.find("]]>") == string::npos) commentDepth++;
757  int comBeg = line.find("<![");
758  line = line.substr(0,comBeg);
759  }
760 
761  // Break lines containing multiple tags into two segments.
762  // (Could be generalized to multiple segments but this is
763  // sufficient to handle at least <tag>info</tag> on same line.
764  size_t firstTagEnd = line.find_first_of(">");
765  size_t secondTagBegin = line.find_first_of("<",firstTagEnd);
766  vector<string> lineVec;
767  if (firstTagEnd != string::npos && secondTagBegin != string::npos) {
768  lineVec.push_back(line.substr(0,secondTagBegin));
769  lineVec.push_back(line.substr(secondTagBegin,
770  line.size()-secondTagBegin));
771  }
772  else {
773  lineVec.push_back(line);
774  }
775 
776  // Loop over segments of current line
777  for (int iVec=0; iVec<int(lineVec.size()); ++iVec) {
778  line = lineVec[iVec];
779 
780  // Clean line to contain only valid characters
781  size_t posBeg = line.find_first_not_of(" \n\t\v\b\r\f\a");
782  size_t posEnd = line.find_last_not_of(" \n\t\v\b\r\f\a");
783  string lineClean = " ";
784  if (posBeg != string::npos && posEnd != string::npos && posBeg
785  < posEnd) {
786  lineClean = line.substr(posBeg, posEnd - posBeg + 1);
787  posBeg = 0;
788  posEnd = lineClean.size();
789  line = lineClean;
790  }
791 
792  // Check for empty line
793  if (lineClean == " " || posBeg >= posEnd) continue;
794 
795  // PZS Jan 2015: Allow multiple open/close tags on a single line.
796  size_t tagBeg = lineClean.find_first_of("<");
797  size_t tagEnd = lineClean.find_first_of(">");
798 
799  while (tagBeg != string::npos && tagBeg < tagEnd) {
800 
801  // Update remainder (non-tag) part of line, for later storage
802  posBeg = tagEnd+1;
803 
804  // Only take the first word of the tag,
805  tag = lineClean.substr(tagBeg + 1, tagEnd - tagBeg - 1);
806  istringstream getfirst(tag);
807  getfirst >> tag;
808 
809  // Prepare for next while iteration:
810  // Look for next tag on line and update posBeg and posEnd.
811  tagBeg = lineClean.find_first_of("<",tagEnd);
812  tagEnd = lineClean.find_first_of(">",tagBeg+1);
813 
814  // Tag present, so handle here
815  if (getfirst) {
816 
817  // Exit condition
818  if (tag == "init") break;
819 
820  // End of header block; keep reading until <init> tag,
821  // but do not store any further information
822  else if (tag == "/header") {
823  read = false;
824  continue;
825 
826  // Opening tag
827  } else if (tag[0] != '/') {
828  keyVec.push_back(tag);
829  newKey = true;
830  continue;
831 
832  // Closing tag that matches current key
833  } else if (tag == "/" + keyVec.back()) {
834  keyVec.pop_back();
835  newKey = true;
836  continue;
837 
838  // Also check for forgotten close tag: next-to-last element
839  } else if (keyVec.size() >= 2
840  && tag == "/" + keyVec[keyVec.size()-2]) {
841  infoPtr->errorMsg("Warning in LHAupLHEF::setInitLHEF:"
842  " corrupt LHEF end tag",keyVec.back());
843  keyVec.pop_back();
844  keyVec.pop_back();
845  newKey = true;
846  continue;
847  }
848 
849  } // if (getfirst)
850 
851  } // Loop over tags
852 
853  // Exit condition
854  if (tag == "init") break;
855 
856  // At this point the (rest of) the line is not a tag;
857  // If no longer reading anything, skip.
858  if (!read) continue;
859 
860  // Check for key change
861  if (newKey) {
862  if (keyVec.empty()) key = "base";
863  else key = keyVec[0];
864  for (size_t i = 1; i < keyVec.size(); i++)
865  key += "." + keyVec[i];
866  newKey = false;
867  }
868 
869  // Check if anything remains to store of this line
870  posBeg = line.find_first_not_of(" \n\t\v\b\r\f\a",posBeg);
871  if (posBeg == string::npos || posBeg > posEnd) continue;
872 
873  // Append information to local storage
874  headerMap[key] += line.substr(posBeg,posEnd - posBeg + 1) + "\n";
875 
876  } // Loop over line segments
877 
878  // Exit condition
879  if (tag == "init") break;
880 
881  } // while (true)
882 
883  // Copy information to info using LHAup::setInfoHeader
884  for (map < string, string >::iterator it = headerMap.begin();
885  it != headerMap.end(); it++)
886  setInfoHeader(it->first, it->second);
887 
888  } // if (readHead == true && tag == headerTag)
889 
890  // Extract beam and strategy info, and store it.
891  int idbmupA, idbmupB;
892  double ebmupA, ebmupB;
893  int pdfgupA, pdfgupB, pdfsupA, pdfsupB, idwtup, nprup;
894 
895  idbmupA = reader.heprup.IDBMUP.first;
896  idbmupB = reader.heprup.IDBMUP.second;
897  ebmupA = reader.heprup.EBMUP.first;
898  ebmupB = reader.heprup.EBMUP.second;
899  pdfgupA = reader.heprup.PDFGUP.first;
900  pdfgupB = reader.heprup.PDFGUP.first;
901  pdfsupA = reader.heprup.PDFSUP.first;
902  pdfsupB = reader.heprup.PDFSUP.second;
903  idwtup = reader.heprup.IDWTUP;
904  nprup = reader.heprup.NPRUP;
905 
906  setBeamA(idbmupA, ebmupA, pdfgupA, pdfsupA);
907  setBeamB(idbmupB, ebmupB, pdfgupB, pdfsupB);
908  setStrategy(idwtup);
909 
910  // Read in process info, one process at a time, and store it.
911  double xsecup, xerrup, xmaxup;
912  xSecSumSave = 0.;
913  xErrSumSave = 0.;
914  int lprup;
915  infoPtr->sigmaLHEFSave.resize(0);
916  for (int ip = 0; ip < nprup; ++ip) {
917  xsecup = reader.heprup.XSECUP[ip];
918  xerrup = reader.heprup.XERRUP[ip];
919  xmaxup = reader.heprup.XMAXUP[ip];
920  lprup = reader.heprup.LPRUP[ip];
921  addProcess(lprup, xsecup, xerrup, xmaxup);
922  xSecSumSave += xsecup;
923  xErrSumSave += pow(xerrup,2);
924  infoPtr->sigmaLHEFSave.push_back(xsecup);
925  }
926  xErrSumSave = sqrt(xErrSumSave);
927 
928  // Now extract LHEF 2.0/3.0 novelties.
929  infoPtr->setLHEF3InitInfo();
930  if (reader.version > 1) {
931  infoPtr->setLHEF3InitInfo( reader.version,
932  &reader.heprup.initrwgt, &(reader.heprup.generators),
933  &(reader.heprup.weightgroups), &(reader.heprup.weights),
934  reader.headerBlock);
935  }
936 
937  // Reading worked.
938  return true;
939 }
940 
941 //--------------------------------------------------------------------------
942 
943 // Routine for doing the job of reading and setting info on next event.
944 
945 bool LHAupLHEF::setNewEventLHEF() {
946 
947  // Done if the reader finished preemptively.
948  if(!reader.readEvent()) return false;
949 
950  // Extract process info and store it.
951  nupSave = reader.hepeup.NUP;
952  idprupSave = reader.hepeup.IDPRUP;
953  xwgtupSave = reader.hepeup.XWGTUP;
954  scalupSave = reader.hepeup.SCALUP;
955  aqedupSave = reader.hepeup.AQEDUP;
956  aqcdupSave = reader.hepeup.AQCDUP;
957 
958  // Reset particlesSave vector, add slot-0 empty particle.
959  particlesSave.clear();
960  particlesSave.push_back( Pythia8::LHAParticle() );
961 
962  // Read in particle info one by one, and store it.
963  // Note unusual C++ loop range, to better reflect LHA/Fortran standard.
964  // (Recall that process(...) above added empty particle at index 0.)
965  int idup, istup, mothup1, mothup2, icolup1, icolup2;
966  double pup1, pup2, pup3, pup4, pup5, vtimup, spinup;
967  for ( int i = 0; i < reader.hepeup.NUP; ++i ) {
968  // Extract information stored in reader.
969  idup = reader.hepeup.IDUP[i];
970  istup = reader.hepeup.ISTUP[i];
971  mothup1 = reader.hepeup.MOTHUP[i].first;
972  mothup2 = reader.hepeup.MOTHUP[i].second;
973  icolup1 = reader.hepeup.ICOLUP[i].first;
974  icolup2 = reader.hepeup.ICOLUP[i].second;
975  pup1 = reader.hepeup.PUP[i][0];
976  pup2 = reader.hepeup.PUP[i][1];
977  pup3 = reader.hepeup.PUP[i][2];
978  pup4 = reader.hepeup.PUP[i][3];
979  pup5 = reader.hepeup.PUP[i][4];
980  vtimup = reader.hepeup.VTIMUP[i];
981  spinup = reader.hepeup.SPINUP[i];
982  particlesSave.push_back( Pythia8::LHAParticle( idup,istup,mothup1,mothup2,
983  icolup1, icolup2, pup1, pup2, pup3, pup4, pup5, vtimup, spinup, -1.) );
984  }
985 
986  // Flavour and x values of hard-process initiators.
987  id1InSave = particlesSave[1].idPart;
988  id2InSave = particlesSave[2].idPart;
989  x1InSave = (eBeamA() > 0.) ? particlesSave[1].ePart / eBeamA() : 0.;
990  x2InSave = (eBeamB() > 0.) ? particlesSave[2].ePart / eBeamB() : 0.;
991 
992  // Parse event comments and look for optional info on the way.
993  std::string line, tag;
994  std::stringstream ss(reader.eventComments);
995  getPDFSave = false;
996  getScale = (setScalesFromLHEF && reader.version == 1) ? false : true;
997  getScaleShowers = false;
998  while (getline(ss, line)) {
999  istringstream getinfo(line);
1000  getinfo >> tag;
1001  if (!getinfo) break;
1002  // Extract PDF info if present.
1003  if (tag == "#pdf" && !getPDFSave) {
1004  getinfo >> id1pdfInSave >> id2pdfInSave >> x1pdfInSave >> x2pdfInSave
1005  >> scalePDFInSave >> pdf1InSave >> pdf2InSave;
1006  if (!getinfo) return false;
1007  getPDFSave = true;
1008  // Extract shower scales info if present.
1009  } else if (tag == "#scaleShowers") {
1010  getinfo >> scaleShowersInSave[0] >> scaleShowersInSave[1];
1011  if (!getinfo) return false;
1012  getScaleShowers = true;
1013  // Extract scale info if present.
1014  } else if (tag == "#" && !getScale) {
1015  double scaleIn = 0;
1016  for (int i = 3; i < int(particlesSave.size()); ++i)
1017  if (particlesSave[i].statusPart == 1) {
1018  if ( !(getinfo >> scaleIn) ) return false;
1019  particlesSave[i].scalePart = scaleIn;
1020  }
1021  if (!getinfo) return false;
1022  getScale = true;
1023  }
1024  }
1025 
1026  // Set production scales from <scales> tag.
1027  if ( setScalesFromLHEF && reader.version > 1 ){
1028  for ( map<string,double>::const_iterator
1029  it = reader.hepeup.scalesSave.attributes.begin();
1030  it != reader.hepeup.scalesSave.attributes.end(); ++it ) {
1031  if ( it->first.find_last_of("_") != string::npos) {
1032  unsigned iFound = it->first.find_last_of("_") + 1;
1033  int iPos = atoi(it->first.substr(iFound).c_str());
1034  // Only set production scales of final particles.
1035  if ( iPos < int(particlesSave.size())
1036  && particlesSave[iPos].statusPart == 1)
1037  particlesSave[iPos].scalePart = it->second;
1038  }
1039  }
1040  }
1041 
1042  // Need id and x values even when no PDF info. Rest empty.
1043  if (!getPDFSave) {
1044  id1pdfInSave = id1InSave;
1045  id2pdfInSave = id2InSave;
1046  x1pdfInSave = x1InSave;
1047  x2pdfInSave = x2InSave;
1048  scalePDFInSave = 0.;
1049  pdf1InSave = 0.;
1050  pdf2InSave = 0.;
1051  }
1052 
1053  // Now extract LHEF 2.0/3.0 novelties.
1054  infoPtr->setLHEF3EventInfo();
1055  // Set everything for 2.0 and 3.0
1056  if (reader.version > 1) {
1057  infoPtr->setLHEF3EventInfo( &reader.hepeup.attributes,
1058  &reader.hepeup.weights_detailed, &reader.hepeup.weights_compressed,
1059  &reader.hepeup.scalesSave, &reader.hepeup.weightsSave,
1060  &reader.hepeup.rwgtSave, reader.weights_detailed_vector(),
1061  reader.weightnames_detailed_vector(),
1062  reader.eventComments, reader.hepeup.XWGTUP);
1063  // Try to at least set the event attributes for 1.0
1064  } else {
1065  infoPtr->setLHEF3EventInfo( &reader.hepeup.attributes, 0, 0, 0, 0, 0,
1066  vector<double>(), vector<string>(), "", reader.hepeup.XWGTUP);
1067  }
1068 
1069  // Reading worked.
1070  return true;
1071 
1072 }
1073 
1074 //==========================================================================
1075 
1076 // LHAupPlugin class.
1077 
1078 //--------------------------------------------------------------------------
1079 
1080 // Constructor.
1081 
1082 LHAupPlugin::LHAupPlugin(string nameIn, Pythia *pythiaPtr) :
1083  LHAup(), lhaPtr(nullptr), libPtr(nullptr), name(nameIn) {
1084 
1085  // Load the plugin library if needed.
1086  if (infoPtr != nullptr)
1087  libPtr = const_cast<Info&>(pythiaPtr->info).plugin(name);
1088  else libPtr = make_shared<Plugin>(name);
1089  if (!libPtr->isLoaded()) return;
1090 
1091  // Create a new LHAup.
1092  NewLHAup* newLHAup = (NewLHAup*)libPtr->symbol("newLHAup");
1093  if (!newLHAup) return;
1094  lhaPtr = newLHAup(pythiaPtr);
1095 
1096 }
1097 
1098 //--------------------------------------------------------------------------
1099 
1100 // Destructor.
1101 
1102 LHAupPlugin::~LHAupPlugin() {
1103 
1104  // Delete the LHAup pointer.
1105  if (lhaPtr == nullptr || !libPtr->isLoaded()) return;
1106  DeleteLHAup* deleteLHAup =
1107  (DeleteLHAup*)libPtr->symbol("deleteLHAup");
1108  if (deleteLHAup) deleteLHAup(lhaPtr);
1109 
1110 }
1111 
1112 //==========================================================================
1113 
1114 // LHAupFromPYTHIA8 class.
1115 
1116 //--------------------------------------------------------------------------
1117 
1118 // Read in initialization information from PYTHIA 8.
1119 
1120 bool LHAupFromPYTHIA8::setInit() {
1121 
1122  // Read in beam from Info class. Parton density left empty.
1123  int idbmupA = infoPtr->idA();
1124  int idbmupB = infoPtr->idB();
1125  double ebmupA = infoPtr->eA();
1126  double ebmupB = infoPtr->eB();
1127  int pdfgupA = 0;
1128  int pdfgupB = 0;
1129  int pdfsupA = 0;
1130  int pdfsupB = 0;
1131  setBeamA(idbmupA, ebmupA, pdfgupA, pdfsupA);
1132  setBeamB(idbmupB, ebmupB, pdfgupB, pdfsupB);
1133 
1134  // Currently only one allowed strategy.
1135  int idwtup = 3;
1136  setStrategy(idwtup);
1137 
1138  // Only one process with dummy information. (Can overwrite at the end.)
1139  int lprup = 9999;
1140  double xsecup = 1.;
1141  double xerrup = 0.;
1142  double xmaxup = 1.;
1143  addProcess(lprup, xsecup, xerrup, xmaxup);
1144 
1145  // Done.
1146  return true;
1147 
1148 }
1149 
1150 //--------------------------------------------------------------------------
1151 
1152 // Read in event information from PYTHIA 8.
1153 
1154 bool LHAupFromPYTHIA8::setEvent( int) {
1155 
1156  // Read process information from Info class, and store it.
1157  // Note: renormalization scale here, factorization further down.
1158  // For now always convert to process 9999, instead of infoPtr->code().
1159  int idprup = 9999;
1160  double xwgtup = infoPtr->weight();
1161  double scalup = infoPtr->QRen();
1162  double aqedup = infoPtr->alphaEM();
1163  double aqcdup = infoPtr->alphaS();
1164  setProcess(idprup, xwgtup, scalup, aqedup, aqcdup);
1165 
1166  // Read in particle info one by one, excluding zero and beams, and store it.
1167  // Note unusual C++ loop range, to better reflect LHA/Fortran standard.
1168  int nup = processPtr->size() - 3;
1169  int n21 = 0;
1170  int idup, statusup, istup, mothup1, mothup2, icolup1, icolup2;
1171  double pup1, pup2, pup3, pup4, pup5, vtimup, spinup;
1172  for (int ip = 1; ip <= nup; ++ip) {
1173  Particle& particle = (*processPtr)[ip + 2];
1174  idup = particle.id();
1175  // Convert from PYTHIA8 to LHA status codes.
1176  statusup = particle.status();
1177  if (statusup == -21) istup = -1;
1178  else if (statusup < 0) istup = 2;
1179  else istup = 1;
1180  mothup1 = max(0, particle.mother1() - 2);
1181  mothup2 = max(0, particle.mother2() - 2);
1182  icolup1 = particle.col();
1183  icolup2 = particle.acol();
1184  pup1 = particle.px();
1185  pup2 = particle.py();
1186  pup3 = particle.pz();
1187  pup4 = particle.e();
1188  pup5 = particle.m();
1189  vtimup = particle.tau();
1190  spinup = particle.pol();
1191  addParticle(idup, istup, mothup1, mothup2, icolup1, icolup2,
1192  pup1, pup2, pup3, pup4, pup5, vtimup, spinup, -1.) ;
1193  if (statusup == -21) ++n21;
1194  }
1195 
1196  // Extract hard-process initiator information from Info class, and store it.
1197  int id1up = infoPtr->id1();
1198  int id2up = infoPtr->id2();
1199  double x1up = infoPtr->x1();
1200  double x2up = infoPtr->x2();
1201  setIdX( id1up, id2up, x1up, x2up);
1202 
1203  // Also extract pdf information from Info class, and store it.
1204  int id1pdfup = infoPtr->id1pdf();
1205  int id2pdfup = infoPtr->id2pdf();
1206  double x1pdfup = infoPtr->x1pdf();
1207  double x2pdfup = infoPtr->x2pdf();
1208  double scalePDFup = infoPtr->QFac();
1209  double pdf1up = infoPtr->pdf1();
1210  double pdf2up = infoPtr->pdf2();
1211  setPdf(id1pdfup, id2pdfup, x1pdfup, x2pdfup, scalePDFup, pdf1up,
1212  pdf2up, true);
1213 
1214  // Extract parton shower scales from for DPS processes.
1215  if (n21 == 4) setScaleShowers( processPtr->scale(),
1216  processPtr->scaleSecond() );
1217 
1218  // Done.
1219  return true;
1220 
1221 }
1222 
1223 //--------------------------------------------------------------------------
1224 
1225 // Update cross-section information at the end of the run.
1226 
1227 bool LHAupFromPYTHIA8::updateSigma() {
1228 
1229  // Read out information from PYTHIA 8 and send it in to LHA.
1230  double sigGen = CONVERTMB2PB * infoPtr->sigmaGen();
1231  double sigErr = CONVERTMB2PB * infoPtr->sigmaErr();
1232  setXSec(0, sigGen);
1233  setXErr(0, sigErr);
1234 
1235  // Done.
1236  return true;
1237 
1238 }
1239 
1240 //==========================================================================
1241 
1242 // LHEF3FromPythia8 class.
1243 
1244 //--------------------------------------------------------------------------
1245 
1246 // Function to open the output file stream.
1247 
1248 bool LHEF3FromPythia8::openLHEF(string fileNameIn) {
1249 
1250  // Open file for writing. Reset it to be empty.
1251  fileName = fileNameIn;
1252  const char* cstring = fileName.c_str();
1253  osLHEF.open(cstring, ios::out | ios::trunc);
1254  if (!osLHEF) {
1255  cout << "Error in LHAup::openLHEF: could not open file "
1256  << fileName << endl;
1257  return false;
1258  }
1259 
1260  // Done.
1261  return true;
1262 }
1263 
1264 //--------------------------------------------------------------------------
1265 
1266 // Routine for reading, setting and printing the initialisation info.
1267 
1268 bool LHEF3FromPythia8::setInit() {
1269 
1270  // Start with clean writer.
1271  writer.headerStream.str("");
1272  writer.initStream.str("");
1273  writer.headerStream.clear();
1274  writer.initStream.clear();
1275 
1276  // PDG id's of beam particles. (first/second is in +/-z direction).
1277  heprup.IDBMUP = make_pair(infoPtr->idA(), infoPtr->idB());
1278 
1279  // Energy of beam particles given in GeV.
1280  heprup.EBMUP = make_pair(infoPtr->eA(),infoPtr->eB());
1281 
1282  // The author group for the PDF used for the beams according to the
1283  // PDFLib specification.
1284  heprup.PDFGUP = make_pair(0,0);
1285 
1286  // The id number the PDF used for the beams according to the
1287  // PDFLib specification.
1288  heprup.PDFSUP = make_pair(0,0);
1289 
1290  // Master switch indicating how the ME generator envisages the
1291  // events weights should be interpreted according to the Les Houches
1292  // accord.
1293  heprup.IDWTUP = -4;
1294 
1295  // The number of different subprocesses in this file.
1296  heprup.NPRUP = 1;
1297 
1298  // The cross sections for the different subprocesses in pb.
1299  vector<double> XSECUP;
1300  for ( int i=0; i < heprup.NPRUP; ++i)
1301  XSECUP.push_back(CONVERTMB2PB * (infoPtr->sigmaGen()));
1302  heprup.XSECUP = XSECUP;
1303 
1304  // The statistical error in the cross sections for the different
1305  // subprocesses in pb.
1306  vector<double> XERRUP;
1307  for ( int i=0; i < heprup.NPRUP; ++i)
1308  XERRUP.push_back(CONVERTMB2PB * infoPtr->sigmaErr());
1309  heprup.XERRUP = XERRUP;
1310 
1311  // The maximum event weights (in HEPEUP::XWGTUP) for different
1312  vector<double> XMAXUP;
1313  for ( int i=0; i < heprup.NPRUP; ++i) XMAXUP.push_back(0.0);
1314  heprup.XMAXUP = XMAXUP;
1315 
1316  // The subprocess code for the different subprocesses.
1317  vector<int> LPRUP;
1318  for ( int i=0; i < heprup.NPRUP; ++i) LPRUP.push_back(9999+i);
1319  heprup.LPRUP = LPRUP;
1320 
1321  // Contents of the LHAinitrwgt tag
1322  if (infoPtr->initrwgt )
1323  heprup.initrwgt = *(infoPtr->initrwgt);
1324 
1325  // Contents of the LHAgenerator tags.
1326  if (infoPtr->generators)
1327  heprup.generators = *(infoPtr->generators);
1328 
1329  // A map of the LHAweightgroup tags, indexed by name.
1330  if (infoPtr->weightgroups)
1331  heprup.weightgroups = *(infoPtr->weightgroups);
1332 
1333  // A map of the LHAweight tags, indexed by name.
1334  if (infoPtr->init_weights)
1335  heprup.weights = *(infoPtr->init_weights);
1336 
1337  // Get init information.
1338  writer.version = 3;
1339 
1340  string line, tag;
1341 
1342  // Not implemented yet:
1343  // Write header block of input LHEF
1344  // Write header comments of input LHEF
1345 
1346  // Print Pythia settings
1347  stringstream setout;
1348  settingsPtr->writeFile(setout, true);
1349  while ( getline(setout,line) )
1350  writer.headerBlock() << line << "\n";
1351 
1352  // Not implemented yet:
1353  // Write init comments of input LHEF.
1354 
1355  writer.heprup = heprup;
1356  writer.init();
1357 
1358  // Done
1359  return true;
1360 }
1361 
1362 //--------------------------------------------------------------------------
1363 
1364 // Routine for reading, setting and printing the next event.
1365 
1366 bool LHEF3FromPythia8::setEvent(int) {
1367 
1368  Event event = *eventPtr;
1369 
1370  // Begin filling Les Houches blocks.
1371  hepeup.clear();
1372  hepeup.resize(0);
1373 
1374  // The number of particle entries in the current event.
1375  hepeup.NUP = 2;
1376  for ( int i = 0; i < int(event.size()); ++i) {
1377  if ( event[i].status() == -22) ++hepeup.NUP;
1378  if ( event[i].isFinal()) ++hepeup.NUP;
1379  }
1380 
1381  // The subprocess code for this event (as given in LPRUP).
1382  hepeup.IDPRUP = 9999;
1383 
1384  // The weight for this event.
1385  hepeup.XWGTUP = infoPtr->weight();
1386 
1387  // The PDF weights for the two incoming partons. Note that this
1388  // variable is not present in the current LesHouches accord
1389  // (<A HREF="http://arxiv.org/abs/hep-ph/0109068">hep-ph/0109068</A>),
1390  // hopefully it will be present in a future accord.
1391  hepeup.XPDWUP = make_pair(0,0);
1392 
1393  // The scale in GeV used in the calculation of the PDF's in this
1394  // event.
1395  hepeup.SCALUP = eventPtr->scale();
1396 
1397  // The value of the QED coupling used in this event.
1398  hepeup.AQEDUP = infoPtr->alphaEM();
1399 
1400  // The value of the QCD coupling used in this event.
1401  hepeup.AQCDUP = infoPtr->alphaS();
1402 
1403  // Find incoming particles.
1404  int in1, in2;
1405  in1 = in2 = 0;
1406  for ( int i = 0; i < int( event.size()); ++i) {
1407  if ( event[i].mother1() == 1 && in1 == 0) in1 = i;
1408  if ( event[i].mother1() == 2 && in2 == 0) in2 = i;
1409  }
1410 
1411  // Find resonances in hard process.
1412  vector<int> hardResonances;
1413  for ( int i = 0; i < int(event.size()); ++i) {
1414  if ( event[i].status() != -22) continue;
1415  if ( event[i].mother1() != 3) continue;
1416  if ( event[i].mother2() != 4) continue;
1417  hardResonances.push_back(i);
1418  }
1419 
1420  // Find resonances and decay products after showering.
1421  vector<int> evolvedResonances;
1422  vector<pair<int,int> > evolvedDecayProducts;
1423  for ( int j = 0; j < int(hardResonances.size()); ++j) {
1424  for ( int i = int(event.size())-1; i > 0; --i) {
1425  if ( i == hardResonances[j]
1426  || (event[i].mother1() == event[i].mother2()
1427  && event[i].isAncestor(hardResonances[j])) ) {
1428  evolvedResonances.push_back(i);
1429  evolvedDecayProducts.push_back(
1430  make_pair(event[i].daughter1(), event[i].daughter2()) );
1431  break;
1432  }
1433  }
1434  }
1435 
1436  // Event for bookkeeping of resonances.
1437  Event now = Event();
1438  now.init("(dummy event)", particleDataPtr);
1439  now.reset();
1440 
1441  // The PDG id's for the particle entries in this event.
1442  // The status codes for the particle entries in this event.
1443  // Indices for the first and last mother for the particle entries in
1444  // this event.
1445  // The colour-line indices (first(second) is (anti)colour) for the
1446  // particle entries in this event.
1447  // Lab frame momentum (Px, Py, Pz, E and M in GeV) for the particle
1448  // entries in this event.
1449  // Invariant lifetime (c*tau, distance from production to decay in
1450  // mm) for the particle entries in this event.
1451  // Spin info for the particle entries in this event given as the
1452  // cosine of the angle between the spin vector of a particle and the
1453  // 3-momentum of the decaying particle, specified in the lab frame.
1454  hepeup.IDUP.push_back(event[in1].id());
1455  hepeup.IDUP.push_back(event[in2].id());
1456  hepeup.ISTUP.push_back(-1);
1457  hepeup.ISTUP.push_back(-1);
1458  hepeup.MOTHUP.push_back(make_pair(0,0));
1459  hepeup.MOTHUP.push_back(make_pair(0,0));
1460  hepeup.ICOLUP.push_back(make_pair(event[in1].col(),event[in1].acol()));
1461  hepeup.ICOLUP.push_back(make_pair(event[in2].col(),event[in2].acol()));
1462  vector <double> p;
1463  p.push_back(0.0);
1464  p.push_back(0.0);
1465  p.push_back(event[in1].pz());
1466  p.push_back(event[in1].e());
1467  p.push_back(event[in1].m());
1468  hepeup.PUP.push_back(p);
1469  p.resize(0);
1470  p.push_back(0.0);
1471  p.push_back(0.0);
1472  p.push_back(event[in2].pz());
1473  p.push_back(event[in2].e());
1474  p.push_back(event[in2].m());
1475  hepeup.PUP.push_back(p);
1476  p.resize(0);
1477  hepeup.VTIMUP.push_back(event[in1].tau());
1478  hepeup.VTIMUP.push_back(event[in2].tau());
1479  hepeup.SPINUP.push_back(event[in1].pol());
1480  hepeup.SPINUP.push_back(event[in2].pol());
1481 
1482  now.append(event[in1]);
1483  now.append(event[in2]);
1484 
1485  // Attach resonances
1486  for ( int j = 0; j < int(evolvedResonances.size()); ++j) {
1487  int i = evolvedResonances[j];
1488  hepeup.IDUP.push_back(event[i].id());
1489  hepeup.ISTUP.push_back(2);
1490  hepeup.MOTHUP.push_back(make_pair(1,2));
1491  hepeup.ICOLUP.push_back(make_pair(event[i].col(),event[i].acol()));
1492  p.push_back(event[i].px());
1493  p.push_back(event[i].py());
1494  p.push_back(event[i].pz());
1495  p.push_back(event[i].e());
1496  p.push_back(event[i].m());
1497  hepeup.PUP.push_back(p);
1498  p.resize(0);
1499  hepeup.VTIMUP.push_back(event[i].tau());
1500  hepeup.SPINUP.push_back(event[i].pol());
1501  now.append(event[i]);
1502  now.back().statusPos();
1503  }
1504 
1505  // Loop through event and attach remaining decays
1506  vector<int> iSkip;
1507  int iDec = 0;
1508  do {
1509 
1510  if ( now[iDec].isFinal() && now[iDec].canDecay()
1511  && now[iDec].mayDecay() && now[iDec].isResonance() ) {
1512 
1513  int iD1 = now[iDec].daughter1();
1514  int iD2 = now[iDec].daughter2();
1515 
1516  // Done if no daughters exist.
1517  if ( iD1 == 0 || iD2 == 0 ) continue;
1518 
1519  // Attach daughters.
1520  for ( int k = iD1; k <= iD2; ++k ) {
1521  Particle partNow = event[k];
1522  hepeup.IDUP.push_back(partNow.id());
1523  hepeup.MOTHUP.push_back(make_pair(iDec,iDec));
1524  hepeup.ICOLUP.push_back(make_pair(partNow.col(),partNow.acol()));
1525  p.push_back(partNow.px());
1526  p.push_back(partNow.py());
1527  p.push_back(partNow.pz());
1528  p.push_back(partNow.e());
1529  p.push_back(partNow.m());
1530  hepeup.PUP.push_back(p);
1531  p.resize(0);
1532  hepeup.VTIMUP.push_back(partNow.tau());
1533  hepeup.SPINUP.push_back(partNow.pol());
1534  now.append(partNow);
1535  if ( partNow.canDecay() && partNow.mayDecay() && partNow.isResonance()){
1536  now.back().statusPos();
1537  hepeup.ISTUP.push_back(2);
1538  } else
1539  hepeup.ISTUP.push_back(1);
1540 
1541  iSkip.push_back(k);
1542  }
1543 
1544  // End of loop over all entries.
1545  }
1546  } while (++iDec < now.size());
1547 
1548  // Attach final state particles
1549  for ( int i = 0; i < int(event.size()); ++i) {
1550  if (!event[i].isFinal()) continue;
1551  // Skip resonance decay products.
1552  bool skip = false;
1553  for ( int j = 0; j < int(evolvedDecayProducts.size()); ++j) {
1554  skip = ( i >= evolvedDecayProducts[j].first
1555  && i <= evolvedDecayProducts[j].second);
1556  }
1557  if (skip) continue;
1558  for ( int j = 0; j < int(iSkip.size()); ++j) {
1559  skip = ( i == iSkip[j] );
1560  }
1561  if (skip) continue;
1562 
1563  hepeup.IDUP.push_back(event[i].id());
1564  hepeup.ISTUP.push_back(1);
1565  hepeup.MOTHUP.push_back(make_pair(1,2));
1566  hepeup.ICOLUP.push_back(make_pair(event[i].col(),event[i].acol()));
1567  p.push_back(event[i].px());
1568  p.push_back(event[i].py());
1569  p.push_back(event[i].pz());
1570  p.push_back(event[i].e());
1571  p.push_back(event[i].m());
1572  hepeup.PUP.push_back(p);
1573  p.resize(0);
1574  hepeup.VTIMUP.push_back(event[i].tau());
1575  hepeup.SPINUP.push_back(event[i].pol());
1576  now.append(event[i]);
1577  }
1578 
1579  // A pointer to the current HEPRUP object.
1580  hepeup.heprup = &heprup;
1581 
1582  // The weights associated with this event, as given by the LHAwgt tags.
1583  if (infoPtr->weights_detailed)
1584  hepeup.weights_detailed = *(infoPtr->weights_detailed);
1585 
1586  // The weights associated with this event, as given by the LHAweights tags.
1587  if (infoPtr->weights_compressed)
1588  hepeup.weights_compressed = *(infoPtr->weights_compressed);
1589 
1590  // Contents of the LHAscales tag
1591  if (infoPtr->scales)
1592  hepeup.scalesSave = *(infoPtr->scales);
1593 
1594  // Contents of the LHAweights tag (compressed format)
1595  if (infoPtr->weights)
1596  hepeup.weightsSave = *(infoPtr->weights);
1597 
1598  // Contents of the LHArwgt tag (detailed format)
1599  if (infoPtr->rwgt)
1600  hepeup.rwgtSave = *(infoPtr->rwgt);
1601 
1602  // Any other attributes.
1603  if (infoPtr->eventAttributes)
1604  hepeup.attributes = *(infoPtr->eventAttributes);
1605 
1606  // Not implemented yet:
1607  // Write event comments of input LHEF.
1608 
1609  writer.hepeup = hepeup;
1610  if (writeToFile) writer.writeEvent(&hepeup,pDigits);
1611 
1612  return true;
1613 
1614 }
1615 
1616 //--------------------------------------------------------------------------
1617 
1618 // Write end of a Les Houches Event File and close it.
1619 
1620 bool LHEF3FromPythia8::closeLHEF(bool updateInit) {
1621 
1622  // Write an end to the file.
1623  osLHEF << "</LesHouchesEvents>" << endl;
1624  osLHEF.close();
1625 
1626  // Optionally update the cross section information.
1627  if (updateInit) {
1628  const char* cstring = fileName.c_str();
1629  osLHEF.open(cstring, ios::in | ios::out);
1630 
1631  setInit();
1632  osLHEF.close();
1633  }
1634 
1635  // Done.
1636  return true;
1637 
1638 }
1639 
1640 //==========================================================================
1641 
1642 } // end namespace Pythia8
Definition: AgUStep.h:26