StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Settings.cc
1 // Settings.cc is a part of the PYTHIA event generator.
2 // Copyright (C) 2018 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 Settings class.
7 
8 #include "Pythia8/Settings.h"
9 
10 // Allow string and character manipulation.
11 #include <cctype>
12 
13 namespace Pythia8 {
14 
15 //==========================================================================
16 
17 // Settings class.
18 // This class contains flags, modes, parms and words used in generation.
19 
20 //--------------------------------------------------------------------------
21 
22 // Read in database from specific file.
23 
24 bool Settings::init(string startFile, bool append) {
25 
26  // Don't initialize if it has already been done and not in append mode.
27  if (isInit && !append) return true;
28  int nError = 0;
29 
30  // Reset readString history.
31  if (!append) {
32  readStringHistory.resize(0);
33  readStringSubrun.clear();
34  }
35 
36  // List of files to be checked. Start with input file.
37  vector<string> files;
38  files.push_back(startFile);
39 
40  // If nontrivial startfile path, then use that for other files as well.
41  string pathName = "";
42  if (startFile.rfind("/") != string::npos)
43  pathName = startFile.substr(0, startFile.rfind("/") + 1);
44 
45  // Loop over files. Open them for read.
46  for (int i = 0; i < int(files.size()); ++i) {
47  const char* cstring = files[i].c_str();
48  ifstream is(cstring);
49 
50  // Check that instream is OK.
51  if (!is.good()) {
52  cout << "\n PYTHIA Error: settings file " << files[i]
53  << " not found" << endl;
54  return false;
55  }
56 
57  // Read in one line at a time.
58  string line;
59  while ( getline(is, line) ) {
60 
61  // Get first word of a line, to interpret it as tag. Remove "more".
62  istringstream getfirst(line);
63  string tag;
64  getfirst >> tag;
65  if (tag.find("more") != string::npos) tag.erase( tag.find("more"), 4);
66 
67  // Skip ahead if not interesting. Only look for new files in startfile.
68  if (tag != "<flag" && tag != "<flagfix" && tag != "<mode"
69  && tag != "<modeopen" && tag != "<modepick" && tag != "<modefix"
70  && tag != "<parm" && tag != "<parmfix" && tag != "<word"
71  && tag != "<wordfix" && tag != "<fvec" && tag != "<fvecfix"
72  && tag != "<mvec" && tag != "<mvecfix"
73  && tag != "<pvec" && tag != "<pvecfix"
74  && tag != "<wvec" && tag != "<wvecfix" && tag != "<aidx") continue;
75 
76  // Read and append continuation line(s) if line does not contain >.
77  while (line.find(">") == string::npos) {
78  string addLine;
79  getline(is, addLine);
80  line += " " + addLine;
81  }
82 
83  // Remove extra blanks before an = sign.
84  while (line.find(" =") != string::npos) line.erase( line.find(" ="), 1);
85 
86  // Add file also to be read.
87  if (tag == "<aidx") {
88  string name = attributeValue( line, "href");
89  if (name == "") {
90  cout << " PYTHIA Error: failed to find name attribute in line "
91  << line << endl;
92  ++nError;
93  continue;
94  }
95  files.push_back(pathName + name + ".xml");
96  continue;
97  }
98 
99  // Find name attribute.
100  string name = attributeValue( line, "name=");
101  if (name == "") {
102  cout << " PYTHIA Error: failed to find name attribute in line "
103  << line << endl;
104  ++nError;
105  continue;
106  }
107 
108  // Check that default value attribute present, and whether max and min.
109  if (line.find("default=") == string::npos) {
110  cout << " PYTHIA Error: failed to find default value token in line "
111  << line << endl;
112  ++nError;
113  continue;
114  }
115  bool hasMin = (line.find("min=") != string::npos);
116  bool hasMax = (line.find("max=") != string::npos);
117 
118  // Check for occurence of a bool and add to flag map.
119  if (tag == "<flag" || tag == "<flagfix") {
120  bool value = boolAttributeValue( line, "default=");
121  addFlag( name, value);
122 
123  // Check for occurence of an int and add to mode map.
124  } else if (tag == "<mode" || tag == "<modeopen"
125  || tag == "<modepick" || tag == "<modefix") {
126  int value = intAttributeValue( line, "default=");
127  int minVal = intAttributeValue( line, "min=");
128  int maxVal = intAttributeValue( line, "max=");
129  // Enforce check that only allowed options are accepted.
130  bool optOnly = false;
131  if (tag == "<modepick" && hasMin && hasMax) optOnly = true;
132  if (tag == "<modefix") {
133  hasMin = true;
134  hasMax = true;
135  minVal = value;
136  maxVal = value;
137  optOnly = true;
138  }
139  addMode( name, value, hasMin, hasMax, minVal, maxVal, optOnly);
140 
141  // Check for occurence of a double and add to parm map.
142  } else if (tag == "<parm" || tag == "<parmfix") {
143  double value = doubleAttributeValue( line, "default=");
144  double minVal = doubleAttributeValue( line, "min=");
145  double maxVal = doubleAttributeValue( line, "max=");
146  addParm( name, value, hasMin, hasMax, minVal, maxVal);
147 
148  // Check for occurence of a string and add to word map.
149  } else if (tag == "<word" || tag == "<wordfix") {
150  string value = attributeValue( line, "default=");
151  addWord( name, value);
152 
153  // Check for occurence of a bool vector and add to fvec map.
154  } else if (tag == "<fvec" || tag == "<fvecfix") {
155  vector<bool> value = boolVectorAttributeValue( line, "default=");
156  addFVec( name, value);
157 
158  // Check for occurence of an int vector and add to mvec map.
159  } else if (tag == "<mvec" || tag == "<mvecfix") {
160  vector<int> value = intVectorAttributeValue( line, "default=");
161  int minVal = intAttributeValue( line, "min=");
162  int maxVal = intAttributeValue( line, "max=");
163  addMVec( name, value, hasMin, hasMax, minVal, maxVal);
164 
165  // Check for occurence of a double vector and add to pvec map.
166  } else if (tag == "<pvec" || tag == "<pvecfix") {
167  vector<double> value = doubleVectorAttributeValue( line, "default=");
168  double minVal = doubleAttributeValue( line, "min=");
169  double maxVal = doubleAttributeValue( line, "max=");
170  addPVec( name, value, hasMin, hasMax, minVal, maxVal);
171 
172  // Check for occurence of a string vector and add to wvec map.
173  } else if (tag == "<wvec" || tag == "<wvecfix") {
174  vector<string> value = stringVectorAttributeValue( line, "default=");
175  addWVec( name, value);
176  }
177 
178  // End of loop over lines in input file and loop over files.
179  };
180  };
181 
182  // Set up default e+e- and pp tunes, if positive.
183  int eeTune = mode("Tune:ee");
184  if (eeTune > 0) initTuneEE( eeTune);
185  int ppTune = mode("Tune:pp");
186  if (ppTune > 0) initTunePP( ppTune);
187 
188  // Done.
189  if (nError > 0) return false;
190  isInit = true;
191  return true;
192 
193 }
194 
195 
196 //--------------------------------------------------------------------------
197 
198 // Read in database from specific stream.
199 
200 bool Settings::init(istream& is, bool append) {
201 
202  // Don't initialize if it has already been done and not in append mode.
203  if (isInit && !append) return true;
204  int nError = 0;
205 
206  // Check that instream is OK.
207  if (!is.good()) {
208  cout << "\n PYTHIA Error: settings stream not found " << endl;
209  return false;
210  }
211 
212  // Reset readString history.
213  if (!append) {
214  readStringHistory.resize(0);
215  readStringSubrun.clear();
216  }
217 
218  // Read in one line at a time.
219  string line;
220  while ( getline(is, line) ) {
221 
222  // Get first word of a line, to interpret it as tag. Remove "more".
223  istringstream getfirst(line);
224  string tag;
225  getfirst >> tag;
226  if (tag.find("more") != string::npos) tag.erase( tag.find("more"), 4);
227 
228  // Skip ahead if not interesting. Only look for new files in startfile.
229  if (tag != "<flag" && tag != "<flagfix" && tag != "<mode"
230  && tag != "<modeopen" && tag != "<modepick" && tag != "<modefix"
231  && tag != "<parm" && tag != "<parmfix" && tag != "<word"
232  && tag != "<wordfix" && tag != "<fvec" && tag != "<fvecfix"
233  && tag != "<mvec" && tag != "<mvecfix"
234  && tag != "<pvec" && tag != "<pvecfix"
235  && tag != "<wvec" && tag != "<wvecfix" && tag != "<aidx") continue;
236 
237  // Read and append continuation line(s) if line does not contain >.
238  while (line.find(">") == string::npos) {
239  string addLine;
240  getline(is, addLine);
241  line += " " + addLine;
242  }
243 
244  // Remove extra blanks before an = sign.
245  while (line.find(" =") != string::npos) line.erase( line.find(" ="), 1);
246 
247  // Find name attribute.
248  string name = attributeValue( line, "name=");
249  if (name == "") {
250  cout << " PYTHIA Error: failed to find name attribute in line "
251  << line << endl;
252  ++nError;
253  continue;
254  }
255 
256  // Check that default value attribute present, and whether max and min.
257  if (line.find("default=") == string::npos) {
258  cout << " PYTHIA Error: failed to find default value token in line "
259  << line << endl;
260  ++nError;
261  continue;
262  }
263  bool hasMin = (line.find("min=") != string::npos);
264  bool hasMax = (line.find("max=") != string::npos);
265 
266  // Check for occurence of a bool and add to flag map.
267  if (tag == "<flag" || tag == "<flagfix") {
268  bool value = boolAttributeValue( line, "default=");
269  addFlag( name, value);
270 
271  // Check for occurence of an int and add to mode map.
272  } else if (tag == "<mode" || tag == "<modeopen"
273  || tag == "<modepick" || tag == "<modefix") {
274  int value = intAttributeValue( line, "default=");
275  int minVal = intAttributeValue( line, "min=");
276  int maxVal = intAttributeValue( line, "max=");
277 
278  // Enforce check that only allowed options are accepted.
279  bool optOnly = false;
280  if (tag == "<modepick" && hasMin && hasMax) optOnly = true;
281  if (tag == "<modefix") {
282  hasMin = true;
283  hasMax = true;
284  minVal = value;
285  maxVal = value;
286  optOnly = true;
287  }
288  addMode( name, value, hasMin, hasMax, minVal, maxVal, optOnly);
289 
290  // Check for occurence of a double and add to parm map.
291  } else if (tag == "<parm" || tag == "<parmfix") {
292  double value = doubleAttributeValue( line, "default=");
293  double minVal = doubleAttributeValue( line, "min=");
294  double maxVal = doubleAttributeValue( line, "max=");
295  addParm( name, value, hasMin, hasMax, minVal, maxVal);
296 
297  // Check for occurence of a string and add to word map.
298  } else if (tag == "<word" || tag == "<wordfix") {
299  string value = attributeValue( line, "default=");
300  addWord( name, value);
301 
302  // Check for occurence of a bool vector and add to fvec map.
303  } else if (tag == "<fvec" || tag == "<fvecfix") {
304  vector<bool> value = boolVectorAttributeValue( line, "default=");
305  addFVec( name, value);
306 
307  // Check for occurence of an int vector and add to mvec map.
308  } else if (tag == "<mvec" || tag == "<mvecfix") {
309  vector<int> value = intVectorAttributeValue( line, "default=");
310  int minVal = intAttributeValue( line, "min=");
311  int maxVal = intAttributeValue( line, "max=");
312  addMVec( name, value, hasMin, hasMax, minVal, maxVal);
313 
314  // Check for occurence of a double vector and add to pvec map.
315  } else if (tag == "<pvec" || tag == "<pvecfix") {
316  vector<double> value = doubleVectorAttributeValue( line, "default=");
317  double minVal = doubleAttributeValue( line, "min=");
318  double maxVal = doubleAttributeValue( line, "max=");
319  addPVec( name, value, hasMin, hasMax, minVal, maxVal);
320 
321  // Check for occurence of a string vector and add to word map.
322  } else if (tag == "<wvec" || tag == "<wvecfix") {
323  vector<string> value = stringVectorAttributeValue( line, "default=");
324  addWVec( name, value);
325  }
326 
327  // End of loop over lines in input file and loop over files.
328  };
329 
330  // Set up default e+e- and pp tunes, if positive.
331  int eeTune = mode("Tune:ee");
332  if (eeTune > 0) initTuneEE( eeTune);
333  int ppTune = mode("Tune:pp");
334  if (ppTune > 0) initTunePP( ppTune);
335 
336  // Done.
337  if (nError > 0) return false;
338  isInit = true;
339  return true;
340 
341 }
342 
343 //--------------------------------------------------------------------------
344 
345 // Overwrite existing database by reading from specific file.
346 
347 bool Settings::reInit(string startFile) {
348 
349  // Reset maps to empty.
350  flags.clear();
351  modes.clear();
352  parms.clear();
353  words.clear();
354  fvecs.clear();
355  mvecs.clear();
356  pvecs.clear();
357  wvecs.clear();
358 
359  // Then let normal init do the rest.
360  isInit = false;
361  return init(startFile, false);
362 
363 }
364 
365 //--------------------------------------------------------------------------
366 
367 // Read in updates from a character string, like a line of a file.
368 // Is used by readString (and readFile) in Pythia.
369 
370 bool Settings::readString(string line, bool warn) {
371 
372  // If empty line then done.
373  if (line.find_first_not_of(" \n\t\v\b\r\f\a") == string::npos) return true;
374 
375  // If unfinished line then add new to existing, else use input line as is.
376  string lineNow = (lineSaved) ? savedLine + line : line;
377  lineSaved = false;
378 
379  // If first character is not a letter, then taken to be a comment line.
380  int firstChar = lineNow.find_first_not_of(" \n\t\v\b\r\f\a");
381  if (!isalpha(lineNow[firstChar])) return true;
382 
383  // Replace an equal sign by a blank to make parsing simpler, except after {.
384  size_t iBrace = (lineNow.find_first_of("{") == string::npos) ? lineNow.size()
385  : lineNow.find_first_of("{");
386  while (lineNow.find("=") != string::npos
387  && lineNow.find_first_of("=") < iBrace) {
388  int firstEqual = lineNow.find_first_of("=");
389  lineNow.replace(firstEqual, 1, " ");
390  }
391 
392  // Get first word of a line.
393  istringstream splitLine(lineNow);
394  string name;
395  splitLine >> name;
396 
397  // Replace two colons by one (:: -> :) to allow for such mistakes.
398  while (name.find("::") != string::npos) {
399  int firstColonColon = name.find_first_of("::");
400  name.replace(firstColonColon, 2, ":");
401  }
402 
403  // Check whether this is in the database.
404  int inDataBase = 0;
405  if (isFlag(name)) inDataBase = 1;
406  else if (isMode(name)) inDataBase = 2;
407  else if (isParm(name)) inDataBase = 3;
408  else if (isWord(name)) inDataBase = 4;
409  else if (isFVec(name)) inDataBase = 5;
410  else if (isMVec(name)) inDataBase = 6;
411  else if (isPVec(name)) inDataBase = 7;
412  else if (isWVec(name)) inDataBase = 8;
413 
414  // For backwards compatibility: old (parts of) names mapped onto new ones.
415  // This code currently has no use, but is partly preserved for the day
416  // it may be needed again.
417  /*
418  if (inDataBase == 0) {
419  bool retry = false;
420  string nameLower = toLower(name);
421  if (!retry && nameLower.find("minbias") != string::npos) {
422  int firstMB = nameLower.find_first_of("minbias");
423  name.replace(firstMB, 7, "nonDiffractive");
424  retry = true;
425  }
426  if (retry) {
427  if (isFlag(name)) inDataBase = 1;
428  else if (isMode(name)) inDataBase = 2;
429  else if (isParm(name)) inDataBase = 3;
430  else if (isWord(name)) inDataBase = 4;
431  else if (isFVec(name)) inDataBase = 5;
432  else if (isMVec(name)) inDataBase = 6;
433  else if (isPVec(name)) inDataBase = 7;
434  else if (isWVec(name)) inDataBase = 8;
435  }
436  }
437  */
438 
439  // Warn and done if not in database.
440  if (inDataBase == 0) {
441  if (warn) cout << "\n PYTHIA Error: input string not found in settings"
442  << " databases::\n " << line << endl;
443  readingFailedSave = true;
444  return false;
445  }
446 
447  // Find value. Warn if none found.
448  string valueString;
449  splitLine >> valueString;
450  if (!splitLine) {
451  if (warn) cout << "\n PYTHIA Error: variable recognized, but its value"
452  << " not meaningful:\n " << line << endl;
453  readingFailedSave = true;
454  return false;
455  }
456 
457  // If value is a ? then echo the current value.
458  if (valueString == "?") {
459  cout << output(name);
460  return true;
461  }
462 
463  // Check for FORCE= statements (to ignore min/max values)
464  bool force = false;
465  if (valueString.find("force") != string::npos) {
466  force = true;
467  // Read value from next word
468  splitLine >> valueString;
469  if (!splitLine) {
470  if (warn) cout << "\n PYTHIA Error: variable recognized, but its value"
471  << " not meaningful:\n " << line << endl;
472  readingFailedSave = true;
473  return false;
474  }
475  }
476 
477 
478  // If string begins with { then find matching } and extract contents.
479  if (valueString[0] == '{') {
480  size_t openBrace = lineNow.find_first_of("{");
481  size_t closeBrace = lineNow.find_first_of("}");
482  // If not } on same line then must append next line and try again.
483  if (closeBrace == string::npos) {
484  lineSaved = true;
485  savedLine = lineNow;
486  return true;
487  }
488  valueString = lineNow.substr( openBrace + 1, closeBrace - openBrace - 1);
489  }
490 
491  // Update flag map; allow many ways to say yes.
492  if (inDataBase == 1) {
493  bool value = boolString(valueString);
494  flag(name, value, force);
495 
496  // Update mode map.
497  } else if (inDataBase == 2) {
498  istringstream modeData(valueString);
499  int value;
500  modeData >> value;
501  if (!modeData) {
502  if (warn) cout << "\n PYTHIA Error: variable recognized, but its value"
503  << " not meaningful:\n " << line << endl;
504  readingFailedSave = true;
505  return false;
506  }
507  if (!mode(name, value, force)) {
508  if (warn) cout << "\n PYTHIA Error: variable recognized, but its value"
509  << " non-existing option:\n " << line << endl;
510  readingFailedSave = true;
511  return false;
512  }
513 
514  // Update parm map.
515  } else if (inDataBase == 3) {
516  istringstream parmData(valueString);
517  double value;
518  parmData >> value;
519  if (!parmData) {
520  if (warn) cout << "\n PYTHIA Error: variable recognized, but its value"
521  << " not meaningful:\n " << line << endl;
522  readingFailedSave = true;
523  return false;
524  }
525  parm(name, value, force);
526 
527  // Update word map.
528  } else if (inDataBase == 4) {
529  word(name, valueString, force);
530 
531  // Update fvec map.
532  } else if (inDataBase == 5) {
533  istringstream fvecData(valueString);
534  vector<bool> value(boolVectorAttributeValue(
535  "value=\"" + valueString + "\"", "value="));
536  if (!fvecData) {
537  if (warn) cout << "\n PYTHIA Error: variable recognized, but its value"
538  << " not meaningful:\n " << line << endl;
539  readingFailedSave = true;
540  return false;
541  }
542  fvec(name, value, force);
543 
544  // Update mvec map.
545  } else if (inDataBase == 6) {
546  istringstream mvecData(valueString);
547  vector<int> value(intVectorAttributeValue(
548  "value=\"" + valueString + "\"", "value="));
549  if (!mvecData) {
550  if (warn) cout << "\n PYTHIA Error: variable recognized, but its value"
551  << " not meaningful:\n " << line << endl;
552  readingFailedSave = true;
553  return false;
554  }
555  mvec(name, value, force);
556 
557  // Update pvec map.
558  } else if (inDataBase == 7) {
559  istringstream pvecData(valueString);
560  vector<double> value(doubleVectorAttributeValue(
561  "value=\"" + valueString + "\"", "value="));
562  if (!pvecData) {
563  if (warn) cout << "\n PYTHIA Error: variable recognized, but its value"
564  << " not meaningful:\n " << line << endl;
565  readingFailedSave = true;
566  return false;
567  }
568  pvec(name, value, force);
569 
570  // Update wvec map.
571  } else if (inDataBase == 8) {
572  istringstream wvecData(valueString);
573  vector<string> value(stringVectorAttributeValue(
574  "value=\"" + valueString + "\"", "value="));
575  if (!wvecData) {
576  if (warn) cout << "\n PYTHIA Error: variable recognized, but its value"
577  << " not meaningful:\n " << line << endl;
578  readingFailedSave = true;
579  return false;
580  }
581  wvec(name, value, force);
582  }
583 
584  // Store history of valid readString statements
585  readStringHistory.push_back(lineNow);
586  int subrun = max(-1,mode("Main:subrun"));
587  if (readStringSubrun.find(subrun) == readStringSubrun.end())
588  readStringSubrun[subrun] = vector<string>();
589  readStringSubrun[subrun].push_back(lineNow);
590 
591  // Done.
592  return true;
593 }
594 
595 //--------------------------------------------------------------------------
596 
597 // Write updates or everything to user-defined file.
598 
599 bool Settings::writeFile(string toFile, bool writeAll) {
600 
601  // Open file for writing.
602  const char* cstring = toFile.c_str();
603  ofstream os(cstring);
604  if (!os) {
605  infoPtr->errorMsg("Error in Settings::writeFile:"
606  " could not open file", toFile);
607  return false;
608  }
609 
610  // Hand over real work to next method.
611  return writeFile( os, writeAll);
612 
613 }
614 
615 //--------------------------------------------------------------------------
616 
617 // Write updates or everything to user-defined stream (or file).
618 
619 bool Settings::writeFile(ostream& os, bool writeAll) {
620 
621  // Write simple header as comment.
622  if (writeAll) os << "! List of all current PYTHIA ";
623  else os << "! List of all modified PYTHIA ";
624  os << fixed << setprecision(3) << parm("Pythia:versionNumber")
625  << " settings.\n";
626 
627  // Iterators for the flag, mode and parm tables.
628  map<string, Flag>::iterator flagEntry = flags.begin();
629  map<string, Mode>::iterator modeEntry = modes.begin();
630  map<string, Parm>::iterator parmEntry = parms.begin();
631  map<string, Word>::iterator wordEntry = words.begin();
632  map<string, FVec>::iterator fvecEntry = fvecs.begin();
633  map<string, MVec>::iterator mvecEntry = mvecs.begin();
634  map<string, PVec>::iterator pvecEntry = pvecs.begin();
635  map<string, WVec>::iterator wvecEntry = wvecs.begin();
636 
637  // Loop while there is something left to do.
638  while (flagEntry != flags.end() || modeEntry != modes.end()
639  || parmEntry != parms.end() || wordEntry != words.end()
640  || fvecEntry != fvecs.end() || mvecEntry != mvecs.end()
641  || pvecEntry != pvecs.end() || wvecEntry != wvecs.end() ) {
642 
643  // Check if a flag is next in lexigraphical order; if so print it.
644  if ( flagEntry != flags.end()
645  && ( modeEntry == modes.end() || flagEntry->first < modeEntry->first )
646  && ( parmEntry == parms.end() || flagEntry->first < parmEntry->first )
647  && ( wordEntry == words.end() || flagEntry->first < wordEntry->first )
648  && ( fvecEntry == fvecs.end() || flagEntry->first < fvecEntry->first )
649  && ( mvecEntry == mvecs.end() || flagEntry->first < mvecEntry->first )
650  && ( pvecEntry == pvecs.end() || flagEntry->first < pvecEntry->first )
651  && ( wvecEntry == wvecs.end() || flagEntry->first < wvecEntry->first )
652  ) {
653  string state[2] = {"off", "on"};
654  bool valNow = flagEntry->second.valNow;
655  bool valDefault = flagEntry->second.valDefault;
656  if ( writeAll || valNow != valDefault )
657  os << flagEntry->second.name << " = " << state[valNow] << "\n";
658  ++flagEntry;
659 
660  // Else check if mode is next, and if so print it.
661  } else if ( modeEntry != modes.end()
662  && ( parmEntry == parms.end() || modeEntry->first < parmEntry->first )
663  && ( wordEntry == words.end() || modeEntry->first < wordEntry->first )
664  && ( fvecEntry == fvecs.end() || modeEntry->first < fvecEntry->first )
665  && ( mvecEntry == mvecs.end() || modeEntry->first < mvecEntry->first )
666  && ( pvecEntry == pvecs.end() || modeEntry->first < pvecEntry->first )
667  && ( wvecEntry == wvecs.end() || modeEntry->first < wvecEntry->first )
668  ) {
669  int valNow = modeEntry->second.valNow;
670  int valDefault = modeEntry->second.valDefault;
671  if ( writeAll || valNow != valDefault )
672  os << modeEntry->second.name << " = " << valNow << "\n";
673  ++modeEntry;
674 
675  // Else check if parm is next, and if so print it; fixed or scientific.
676  } else if ( parmEntry != parms.end()
677  && ( wordEntry == words.end() || parmEntry->first < wordEntry->first )
678  && ( fvecEntry == fvecs.end() || parmEntry->first < fvecEntry->first )
679  && ( mvecEntry == mvecs.end() || parmEntry->first < mvecEntry->first )
680  && ( pvecEntry == pvecs.end() || parmEntry->first < pvecEntry->first )
681  && ( wvecEntry == wvecs.end() || parmEntry->first < wvecEntry->first )
682  ) {
683  double valNow = parmEntry->second.valNow;
684  double valDefault = parmEntry->second.valDefault;
685  if ( writeAll || valNow != valDefault ) {
686  os << parmEntry->second.name << " = ";
687  if ( valNow == 0. ) os << fixed << setprecision(1);
688  else if ( abs(valNow) < 0.001 ) os << scientific << setprecision(4);
689  else if ( abs(valNow) < 0.1 ) os << fixed << setprecision(7);
690  else if ( abs(valNow) < 1000. ) os << fixed << setprecision(5);
691  else if ( abs(valNow) < 1000000. ) os << fixed << setprecision(3);
692  else os << scientific << setprecision(4);
693  os << valNow << "\n";
694  }
695  ++parmEntry;
696 
697  // Else check if word is next, and if so print it.
698  } else if ( wordEntry != words.end()
699  && ( fvecEntry == fvecs.end() || wordEntry->first < fvecEntry->first )
700  && ( mvecEntry == mvecs.end() || wordEntry->first < mvecEntry->first )
701  && ( pvecEntry == pvecs.end() || wordEntry->first < pvecEntry->first )
702  && ( wvecEntry == wvecs.end() || wordEntry->first < wvecEntry->first )
703  ) {
704  string valNow = wordEntry->second.valNow;
705  string valDefault = wordEntry->second.valDefault;
706  if ( writeAll || valNow != valDefault )
707  os << wordEntry->second.name << " = " << valNow << "\n";
708  ++wordEntry;
709 
710  // Else check if fvec is next, and if so print it.
711  } else if ( fvecEntry != fvecs.end()
712  && ( mvecEntry == mvecs.end() || fvecEntry->first < mvecEntry->first )
713  && ( pvecEntry == pvecs.end() || fvecEntry->first < pvecEntry->first )
714  && ( wvecEntry == wvecs.end() || fvecEntry->first < wvecEntry->first )
715  ) {
716  string state[2] = {"off", "on"};
717  vector<bool> valNow = fvecEntry->second.valNow;
718  vector<bool> valDefault = fvecEntry->second.valDefault;
719  if ( writeAll || valNow != valDefault ) {
720  os << fvecEntry->second.name << " = ";
721  for (vector<bool>::iterator val = valNow.begin();
722  val != --valNow.end(); ++val) os << state[*val] << ",";
723  os << *(--valNow.end()) << "\n";
724  }
725  ++fvecEntry;
726 
727  // Else check if mvec is next, and if so print it.
728  } else if ( mvecEntry != mvecs.end()
729  && ( pvecEntry == pvecs.end() || mvecEntry->first < pvecEntry->first )
730  && ( wvecEntry == wvecs.end() || mvecEntry->first < wvecEntry->first )
731  ) {
732  vector<int> valNow = mvecEntry->second.valNow;
733  vector<int> valDefault = mvecEntry->second.valDefault;
734  if ( writeAll || valNow != valDefault ) {
735  os << mvecEntry->second.name << " = ";
736  for (vector<int>::iterator val = valNow.begin();
737  val != --valNow.end(); ++val) os << *val << ",";
738  os << *(--valNow.end()) << "\n";
739  }
740  ++mvecEntry;
741 
742  // Else check if pvec is next; print fixed or scientific.
743  } else if ( pvecEntry != pvecs.end()
744  && ( wvecEntry == wvecs.end() || pvecEntry->first < wvecEntry->first )
745  ) {
746  vector<double> valNow = pvecEntry->second.valNow;
747  vector<double> valDefault = pvecEntry->second.valDefault;
748  if ( writeAll || valNow != valDefault ) {
749  os << pvecEntry->second.name << " = ";
750  for (vector<double>::iterator val = valNow.begin();
751  val != --valNow.end(); ++val) {
752  if ( *val == 0. ) os << fixed << setprecision(1);
753  else if ( abs(*val) < 0.001 ) os << scientific << setprecision(4);
754  else if ( abs(*val) < 0.1 ) os << fixed << setprecision(7);
755  else if ( abs(*val) < 1000. ) os << fixed << setprecision(5);
756  else if ( abs(*val) < 1000000. ) os << fixed << setprecision(3);
757  else os << scientific << setprecision(4);
758  os << *val << ",";
759  } os << *(--valNow.end()) << "\n";
760  }
761  ++pvecEntry;
762 
763  // Else print wvec.
764  } else {
765  vector<string> valNow = wvecEntry->second.valNow;
766  vector<string> valDefault = wvecEntry->second.valDefault;
767  if ( writeAll || valNow != valDefault ) {
768  os << wvecEntry->second.name << " = ";
769  for (vector<string>::iterator val = valNow.begin();
770  val != --valNow.end(); ++val) os << *val << ",";
771  os << *(--valNow.end()) << "\n";
772  }
773  ++wvecEntry;
774  }
775  } ;
776 
777  // Done.
778  return true;
779 }
780 
781 //--------------------------------------------------------------------------
782 
783  // Write updates or everything to user-defined stream (or file).
784 
785 bool Settings::writeFileXML(ostream& os) {
786 
787  // Iterators for the flag, mode and parm tables.
788  map<string, Flag>::iterator flagEntry = flags.begin();
789  map<string, Mode>::iterator modeEntry = modes.begin();
790  map<string, Parm>::iterator parmEntry = parms.begin();
791  map<string, Word>::iterator wordEntry = words.begin();
792  map<string, FVec>::iterator fvecEntry = fvecs.begin();
793  map<string, MVec>::iterator mvecEntry = mvecs.begin();
794  map<string, PVec>::iterator pvecEntry = pvecs.begin();
795  map<string, WVec>::iterator wvecEntry = wvecs.begin();
796 
797  // Loop while there is something left to do.
798  while (flagEntry != flags.end() || modeEntry != modes.end()
799  || parmEntry != parms.end() || wordEntry != words.end()
800  || fvecEntry != fvecs.end() || mvecEntry != mvecs.end()
801  || pvecEntry != pvecs.end() || wvecEntry != wvecs.end() ) {
802 
803  // Check if a flag is next in lexigraphical order; if so print it.
804  if ( flagEntry != flags.end()
805  && ( modeEntry == modes.end() || flagEntry->first < modeEntry->first )
806  && ( parmEntry == parms.end() || flagEntry->first < parmEntry->first )
807  && ( wordEntry == words.end() || flagEntry->first < wordEntry->first )
808  && ( fvecEntry == fvecs.end() || flagEntry->first < fvecEntry->first )
809  && ( mvecEntry == mvecs.end() || flagEntry->first < mvecEntry->first )
810  && ( pvecEntry == pvecs.end() || flagEntry->first < pvecEntry->first )
811  && ( wvecEntry == wvecs.end() || flagEntry->first < wvecEntry->first )
812  ) {
813  string state[2] = {"off", "on"};
814  bool valDefault = flagEntry->second.valDefault;
815  os << "<flag name=\"" << flagEntry->second.name << "\" default=\""
816  << state[valDefault] << "\"></flag>" << endl;
817  ++flagEntry;
818 
819  // Else check if mode is next, and if so print it.
820  } else if ( modeEntry != modes.end()
821  && ( parmEntry == parms.end() || modeEntry->first < parmEntry->first )
822  && ( wordEntry == words.end() || modeEntry->first < wordEntry->first )
823  && ( fvecEntry == fvecs.end() || modeEntry->first < fvecEntry->first )
824  && ( mvecEntry == mvecs.end() || modeEntry->first < mvecEntry->first )
825  && ( pvecEntry == pvecs.end() || modeEntry->first < pvecEntry->first )
826  && ( wvecEntry == wvecs.end() || modeEntry->first < wvecEntry->first )
827  ) {
828  int valDefault = modeEntry->second.valDefault;
829  os << "<mode name=\"" << modeEntry->second.name << "\" default=\""
830  << valDefault << "\">";
831  if (modeEntry->second.hasMin ) os << " min=\""
832  << modeEntry->second.valMin << "\"";
833  if (modeEntry->second.hasMax ) os << " max=\""
834  << modeEntry->second.valMax << "\"";
835  os << "</mode>" << endl;
836  ++modeEntry;
837 
838  // Else check if parm is next, and if so print it; fixed or scientific.
839  } else if ( parmEntry != parms.end()
840  && ( wordEntry == words.end() || parmEntry->first < wordEntry->first )
841  && ( fvecEntry == fvecs.end() || parmEntry->first < fvecEntry->first )
842  && ( mvecEntry == mvecs.end() || parmEntry->first < mvecEntry->first )
843  && ( pvecEntry == pvecs.end() || parmEntry->first < pvecEntry->first )
844  && ( wvecEntry == wvecs.end() || parmEntry->first < wvecEntry->first )
845  ) {
846  double valDefault = parmEntry->second.valDefault;
847  os << "<parm name=\"" << parmEntry->second.name << "\" default=\"";
848  if ( valDefault == 0. ) os << fixed << setprecision(1);
849  else if ( abs(valDefault) < 0.001 ) os << scientific << setprecision(4);
850  else if ( abs(valDefault) < 0.1 ) os << fixed << setprecision(7);
851  else if ( abs(valDefault) < 1000. ) os << fixed << setprecision(5);
852  else if ( abs(valDefault) < 1000000. ) os << fixed << setprecision(3);
853  else os << scientific << setprecision(4);
854  os << valDefault << "\">";
855  if (parmEntry->second.hasMin) {
856  os << " min=\"";
857  valDefault = parmEntry->second.valMin;
858  if ( valDefault == 0. ) os << fixed << setprecision(1);
859  else if ( abs(valDefault) < 0.001) os << scientific << setprecision(4);
860  else if ( abs(valDefault) < 0.1 ) os << fixed << setprecision(7);
861  else if ( abs(valDefault) < 1000. ) os << fixed << setprecision(5);
862  else if ( abs(valDefault) < 1000000. ) os << fixed << setprecision(3);
863  else os << scientific << setprecision(4);
864  os << valDefault << "\">";
865  }
866  if (parmEntry->second.hasMax) {
867  os << " max=\"";
868  valDefault = parmEntry->second.valMax;
869  if ( valDefault == 0. ) os << fixed << setprecision(1);
870  else if ( abs(valDefault) < 0.001) os << scientific << setprecision(4);
871  else if ( abs(valDefault) < 0.1 ) os << fixed << setprecision(7);
872  else if ( abs(valDefault) < 1000. ) os << fixed << setprecision(5);
873  else if ( abs(valDefault) < 1000000. ) os << fixed << setprecision(3);
874  else os << scientific << setprecision(4);
875  os << valDefault << "\">";
876  }
877  os << "</parm>" << endl;
878  ++parmEntry;
879 
880  // Else check if word is next, and if so print it.
881  } else if ( wordEntry != words.end()
882  && ( fvecEntry == fvecs.end() || wordEntry->first < fvecEntry->first )
883  && ( mvecEntry == mvecs.end() || wordEntry->first < mvecEntry->first )
884  && ( pvecEntry == pvecs.end() || wordEntry->first < pvecEntry->first )
885  && ( wvecEntry == wvecs.end() || wordEntry->first < wvecEntry->first )
886  ) {
887  string valDefault = wordEntry->second.valDefault;
888  os << "<word name=\"" << wordEntry->second.name << "\" default=\""
889  << valDefault << "\"></word>" << endl;
890  ++wordEntry;
891 
892  // Else check if fvec is next, and if so print it.
893  } else if ( fvecEntry != fvecs.end()
894  && ( mvecEntry == mvecs.end() || fvecEntry->first < mvecEntry->first )
895  && ( pvecEntry == pvecs.end() || fvecEntry->first < pvecEntry->first )
896  && ( wvecEntry == wvecs.end() || fvecEntry->first < wvecEntry->first )
897  ) {
898  string state[2] = {"off", "on"};
899  vector<bool> valDefault = fvecEntry->second.valDefault;
900  os << "<fvec name=\"" << fvecEntry->second.name << "\" default=\"";
901  for (vector<bool>::iterator val = valDefault.begin();
902  val != --valDefault.end(); ++val) os << state[*val] << ",";
903  os << state[*(--valDefault.end())] << "\"></fvec>" << endl;
904  ++fvecEntry;
905 
906  // Else check if mvec is next, and if so print it.
907  } else if ( mvecEntry != mvecs.end()
908  && ( pvecEntry == pvecs.end() || mvecEntry->first < pvecEntry->first )
909  && ( wvecEntry == wvecs.end() || mvecEntry->first < wvecEntry->first )
910  ) {
911  vector<int> valDefault = mvecEntry->second.valDefault;
912  os << "<mvec name=\"" << mvecEntry->second.name << "\" default=\"";
913  for (vector<int>::iterator val = valDefault.begin();
914  val != --valDefault.end(); ++val) os << *val << ",";
915  os << *(--valDefault.end()) << "\">";
916  if (mvecEntry->second.hasMin ) os << " min=\""
917  << mvecEntry->second.valMin << "\"";
918  if (mvecEntry->second.hasMax ) os << " max=\""
919  << mvecEntry->second.valMax << "\"";
920  os << "</mvec>" << endl;
921  ++mvecEntry;
922 
923  // Else check if pvec is next; print fixed or scientific.
924  } else if ( pvecEntry != pvecs.end()
925  && ( wvecEntry == wvecs.end() || pvecEntry->first < wvecEntry->first )
926  ) {
927  vector<double> valDefault = pvecEntry->second.valDefault;
928  os << "<pvec name=\"" << pvecEntry->second.name << "\" default=\"";
929  for (vector<double>::iterator val = valDefault.begin();
930  val != --valDefault.end(); ++val) {
931  if ( *val == 0. ) os << fixed << setprecision(1);
932  else if ( abs(*val) < 0.001 ) os << scientific << setprecision(4);
933  else if ( abs(*val) < 0.1 ) os << fixed << setprecision(7);
934  else if ( abs(*val) < 1000. ) os << fixed << setprecision(5);
935  else if ( abs(*val) < 1000000. ) os << fixed << setprecision(3);
936  else os << scientific << setprecision(4);
937  os << *val << ",";
938  }
939  os << *(--valDefault.end()) << "\">";
940  if (pvecEntry->second.hasMin ) {
941  double valLocal = pvecEntry->second.valMin;
942  os << " min=\"";
943  if ( valLocal == 0. ) os << fixed << setprecision(1);
944  else if ( abs(valLocal) < 0.001 ) os << scientific << setprecision(4);
945  else if ( abs(valLocal) < 0.1 ) os << fixed << setprecision(7);
946  else if ( abs(valLocal) < 1000. ) os << fixed << setprecision(5);
947  else if ( abs(valLocal) < 1000000. ) os << fixed << setprecision(3);
948  else os << scientific << setprecision(4);
949  os << valLocal << "\">";
950  }
951  if (pvecEntry->second.hasMax ) {
952  double valLocal = pvecEntry->second.valMax;
953  os << " max=\"";
954  if ( valLocal == 0. ) os << fixed << setprecision(1);
955  else if ( abs(valLocal) < 0.001 ) os << scientific << setprecision(4);
956  else if ( abs(valLocal) < 0.1 ) os << fixed << setprecision(7);
957  else if ( abs(valLocal) < 1000. ) os << fixed << setprecision(5);
958  else if ( abs(valLocal) < 1000000. ) os << fixed << setprecision(3);
959  else os << scientific << setprecision(4);
960  os << valLocal << "\">";
961  }
962  os << "</pvec>" << endl;
963  ++pvecEntry;
964 
965  // Else print wvec.
966  } else {
967  vector<string> valDefault = wvecEntry->second.valDefault;
968  os << "<wvec name=\"" << wvecEntry->second.name << "\" default=\"";
969  for (vector<string>::iterator val = valDefault.begin();
970  val != --valDefault.end(); ++val) os << *val << ",";
971  os << *(--valDefault.end()) << "\">";
972  os << "</wvec>" << endl;
973  ++wvecEntry;
974  }
975  } ;
976 
977  // Done.
978  return true;
979 }
980 
981 //--------------------------------------------------------------------------
982 
983 // Print out table of database in lexigraphical order.
984 
985 void Settings::list(bool doListAll, bool doListString, string match) {
986 
987  // Table header; output for bool as off/on.
988  if (doListAll)
989  cout << "\n *------- PYTHIA Flag + Mode + Parm + Word + FVec + MVec "
990  << "+ PVec + WVec Settings (all) ---------------------------* \n";
991  else if (!doListString)
992  cout << "\n *------- PYTHIA Flag + Mode + Parm + Word + FVec + MVec "
993  << "+ PVec + WVec Settings (changes only) ------------------* \n" ;
994  else
995  cout << "\n *------- PYTHIA Flag + Mode + Parm + Word + FVec + MVec "
996  << "+ PVec + WVec Settings (with requested string) ----------* \n" ;
997  cout << " | "
998  << " | \n"
999  << " | Name | "
1000  << " Now | Default Min Max | \n"
1001  << " | | "
1002  << " | | \n";
1003 
1004  // Convert input string to lowercase for match.
1005  toLowerRep(match);
1006  if (match == "") match = " ";
1007 
1008  // Iterators for the flag, mode and parm tables.
1009  map<string, Flag>::iterator flagEntry = flags.begin();
1010  map<string, Mode>::iterator modeEntry = modes.begin();
1011  map<string, Parm>::iterator parmEntry = parms.begin();
1012  map<string, Word>::iterator wordEntry = words.begin();
1013  map<string, FVec>::iterator fvecEntry = fvecs.begin();
1014  map<string, MVec>::iterator mvecEntry = mvecs.begin();
1015  map<string, PVec>::iterator pvecEntry = pvecs.begin();
1016  map<string, WVec>::iterator wvecEntry = wvecs.begin();
1017 
1018  // Loop while there is something left to do.
1019  while (flagEntry != flags.end() || modeEntry != modes.end()
1020  || parmEntry != parms.end() || wordEntry != words.end()
1021  || fvecEntry != fvecs.end() || mvecEntry != mvecs.end()
1022  || pvecEntry != pvecs.end() || wvecEntry != wvecs.end() ) {
1023 
1024  // Check if a flag is next in lexigraphical order; if so print it.
1025  if ( flagEntry != flags.end()
1026  && ( modeEntry == modes.end() || flagEntry->first < modeEntry->first )
1027  && ( parmEntry == parms.end() || flagEntry->first < parmEntry->first )
1028  && ( wordEntry == words.end() || flagEntry->first < wordEntry->first )
1029  && ( fvecEntry == fvecs.end() || flagEntry->first < fvecEntry->first )
1030  && ( mvecEntry == mvecs.end() || flagEntry->first < mvecEntry->first )
1031  && ( pvecEntry == pvecs.end() || flagEntry->first < pvecEntry->first )
1032  && ( wvecEntry == wvecs.end() || flagEntry->first < wvecEntry->first )
1033  ) {
1034  string state[2] = {"off", "on"};
1035  bool valNow = flagEntry->second.valNow;
1036  bool valDefault = flagEntry->second.valDefault;
1037  if ( doListAll || (!doListString && valNow != valDefault)
1038  || (doListString && flagEntry->first.find(match) != string::npos) )
1039  cout << " | " << setw(45) << left
1040  << flagEntry->second.name << " | " << setw(24) << right
1041  << state[valNow] << " | " << setw(12) << state[valDefault]
1042  << " | \n";
1043  ++flagEntry;
1044 
1045  // Else check if mode is next, and if so print it.
1046  } else if ( modeEntry != modes.end()
1047  && ( parmEntry == parms.end() || modeEntry->first < parmEntry->first )
1048  && ( wordEntry == words.end() || modeEntry->first < wordEntry->first )
1049  && ( fvecEntry == fvecs.end() || modeEntry->first < fvecEntry->first )
1050  && ( mvecEntry == mvecs.end() || modeEntry->first < mvecEntry->first )
1051  && ( pvecEntry == pvecs.end() || modeEntry->first < pvecEntry->first )
1052  && ( wvecEntry == wvecs.end() || modeEntry->first < wvecEntry->first )
1053  ) {
1054  int valNow = modeEntry->second.valNow;
1055  int valDefault = modeEntry->second.valDefault;
1056  if ( doListAll || (!doListString && valNow != valDefault)
1057  || (doListString && modeEntry->first.find(match) != string::npos) ) {
1058  cout << " | " << setw(45) << left
1059  << modeEntry->second.name << " | " << setw(24) << right
1060  << valNow << " | " << setw(12) << valDefault;
1061  if (modeEntry->second.hasMin)
1062  cout << setw(12) << modeEntry->second.valMin;
1063  else cout << " ";
1064  if (modeEntry->second.hasMax)
1065  cout << setw(12) << modeEntry->second.valMax;
1066  else cout << " ";
1067  cout << " | \n";
1068  }
1069  ++modeEntry;
1070 
1071  // Else check if parm is next, and if so print it; fixed or scientific.
1072  } else if ( parmEntry != parms.end()
1073  && ( wordEntry == words.end() || parmEntry->first < wordEntry->first )
1074  && ( fvecEntry == fvecs.end() || parmEntry->first < fvecEntry->first )
1075  && ( mvecEntry == mvecs.end() || parmEntry->first < mvecEntry->first )
1076  && ( pvecEntry == pvecs.end() || parmEntry->first < pvecEntry->first )
1077  && ( wvecEntry == wvecs.end() || parmEntry->first < wvecEntry->first )
1078  ) {
1079  double valNow = parmEntry->second.valNow;
1080  double valDefault = parmEntry->second.valDefault;
1081  if ( doListAll || (!doListString && valNow != valDefault )
1082  || (doListString && parmEntry->first.find(match) != string::npos) ) {
1083  cout << " | " << setw(45) << left
1084  << parmEntry->second.name << right << " | ";
1085  for (int i = 0; i < 4; ++i) {
1086  if (i == 1) valNow = valDefault;
1087  if (i == 2) valNow = parmEntry->second.valMin;
1088  if (i == 3) valNow = parmEntry->second.valMax;
1089  if ( (i == 2 && !parmEntry->second.hasMin)
1090  || (i == 3 && !parmEntry->second.hasMax) )
1091  cout << " ";
1092  else if ( valNow == 0. )
1093  cout << fixed << setprecision(1) << setw(12) << valNow;
1094  else if ( abs(valNow) < 0.001 )
1095  cout << scientific << setprecision(4) << setw(12) << valNow;
1096  else if ( abs(valNow) < 0.1 )
1097  cout << fixed << setprecision(7) << setw(12) << valNow;
1098  else if ( abs(valNow) < 1000. )
1099  cout << fixed << setprecision(5) << setw(12) << valNow;
1100  else if ( abs(valNow) < 1000000. )
1101  cout << fixed << setprecision(3) << setw(12) << valNow;
1102  else
1103  cout << scientific << setprecision(4) << setw(12) << valNow;
1104  if (i == 0) cout << " | ";
1105  }
1106  cout << " | \n";
1107  }
1108  ++parmEntry;
1109 
1110  // Else check if word is next, and if so print it.
1111  } else if ( wordEntry != words.end()
1112  && ( fvecEntry == fvecs.end() || wordEntry->first < fvecEntry->first )
1113  && ( mvecEntry == mvecs.end() || wordEntry->first < mvecEntry->first )
1114  && ( pvecEntry == pvecs.end() || wordEntry->first < pvecEntry->first )
1115  && ( wvecEntry == wvecs.end() || wordEntry->first < wvecEntry->first )
1116  ) {
1117  string valNow = wordEntry->second.valNow;
1118  string valDefault = wordEntry->second.valDefault;
1119  int blankLeft = max(0, 60 - max(24, int(valNow.length()) )
1120  - max(12, int(valDefault.length()) ) );
1121  string blankPad( blankLeft, ' ');
1122  if ( doListAll || (!doListString && valNow != valDefault)
1123  || (doListString && wordEntry->first.find(match) != string::npos) )
1124  cout << " | " << setw(45) << left
1125  << wordEntry->second.name << " | " << setw(24) << right
1126  << valNow << " | " << setw(12) << valDefault << blankPad
1127  << " | \n";
1128  ++wordEntry;
1129 
1130  // Else check if fvec is next, and if so print it.
1131  } else if ( fvecEntry != fvecs.end()
1132  && ( mvecEntry == mvecs.end() || fvecEntry->first < mvecEntry->first )
1133  && ( pvecEntry == pvecs.end() || fvecEntry->first < pvecEntry->first )
1134  && ( wvecEntry == wvecs.end() || fvecEntry->first < wvecEntry->first )
1135  ) {
1136  string state[2] = {"off", "on"};
1137  vector<bool> valsNow = fvecEntry->second.valNow;
1138  vector<bool> valsDefault = fvecEntry->second.valDefault;
1139  bool valNow(false), valDefault(false);
1140  if ( doListAll || (!doListString && valsNow != valsDefault )
1141  || (doListString && fvecEntry->first.find(match) != string::npos) ) {
1142  for (unsigned int i = 0; i < valsNow.size() || i < valsDefault.size();
1143  ++i) {
1144  if ( i == 0 )
1145  cout << " | " << setw(45) << left
1146  << fvecEntry->second.name << right << " | ";
1147  else
1148  cout << " | " << setw(45) << " " << right << " | ";
1149  for (int j = 0; j < 4; ++j) {
1150  if (i < valsNow.size()) valNow = valsNow[i];
1151  if (i < valsDefault.size()) valDefault = valsDefault[i];
1152  if (j == 1) valNow = valDefault;
1153  if ( (j == 0 && i >= valsNow.size())
1154  || (j == 1 && i >= valsDefault.size()) || (j > 1) )
1155  cout << " ";
1156  else cout << setw(12) << state[valNow];
1157  if (j == 0) cout << " | ";
1158  }
1159  cout << " | \n";
1160  }
1161  }
1162  ++fvecEntry;
1163 
1164  // Else check if mvec is next, and if so print it.
1165  } else if ( mvecEntry != mvecs.end()
1166  && ( pvecEntry == pvecs.end() || mvecEntry->first < pvecEntry->first )
1167  && ( wvecEntry == wvecs.end() || mvecEntry->first < wvecEntry->first )
1168  ) {
1169  vector<int> valsNow = mvecEntry->second.valNow;
1170  vector<int> valsDefault = mvecEntry->second.valDefault;
1171  int valNow(0), valDefault(0);
1172  if ( doListAll || (!doListString && valsNow != valsDefault )
1173  || (doListString && mvecEntry->first.find(match) != string::npos) ) {
1174  for (unsigned int i = 0; i < valsNow.size() || i < valsDefault.size();
1175  ++i) {
1176  if ( i == 0 )
1177  cout << " | " << setw(45) << left
1178  << mvecEntry->second.name << right << " | ";
1179  else
1180  cout << " | " << setw(45) << " " << right << " | ";
1181  for (int j = 0; j < 4; ++j) {
1182  if (i < valsNow.size()) valNow = valsNow[i];
1183  if (i < valsDefault.size()) valDefault = valsDefault[i];
1184  if (j == 1) valNow = valDefault;
1185  if (j == 2) valNow = mvecEntry->second.valMin;
1186  if (j == 3) valNow = mvecEntry->second.valMax;
1187  if ( (j == 0 && i >= valsNow.size())
1188  || (j == 1 && i >= valsDefault.size())
1189  || (j == 2 && !mvecEntry->second.hasMin)
1190  || (j == 3 && !mvecEntry->second.hasMax) )
1191  cout << " ";
1192  else cout << setw(12) << valNow;
1193  if (j == 0) cout << " | ";
1194  }
1195  cout << " | \n";
1196  }
1197  }
1198  ++mvecEntry;
1199 
1200  // Else check if pvec is next; print fixed or scientific.
1201  } else if ( pvecEntry != pvecs.end()
1202  && ( wvecEntry == wvecs.end() || pvecEntry->first < wvecEntry->first )
1203  ) {
1204  vector<double> valsNow = pvecEntry->second.valNow;
1205  vector<double> valsDefault = pvecEntry->second.valDefault;
1206  double valNow(0), valDefault(0);
1207  if ( doListAll || (!doListString && valsNow != valsDefault )
1208  || (doListString && pvecEntry->first.find(match) != string::npos) ) {
1209  for (unsigned int i = 0; i < valsNow.size() || i < valsDefault.size();
1210  ++i) {
1211  if ( i == 0 )
1212  cout << " | " << setw(45) << left
1213  << pvecEntry->second.name << right << " | ";
1214  else
1215  cout << " | " << setw(45) << " " << right << " | ";
1216  for (int j = 0; j < 4; ++j) {
1217  if (i < valsNow.size()) valNow = valsNow[i];
1218  if (i < valsDefault.size()) valDefault = valsDefault[i];
1219  if (j == 1) valNow = valDefault;
1220  if (j == 2) valNow = pvecEntry->second.valMin;
1221  if (j == 3) valNow = pvecEntry->second.valMax;
1222  if ( (j == 0 && i >= valsNow.size())
1223  || (j == 1 && i >= valsDefault.size())
1224  || (j == 2 && !pvecEntry->second.hasMin)
1225  || (j == 3 && !pvecEntry->second.hasMax) )
1226  cout << " ";
1227  else if ( valNow == 0. )
1228  cout << fixed << setprecision(1) << setw(12) << valNow;
1229  else if ( abs(valNow) < 0.001 )
1230  cout << scientific << setprecision(4) << setw(12) << valNow;
1231  else if ( abs(valNow) < 0.1 )
1232  cout << fixed << setprecision(7) << setw(12) << valNow;
1233  else if ( abs(valNow) < 1000. )
1234  cout << fixed << setprecision(5) << setw(12) << valNow;
1235  else if ( abs(valNow) < 1000000. )
1236  cout << fixed << setprecision(3) << setw(12) << valNow;
1237  else
1238  cout << scientific << setprecision(4) << setw(12) << valNow;
1239  if (j == 0) cout << " | ";
1240  }
1241  cout << " | \n";
1242  }
1243  }
1244  ++pvecEntry;
1245 
1246  // Else print wvec.
1247  } else {
1248  vector<string> valsNow = wvecEntry->second.valNow;
1249  vector<string> valsDefault = wvecEntry->second.valDefault;
1250  if ( doListAll || (!doListString && valsNow != valsDefault )
1251  || (doListString && wvecEntry->first.find(match) != string::npos) ) {
1252  for (unsigned int i = 0; i < valsNow.size() || i < valsDefault.size();
1253  ++i) {
1254  if ( i == 0 )
1255  cout << " | " << setw(45) << left
1256  << wvecEntry->second.name << right << " | ";
1257  else
1258  cout << " | " << setw(45) << " " << right << " | ";
1259  string valNow = (i < valsNow.size()) ? valsNow[i] : " ";
1260  string valDefault = (i < valsDefault.size()) ? valsDefault[i] : " ";
1261  int blankLeft = max(0, 60 - max(24, int(valNow.length()) )
1262  - max(12, int(valDefault.length()) ) );
1263  string blankPad( blankLeft, ' ');
1264  cout << setw(24) << right << valNow << " | " << setw(12)
1265  << valDefault << blankPad << " | \n";
1266  }
1267  }
1268  ++wvecEntry;
1269  }
1270  } ;
1271 
1272  // End of loop over database contents.
1273  cout << " | "
1274  << " | \n"
1275  << " *------- End PYTHIA Flag + Mode + Parm + Word + FVec + MVec "
1276  << "+ PVec + WVec Settings -----------------------------* " << endl;
1277 
1278 }
1279 
1280 //--------------------------------------------------------------------------
1281 
1282 // Give back current value(s) as a string, whatever the type.
1283 
1284 string Settings::output(string keyIn, bool fullLine) {
1285 
1286  // Default string echoes input key =.
1287  string outVal = (fullLine) ? " " + keyIn + " = " : "";
1288 
1289  // Identify flag, mode, parm or word, and convert to string.
1290  if (isFlag(keyIn)) {
1291  outVal += (flag(keyIn)) ? "true" : "false";
1292  } else if (isMode(keyIn)) {
1293  ostringstream ostr;
1294  ostr << mode(keyIn);
1295  outVal += ostr.str();
1296  } else if (isParm(keyIn)) {
1297  ostringstream ostr;
1298  ostr << scientific << setprecision(5) << parm(keyIn);
1299  outVal += ostr.str();
1300  } else if (isWord(keyIn)) {
1301  outVal += word(keyIn);
1302 
1303  // Identify fvec, mvec, pvec or wvec, and convert to string.
1304  } else if (isFVec(keyIn)) {
1305  vector<bool> outVec = fvec(keyIn);
1306  for (int i = 0; i < int(outVec.size()); ++i) {
1307  outVal += (outVec[i]) ? "true" : "false";
1308  if (i != int(outVec.size()) - 1) outVal += " ";
1309  }
1310  } else if (isMVec(keyIn)) {
1311  vector<int> outVec = mvec(keyIn);
1312  for (int i = 0; i < int(outVec.size()); ++i) {
1313  ostringstream ostr;
1314  ostr << outVec[i];
1315  outVal += ostr.str();
1316  if (i != int(outVec.size()) - 1) outVal += " ";
1317  }
1318  } else if (isPVec(keyIn)) {
1319  vector<double> outVec = pvec(keyIn);
1320  for (int i = 0; i < int(outVec.size()); ++i) {
1321  ostringstream ostr;
1322  ostr << scientific << setprecision(5) << outVec[i];
1323  outVal += ostr.str();
1324  if (i != int(outVec.size()) - 1) outVal += " ";
1325  }
1326  } else if (isWVec(keyIn)) {
1327  vector<string> outVec = wvec(keyIn);
1328  for (int i = 0; i < int(outVec.size()); ++i) {
1329  outVal += outVec[i];
1330  if (i != int(outVec.size()) - 1) outVal += " ";
1331  }
1332 
1333  // Default value, possible endline and done.
1334  } else outVal += "unknown";
1335  if (fullLine) outVal += "\n";
1336  return outVal;
1337 
1338 }
1339 
1340 //--------------------------------------------------------------------------
1341 
1342 // Reset all values to their defaults.
1343 
1344 void Settings::resetAll() {
1345 
1346  // Loop through the flags table, resetting all entries.
1347  for (map<string, Flag>::iterator flagEntry = flags.begin();
1348  flagEntry != flags.end(); ++flagEntry) {
1349  string name = flagEntry->first;
1350  resetFlag(name);
1351  }
1352 
1353  // Loop through the modes table, resetting all entries.
1354  for (map<string, Mode>::iterator modeEntry = modes.begin();
1355  modeEntry != modes.end(); ++modeEntry) {
1356  string name = modeEntry->first;
1357  resetMode(name);
1358  }
1359 
1360  // Loop through the parms table, resetting all entries.
1361  for (map<string, Parm>::iterator parmEntry = parms.begin();
1362  parmEntry != parms.end(); ++parmEntry) {
1363  string name = parmEntry->first;
1364  resetParm(name);
1365  }
1366 
1367  // Loop through the words table, resetting all entries.
1368  for (map<string, Word>::iterator wordEntry = words.begin();
1369  wordEntry != words.end(); ++wordEntry) {
1370  string name = wordEntry->first;
1371  resetWord(name);
1372  }
1373 
1374  // Loop through the fvecs table, resetting all entries.
1375  for (map<string, FVec>::iterator fvecEntry = fvecs.begin();
1376  fvecEntry != fvecs.end(); ++fvecEntry) {
1377  string name = fvecEntry->first;
1378  resetFVec(name);
1379  }
1380 
1381  // Loop through the mvecs table, resetting all entries.
1382  for (map<string, MVec>::iterator mvecEntry = mvecs.begin();
1383  mvecEntry != mvecs.end(); ++mvecEntry) {
1384  string name = mvecEntry->first;
1385  resetMVec(name);
1386  }
1387 
1388  // Loop through the pvecs table, resetting all entries.
1389  for (map<string, PVec>::iterator pvecEntry = pvecs.begin();
1390  pvecEntry != pvecs.end(); ++pvecEntry) {
1391  string name = pvecEntry->first;
1392  resetPVec(name);
1393  }
1394 
1395  // Loop through the wvecs table, resetting all entries.
1396  for (map<string, WVec>::iterator wvecEntry = wvecs.begin();
1397  wvecEntry != wvecs.end(); ++wvecEntry) {
1398  string name = wvecEntry->first;
1399  resetWVec(name);
1400  }
1401 
1402 }
1403 
1404 //--------------------------------------------------------------------------
1405 
1406 // Give back current value, with check that key exists.
1407 
1408 bool Settings::flag(string keyIn) {
1409  if (isFlag(keyIn)) return flags[toLower(keyIn)].valNow;
1410  infoPtr->errorMsg("Error in Settings::flag: unknown key", keyIn);
1411  return false;
1412 }
1413 
1414 int Settings::mode(string keyIn) {
1415  if (isMode(keyIn)) return modes[toLower(keyIn)].valNow;
1416  infoPtr->errorMsg("Error in Settings::mode: unknown key", keyIn);
1417  return 0;
1418 }
1419 
1420 double Settings::parm(string keyIn) {
1421  if (isParm(keyIn)) return parms[toLower(keyIn)].valNow;
1422  infoPtr->errorMsg("Error in Settings::parm: unknown key", keyIn);
1423  return 0.;
1424 }
1425 
1426 string Settings::word(string keyIn) {
1427  if (isWord(keyIn)) return words[toLower(keyIn)].valNow;
1428  infoPtr->errorMsg("Error in Settings::word: unknown key", keyIn);
1429  return " ";
1430 }
1431 
1432 vector<bool> Settings::fvec(string keyIn) {
1433  if (isFVec(keyIn)) return fvecs[toLower(keyIn)].valNow;
1434  infoPtr->errorMsg("Error in Settings::fvec: unknown key", keyIn);
1435  return vector<bool>(1, false);
1436 }
1437 
1438 vector<int> Settings::mvec(string keyIn) {
1439  if (isMVec(keyIn)) return mvecs[toLower(keyIn)].valNow;
1440  infoPtr->errorMsg("Error in Settings::mvec: unknown key", keyIn);
1441  return vector<int>(1, 0);
1442 }
1443 
1444 vector<double> Settings::pvec(string keyIn) {
1445  if (isPVec(keyIn)) return pvecs[toLower(keyIn)].valNow;
1446  infoPtr->errorMsg("Error in Settings::pvec: unknown key", keyIn);
1447  return vector<double>(1, 0.);
1448 }
1449 
1450 vector<string> Settings::wvec(string keyIn) {
1451  if (isWVec(keyIn)) return wvecs[toLower(keyIn)].valNow;
1452  infoPtr->errorMsg("Error in Settings::wvec: unknown key", keyIn);
1453  return vector<string>(1, " ");
1454 }
1455 
1456 //--------------------------------------------------------------------------
1457 
1458 // Give back default value, with check that key exists.
1459 
1460 bool Settings::flagDefault(string keyIn) {
1461  if (isFlag(keyIn)) return flags[toLower(keyIn)].valDefault;
1462  infoPtr->errorMsg("Error in Settings::flagDefault: unknown key", keyIn);
1463  return false;
1464 }
1465 
1466 int Settings::modeDefault(string keyIn) {
1467  if (isMode(keyIn)) return modes[toLower(keyIn)].valDefault;
1468  infoPtr->errorMsg("Error in Settings::modeDefault: unknown key", keyIn);
1469  return 0;
1470 }
1471 
1472 double Settings::parmDefault(string keyIn) {
1473  if (isParm(keyIn)) return parms[toLower(keyIn)].valDefault;
1474  infoPtr->errorMsg("Error in Settings::parmDefault: unknown key", keyIn);
1475  return 0.;
1476 }
1477 
1478 string Settings::wordDefault(string keyIn) {
1479  if (isWord(keyIn)) return words[toLower(keyIn)].valDefault;
1480  infoPtr->errorMsg("Error in Settings::wordDefault: unknown key", keyIn);
1481  return " ";
1482 }
1483 
1484 vector<bool> Settings::fvecDefault(string keyIn) {
1485  if (isFVec(keyIn)) return fvecs[toLower(keyIn)].valDefault;
1486  infoPtr->errorMsg("Error in Settings::fvecDefault: unknown key", keyIn);
1487  return vector<bool>(1, false);
1488 }
1489 
1490 vector<int> Settings::mvecDefault(string keyIn) {
1491  if (isMVec(keyIn)) return mvecs[toLower(keyIn)].valDefault;
1492  infoPtr->errorMsg("Error in Settings::mvecDefault: unknown key", keyIn);
1493  return vector<int>(1, 0);
1494 }
1495 
1496 vector<double> Settings::pvecDefault(string keyIn) {
1497  if (isPVec(keyIn)) return pvecs[toLower(keyIn)].valDefault;
1498  infoPtr->errorMsg("Error in Settings::pvecDefault: unknown key", keyIn);
1499  return vector<double>(1, 0.);
1500 }
1501 
1502 vector<string> Settings::wvecDefault(string keyIn) {
1503  if (isWVec(keyIn)) return wvecs[toLower(keyIn)].valDefault;
1504  infoPtr->errorMsg("Error in Settings::wvecDefault: unknown key", keyIn);
1505  return vector<string>(1, " ");
1506 }
1507 
1508 //--------------------------------------------------------------------------
1509 
1510 // Get a map of entries whose names contain the string "match".
1511 
1512 map<string, Flag> Settings::getFlagMap(string match) {
1513  // Make the match string lower case. Start with an empty map.
1514  toLowerRep(match);
1515  map<string, Flag> flagMap;
1516  // Loop over the flag map (using iterator).
1517  for (map<string,Flag>::iterator flagEntry = flags.begin();
1518  flagEntry != flags.end(); ++flagEntry)
1519  if (flagEntry->first.find(match) != string::npos)
1520  flagMap[flagEntry->first] = flagEntry->second;
1521  return flagMap;
1522 }
1523 
1524 map<string, Mode> Settings::getModeMap(string match) {
1525  // Make the match string lower case. Start with an empty map.
1526  toLowerRep(match);
1527  map<string, Mode> modeMap;
1528  // Loop over the mode map (using iterator).
1529  for (map<string,Mode>::iterator modeEntry = modes.begin();
1530  modeEntry != modes.end(); ++modeEntry)
1531  if (modeEntry->first.find(match) != string::npos)
1532  modeMap[modeEntry->first] = modeEntry->second;
1533  return modeMap;
1534 }
1535 
1536 map<string, Parm> Settings::getParmMap(string match) {
1537  // Make the match string lower case. Start with an empty map.
1538  toLowerRep(match);
1539  map<string, Parm> parmMap;
1540  // Loop over the parm map (using iterator).
1541  for (map<string,Parm>::iterator parmEntry = parms.begin();
1542  parmEntry != parms.end(); ++parmEntry)
1543  if (parmEntry->first.find(match) != string::npos)
1544  parmMap[parmEntry->first] = parmEntry->second;
1545  return parmMap;
1546 }
1547 
1548 map<string, Word> Settings::getWordMap(string match) {
1549  // Make the match string lower case. Start with an empty map.
1550  toLowerRep(match);
1551  map<string, Word> wordMap;
1552  // Loop over the word map (using iterator).
1553  for (map<string,Word>::iterator wordEntry = words.begin();
1554  wordEntry != words.end(); ++wordEntry)
1555  if (wordEntry->first.find(match) != string::npos)
1556  wordMap[wordEntry->first] = wordEntry->second;
1557  return wordMap;
1558 }
1559 
1560 map<string, FVec> Settings::getFVecMap(string match) {
1561  // Make the match string lower case. Start with an empty map.
1562  toLowerRep(match);
1563  map<string, FVec> fvecMap;
1564  // Loop over the fvec map (using iterator).
1565  for (map<string,FVec>::iterator fvecEntry = fvecs.begin();
1566  fvecEntry != fvecs.end(); ++fvecEntry)
1567  if (fvecEntry->first.find(match) != string::npos)
1568  fvecMap[fvecEntry->first] = fvecEntry->second;
1569  return fvecMap;
1570 }
1571 
1572 map<string, MVec> Settings::getMVecMap(string match) {
1573  // Make the match string lower case. Start with an empty map.
1574  toLowerRep(match);
1575  map<string, MVec> mvecMap;
1576  // Loop over the mvec map (using iterator).
1577  for (map<string,MVec>::iterator mvecEntry = mvecs.begin();
1578  mvecEntry != mvecs.end(); ++mvecEntry)
1579  if (mvecEntry->first.find(match) != string::npos)
1580  mvecMap[mvecEntry->first] = mvecEntry->second;
1581  return mvecMap;
1582 }
1583 
1584 map<string, PVec> Settings::getPVecMap(string match) {
1585  // Make the match string lower case. Start with an empty map.
1586  toLowerRep(match);
1587  map<string, PVec> pvecMap;
1588  // Loop over the pvec map (using iterator).
1589  for (map<string,PVec>::iterator pvecEntry = pvecs.begin();
1590  pvecEntry != pvecs.end(); ++pvecEntry)
1591  if (pvecEntry->first.find(match) != string::npos)
1592  pvecMap[pvecEntry->first] = pvecEntry->second;
1593  return pvecMap;
1594 }
1595 
1596 map<string, WVec> Settings::getWVecMap(string match) {
1597  // Make the match string lower case. Start with an empty map.
1598  toLowerRep(match);
1599  map<string, WVec> wvecMap;
1600  // Loop over the wvec map (using iterator).
1601  for (map<string,WVec>::iterator wvecEntry = wvecs.begin();
1602  wvecEntry != wvecs.end(); ++wvecEntry)
1603  if (wvecEntry->first.find(match) != string::npos)
1604  wvecMap[wvecEntry->first] = wvecEntry->second;
1605  return wvecMap;
1606 }
1607 
1608 //--------------------------------------------------------------------------
1609 
1610 // Change current value. Respect limits unless force==true.
1611 // If key not recognised, add new key if force==true, otherwise ignore.
1612 
1613 void Settings::flag(string keyIn, bool nowIn, bool force) {
1614  string keyLower = toLower(keyIn);
1615  if (isFlag(keyIn)) flags[keyLower].valNow = nowIn;
1616  else if (force) addFlag( keyIn, nowIn);
1617  // Print:quiet triggers a whole set of changes.
1618  if (keyLower == "print:quiet") printQuiet( nowIn);
1619 }
1620 
1621 bool Settings::mode(string keyIn, int nowIn, bool force) {
1622  if (isMode(keyIn)) {
1623  string keyLower = toLower(keyIn);
1624  Mode& modeNow = modes[keyLower];
1625  // For modepick and modefix fail if values are outside range.
1626  if (!force && modeNow.optOnly
1627  && (nowIn < modeNow.valMin || nowIn > modeNow.valMax) ) return false;
1628  if (!force && modeNow.hasMin && nowIn < modeNow.valMin)
1629  modeNow.valNow = modeNow.valMin;
1630  else if (!force && modeNow.hasMax && nowIn > modeNow.valMax)
1631  modeNow.valNow = modeNow.valMax;
1632  else modeNow.valNow = nowIn;
1633  // Tune:ee and Tune:pp each trigger a whole set of changes.
1634  if (keyLower == "tune:ee") initTuneEE( modeNow.valNow);
1635  if (keyLower == "tune:pp") initTunePP( modeNow.valNow);
1636  }
1637  else if (force) {
1638  addMode(keyIn, nowIn, false, false, 0, 0);
1639  }
1640  return true;
1641 }
1642 
1643 void Settings::parm(string keyIn, double nowIn, bool force) {
1644  if (isParm(keyIn)) {
1645  Parm& parmNow = parms[toLower(keyIn)];
1646  if (!force && parmNow.hasMin && nowIn < parmNow.valMin)
1647  parmNow.valNow = parmNow.valMin;
1648  else if (!force && parmNow.hasMax && nowIn > parmNow.valMax)
1649  parmNow.valNow = parmNow.valMax;
1650  else parmNow.valNow = nowIn;
1651  }
1652  else if (force) {
1653  addParm(keyIn, nowIn, false, false, 0., 0.);
1654  }
1655 }
1656 
1657 void Settings::word(string keyIn, string nowIn, bool force) {
1658  if (isWord(keyIn)) words[toLower(keyIn)].valNow = nowIn;
1659  else if (force) addWord(keyIn, nowIn);
1660 }
1661 
1662 void Settings::fvec(string keyIn, vector<bool> nowIn, bool force) {
1663  if (isFVec(keyIn)) {
1664  FVec& fvecNow = fvecs[toLower(keyIn)];
1665  fvecNow.valNow.clear();
1666  for (vector<bool>::iterator now = nowIn.begin();
1667  now != nowIn.end(); now++)
1668  fvecNow.valNow.push_back(*now);
1669  }
1670  else if (force) addFVec(keyIn, nowIn);
1671 }
1672 
1673 void Settings::mvec(string keyIn, vector<int> nowIn, bool force) {
1674  if (isMVec(keyIn)) {
1675  MVec& mvecNow = mvecs[toLower(keyIn)];
1676  mvecNow.valNow.clear();
1677  for (vector<int>::iterator now = nowIn.begin();
1678  now != nowIn.end(); now++) {
1679  if (!force && mvecNow.hasMin && *now < mvecNow.valMin)
1680  mvecNow.valNow.push_back(mvecNow.valMin);
1681  else if (!force && mvecNow.hasMax && *now > mvecNow.valMax)
1682  mvecNow.valNow.push_back(mvecNow.valMax);
1683  else mvecNow.valNow.push_back(*now);
1684  }
1685  }
1686  else if (force) addMVec(keyIn, nowIn, false, false, 0, 0);
1687 }
1688 
1689 void Settings::pvec(string keyIn, vector<double> nowIn, bool force) {
1690  if (isPVec(keyIn)) {
1691  PVec& pvecNow = pvecs[toLower(keyIn)];
1692  pvecNow.valNow.clear();
1693  for (vector<double>::iterator now = nowIn.begin();
1694  now != nowIn.end(); now++) {
1695  if (!force && pvecNow.hasMin && *now < pvecNow.valMin)
1696  pvecNow.valNow.push_back(pvecNow.valMin);
1697  else if (!force && pvecNow.hasMax && *now > pvecNow.valMax)
1698  pvecNow.valNow.push_back(pvecNow.valMax);
1699  else pvecNow.valNow.push_back(*now);
1700  }
1701  }
1702  else if (force) addPVec(keyIn, nowIn, false, false, 0., 0.);
1703 }
1704 
1705 void Settings::wvec(string keyIn, vector<string> nowIn, bool force) {
1706  if (isWVec(keyIn)) {
1707  WVec& wvecNow = wvecs[toLower(keyIn)];
1708  wvecNow.valNow.clear();
1709  for (vector<string>::iterator now = nowIn.begin();
1710  now != nowIn.end(); now++)
1711  wvecNow.valNow.push_back(*now);
1712  }
1713  else if (force) addWVec(keyIn, nowIn);
1714 }
1715 
1716 //--------------------------------------------------------------------------
1717 
1718 // Restore current value to default.
1719 
1720 void Settings::resetFlag(string keyIn) {
1721  if (isFlag(keyIn)) flags[toLower(keyIn)].valNow
1722  = flags[toLower(keyIn)].valDefault ;
1723 }
1724 
1725 void Settings::resetMode(string keyIn) {
1726  string keyLower = toLower(keyIn);
1727  if (isMode(keyIn)) modes[keyLower].valNow
1728  = modes[toLower(keyIn)].valDefault ;
1729  // For Tune:ee and Tune:pp must also restore variables involved in tunes.
1730  if (keyLower == "tune:ee") resetTuneEE();
1731  if (keyLower == "tune:pp") resetTunePP();
1732 }
1733 
1734 void Settings::resetParm(string keyIn) {
1735  if (isParm(keyIn)) parms[toLower(keyIn)].valNow
1736  = parms[toLower(keyIn)].valDefault ;
1737 }
1738 
1739 void Settings::resetWord(string keyIn) {
1740  if (isWord(keyIn)) words[toLower(keyIn)].valNow
1741  = words[toLower(keyIn)].valDefault ;
1742 }
1743 
1744 void Settings::resetFVec(string keyIn) {
1745  if (isFVec(keyIn)) fvecs[toLower(keyIn)].valNow
1746  = fvecs[toLower(keyIn)].valDefault ;
1747 }
1748 
1749 void Settings::resetMVec(string keyIn) {
1750  if (isMVec(keyIn)) mvecs[toLower(keyIn)].valNow
1751  = mvecs[toLower(keyIn)].valDefault ;
1752 }
1753 
1754 void Settings::resetPVec(string keyIn) {
1755  if (isPVec(keyIn)) pvecs[toLower(keyIn)].valNow
1756  = pvecs[toLower(keyIn)].valDefault ;
1757 }
1758 
1759 void Settings::resetWVec(string keyIn) {
1760  if (isWVec(keyIn)) wvecs[toLower(keyIn)].valNow
1761  = wvecs[toLower(keyIn)].valDefault ;
1762 }
1763 
1764 //--------------------------------------------------------------------------
1765 
1766 // Regulate level of printout by overall change of settings.
1767 
1768 void Settings::printQuiet(bool quiet) {
1769 
1770  // Switch off as much output as possible.
1771  if (quiet) {
1772  flag("Init:showProcesses", false );
1773  flag("Init:showMultipartonInteractions", false );
1774  flag("Init:showChangedSettings", false );
1775  flag("Init:showAllSettings", false );
1776  flag("Init:showChangedParticleData", false );
1777  flag("Init:showChangedResonanceData", false );
1778  flag("Init:showAllParticleData", false );
1779  mode("Init:showOneParticleData", 0 );
1780  mode("Next:numberCount", 0 );
1781  mode("Next:numberShowLHA", 0 );
1782  mode("Next:numberShowInfo", 0 );
1783  mode("Next:numberShowProcess", 0 );
1784  mode("Next:numberShowEvent", 0 );
1785 
1786  // Restore ouput settings to default.
1787  } else {
1788  resetFlag("Init:showProcesses");
1789  resetFlag("Init:showMultipartonInteractions");
1790  resetFlag("Init:showChangedSettings");
1791  resetFlag("Init:showAllSettings");
1792  resetFlag("Init:showChangedParticleData");
1793  resetFlag("Init:showChangedResonanceData");
1794  resetFlag("Init:showAllParticleData");
1795  resetMode("Init:showOneParticleData");
1796  resetMode("Next:numberCount");
1797  resetMode("Next:numberShowLHA");
1798  resetMode("Next:numberShowInfo");
1799  resetMode("Next:numberShowProcess");
1800  resetMode("Next:numberShowEvent");
1801  }
1802 
1803 }
1804 
1805 //--------------------------------------------------------------------------
1806 
1807 // Restore all e+e- settings to their original values.
1808 
1809 void Settings::resetTuneEE() {
1810 
1811  // Flavour composition.
1812  resetParm("StringFlav:probStoUD");
1813  resetParm("StringFlav:probQQtoQ");
1814  resetParm("StringFlav:probSQtoQQ");
1815  resetParm("StringFlav:probQQ1toQQ0");
1816  resetParm("StringFlav:mesonUDvector");
1817  resetParm("StringFlav:mesonSvector");
1818  resetParm("StringFlav:mesonCvector");
1819  resetParm("StringFlav:mesonBvector");
1820  resetParm("StringFlav:etaSup");
1821  resetParm("StringFlav:etaPrimeSup");
1822  resetParm("StringFlav:popcornSpair");
1823  resetParm("StringFlav:popcornSmeson");
1824  resetFlag("StringFlav:suppressLeadingB");
1825 
1826  // String breaks: z.
1827  resetParm("StringZ:aLund");
1828  resetParm("StringZ:bLund");
1829  resetParm("StringZ:aExtraSquark");
1830  resetParm("StringZ:aExtraDiquark");
1831  resetParm("StringZ:rFactC");
1832  resetParm("StringZ:rFactB");
1833 
1834  // String breaks: pT.
1835  resetParm("StringPT:sigma");
1836  resetParm("StringPT:enhancedFraction");
1837  resetParm("StringPT:enhancedWidth");
1838 
1839  // FSR: strong coupling, IR cutoff.
1840  resetParm("TimeShower:alphaSvalue");
1841  resetMode("TimeShower:alphaSorder");
1842  resetFlag("TimeShower:alphaSuseCMW");
1843  resetParm("TimeShower:pTmin");
1844  resetParm("TimeShower:pTminChgQ");
1845 
1846 }
1847 
1848 //--------------------------------------------------------------------------
1849 
1850 // Restore all pp settings to their original values.
1851 
1852 void Settings::resetTunePP() {
1853 
1854  // PDF set.
1855  resetWord("PDF:pSet");
1856 
1857  // Hard matrix elements alpha_s value.
1858  resetParm("SigmaProcess:alphaSvalue");
1859 
1860  // Diffraction: cross sections and mass distributions.
1861  resetFlag("SigmaTotal:zeroAXB");
1862  resetFlag("SigmaDiffractive:dampen");
1863  resetParm("SigmaDiffractive:maxXB");
1864  resetParm("SigmaDiffractive:maxAX");
1865  resetParm("SigmaDiffractive:maxXX");
1866  resetParm("Diffraction:largeMassSuppress");
1867 
1868  // FSR: dipoles to beam, spin correlations.
1869  resetFlag("TimeShower:dampenBeamRecoil");
1870  resetFlag("TimeShower:phiPolAsym");
1871 
1872  // ISR: strong coupling, IR cutoff, coherence and spin correlations.
1873  resetParm("SpaceShower:alphaSvalue");
1874  resetMode("SpaceShower:alphaSorder");
1875  resetParm("SpaceShower:alphaSuseCMW");
1876  resetFlag("SpaceShower:samePTasMPI");
1877  resetParm("SpaceShower:pT0Ref");
1878  resetParm("SpaceShower:ecmRef");
1879  resetParm("SpaceShower:ecmPow");
1880  resetParm("SpaceShower:pTmaxFudge");
1881  resetParm("SpaceShower:pTdampFudge");
1882  resetFlag("SpaceShower:rapidityOrder");
1883  resetFlag("SpaceShower:rapidityOrderMPI");
1884  resetFlag("SpaceShower:phiPolAsym");
1885  resetFlag("SpaceShower:phiIntAsym");
1886 
1887  // MPI: strong coupling, IR regularization, energy scaling.
1888  resetParm("MultipartonInteractions:alphaSvalue");
1889  resetParm("MultipartonInteractions:pT0Ref");
1890  resetParm("MultipartonInteractions:ecmRef");
1891  resetParm("MultipartonInteractions:ecmPow");
1892  resetMode("MultipartonInteractions:bProfile");
1893  resetParm("MultipartonInteractions:expPow");
1894  resetParm("MultipartonInteractions:a1");
1895 
1896  // Beam remnant parameters.
1897  resetParm("BeamRemnants:primordialKTsoft");
1898  resetParm("BeamRemnants:primordialKThard");
1899  resetParm("BeamRemnants:halfScaleForKT");
1900  resetParm("BeamRemnants:halfMassForKT");
1901 
1902  // Colour reconnection parameters.
1903  resetMode("ColourReconnection:mode");
1904  resetParm("ColourReconnection:range");
1905 
1906 }
1907 
1908 //--------------------------------------------------------------------------
1909 
1910 // Set the values related to a tune of e+e- data,
1911 // i.e. mainly for final-state radiation and hadronization.
1912 
1913 void Settings::initTuneEE( int eeTune) {
1914 
1915  // Restore all e+e- settings to their original values.
1916  // Is first step for setting up a specific tune.
1917  if (eeTune != 0) resetTuneEE();
1918 
1919  // Old flavour and FSR defaults carried over from very old JETSET tune,
1920  // only with alphaS roughly tuned for "new" pT-ordered shower.
1921  if (eeTune == 1) {
1922  parm("StringFlav:probStoUD", 0.30 );
1923  parm("StringFlav:probQQtoQ", 0.10 );
1924  parm("StringFlav:probSQtoQQ", 0.40 );
1925  parm("StringFlav:probQQ1toQQ0", 0.05 );
1926  parm("StringFlav:mesonUDvector", 1.00 );
1927  parm("StringFlav:mesonSvector", 1.50 );
1928  parm("StringFlav:mesonCvector", 2.50 );
1929  parm("StringFlav:mesonBvector", 3.00 );
1930  parm("StringFlav:etaSup", 1.00 );
1931  parm("StringFlav:etaPrimeSup", 0.40 );
1932  parm("StringFlav:popcornSpair", 0.50 );
1933  parm("StringFlav:popcornSmeson", 0.50 );
1934  flag("StringFlav:suppressLeadingB", false );
1935  parm("StringZ:aLund", 0.30 );
1936  parm("StringZ:bLund", 0.58 );
1937  parm("StringZ:aExtraSquark", 0.00 );
1938  parm("StringZ:aExtraDiquark", 0.50 );
1939  parm("StringZ:rFactC", 1.00 );
1940  parm("StringZ:rFactB", 1.00 );
1941  parm("StringPT:sigma", 0.36 );
1942  parm("StringPT:enhancedFraction", 0.01 );
1943  parm("StringPT:enhancedWidth", 2.0 );
1944  parm("TimeShower:alphaSvalue", 0.137 );
1945  mode("TimeShower:alphaSorder", 1 );
1946  flag("TimeShower:alphaSuseCMW", false );
1947  parm("TimeShower:pTmin", 0.5 );
1948  parm("TimeShower:pTminChgQ", 0.5 );
1949  }
1950 
1951  // Marc Montull's tune to particle composition at LEP1 (August 2007).
1952  else if (eeTune == 2) {
1953  parm("StringFlav:probStoUD", 0.22 );
1954  parm("StringFlav:probQQtoQ", 0.08 );
1955  parm("StringFlav:probSQtoQQ", 0.75 );
1956  parm("StringFlav:probQQ1toQQ0", 0.025 );
1957  parm("StringFlav:mesonUDvector", 0.5 );
1958  parm("StringFlav:mesonSvector", 0.6 );
1959  parm("StringFlav:mesonCvector", 1.5 );
1960  parm("StringFlav:mesonBvector", 2.5 );
1961  parm("StringFlav:etaSup", 0.60 );
1962  parm("StringFlav:etaPrimeSup", 0.15 );
1963  parm("StringFlav:popcornSpair", 1.0 );
1964  parm("StringFlav:popcornSmeson", 1.0 );
1965  flag("StringFlav:suppressLeadingB", false ); // kept fixed
1966  parm("StringZ:aLund", 0.76 );
1967  parm("StringZ:bLund", 0.58 ); // kept fixed
1968  parm("StringZ:aExtraSquark", 0.00 ); // kept fixed
1969  parm("StringZ:aExtraDiquark", 0.50 ); // kept fixed
1970  parm("StringZ:rFactC", 1.00 ); // kept fixed
1971  parm("StringZ:rFactB", 1.00 ); // kept fixed
1972  parm("StringPT:sigma", 0.36 ); // kept fixed
1973  parm("StringPT:enhancedFraction", 0.01 ); // kept fixed
1974  parm("StringPT:enhancedWidth", 2.0 ); // kept fixed
1975  parm("TimeShower:alphaSvalue", 0.137 ); // kept fixed
1976  mode("TimeShower:alphaSorder", 1 ); // kept fixed
1977  flag("TimeShower:alphaSuseCMW", false ); // kept fixed
1978  parm("TimeShower:pTmin", 0.5 ); // kept fixed
1979  parm("TimeShower:pTminChgQ", 0.5 ); // kept fixed
1980  }
1981 
1982  // Full e+e- tune of flavours and FSR to LEP1 data within the
1983  // Rivet + Professor framework, by Hendrik Hoeth (June 2009).
1984  else if (eeTune == 3) {
1985  parm("StringFlav:probStoUD", 0.19 );
1986  parm("StringFlav:probQQtoQ", 0.09 );
1987  parm("StringFlav:probSQtoQQ", 1.00 );
1988  parm("StringFlav:probQQ1toQQ0", 0.027 );
1989  parm("StringFlav:mesonUDvector", 0.62 );
1990  parm("StringFlav:mesonSvector", 0.725 );
1991  parm("StringFlav:mesonCvector", 1.06 );
1992  parm("StringFlav:mesonBvector", 3.0 );
1993  parm("StringFlav:etaSup", 0.63 );
1994  parm("StringFlav:etaPrimeSup", 0.12 );
1995  parm("StringFlav:popcornSpair", 0.5 ); // kept fixed
1996  parm("StringFlav:popcornSmeson", 0.5 ); // kept fixed
1997  flag("StringFlav:suppressLeadingB", false ); // kept fixed
1998  parm("StringZ:aLund", 0.3 ); // kept fixed
1999  parm("StringZ:bLund", 0.8 );
2000  parm("StringZ:aExtraSquark", 0.00 ); // kept fixed
2001  parm("StringZ:aExtraDiquark", 0.50 ); // kept fixed
2002  parm("StringZ:rFactC", 1.00 ); // kept fixed
2003  parm("StringZ:rFactB", 0.67 );
2004  parm("StringPT:sigma", 0.304 );
2005  parm("StringPT:enhancedFraction", 0.01 ); // kept fixed
2006  parm("StringPT:enhancedWidth", 2.0 ); // kept fixed
2007  parm("TimeShower:alphaSvalue", 0.1383);
2008  mode("TimeShower:alphaSorder", 1 ); // kept fixed
2009  flag("TimeShower:alphaSuseCMW", false ); // kept fixed
2010  parm("TimeShower:pTmin", 0.4 ); // kept fixed (near limit)
2011  parm("TimeShower:pTminChgQ", 0.4 ); // kept same as pTmin
2012  }
2013 
2014  // Full e+e- tune of flavours and FSR to LEP1 data, by Peter Skands
2015  // (September 2013). Note use of CMW convention for shower.
2016  else if (eeTune == 4) {
2017  parm("StringFlav:probStoUD", 0.21 );
2018  parm("StringFlav:probQQtoQ", 0.086 );
2019  parm("StringFlav:probSQtoQQ", 1.00 );
2020  parm("StringFlav:probQQ1toQQ0", 0.031 );
2021  parm("StringFlav:mesonUDvector", 0.45 );
2022  parm("StringFlav:mesonSvector", 0.60 );
2023  parm("StringFlav:mesonCvector", 0.95 );
2024  parm("StringFlav:mesonBvector", 3.0 ); // kept fixed
2025  parm("StringFlav:etaSup", 0.65 );
2026  parm("StringFlav:etaPrimeSup", 0.08 );
2027  parm("StringFlav:popcornSpair", 0.5 ); // kept fixed
2028  parm("StringFlav:popcornSmeson", 0.5 ); // kept fixed
2029  flag("StringFlav:suppressLeadingB", false ); // kept fixed
2030  parm("StringZ:aLund", 0.55 );
2031  parm("StringZ:bLund", 1.08 );
2032  parm("StringZ:aExtraSquark", 0.00 ); // kept fixed
2033  parm("StringZ:aExtraDiquark", 1.00 );
2034  parm("StringZ:rFactC", 1.00 ); // kept fixed
2035  parm("StringZ:rFactB", 0.85 );
2036  parm("StringPT:sigma", 0.305 );
2037  parm("StringPT:enhancedFraction", 0.01 ); // kept fixed
2038  parm("StringPT:enhancedWidth", 2.0 ); // kept fixed
2039  parm("TimeShower:alphaSvalue", 0.127 );
2040  mode("TimeShower:alphaSorder", 1 ); // kept fixed
2041  flag("TimeShower:alphaSuseCMW", true );
2042  parm("TimeShower:pTmin", 0.4 );
2043  parm("TimeShower:pTminChgQ", 0.4 ); // kept same as pTmin
2044  }
2045 
2046  // First e+e- tune by Nadine Fischer, using eeTune = 3 for flavour
2047  // composition (September 2013).
2048  else if (eeTune == 5) {
2049  parm("StringFlav:probStoUD", 0.19 ); // kept fixed
2050  parm("StringFlav:probQQtoQ", 0.09 ); // kept fixed
2051  parm("StringFlav:probSQtoQQ", 1.00 ); // kept fixed
2052  parm("StringFlav:probQQ1toQQ0", 0.027 ); // kept fixed
2053  parm("StringFlav:mesonUDvector", 0.62 ); // kept fixed
2054  parm("StringFlav:mesonSvector", 0.725 ); // kept fixed
2055  parm("StringFlav:mesonCvector", 1.06 ); // kept fixed
2056  parm("StringFlav:mesonBvector", 3.0 ); // kept fixed
2057  parm("StringFlav:etaSup", 0.63 ); // kept fixed
2058  parm("StringFlav:etaPrimeSup", 0.12 ); // kept fixed
2059  parm("StringFlav:popcornSpair", 0.5 ); // kept fixed
2060  parm("StringFlav:popcornSmeson", 0.5 ); // kept fixed
2061  flag("StringFlav:suppressLeadingB", false ); // kept fixed
2062  parm("StringZ:aLund", 0.386 );
2063  parm("StringZ:bLund", 0.977 );
2064  parm("StringZ:aExtraSquark", 0.00 ); // kept fixed
2065  parm("StringZ:aExtraDiquark", 0.940 );
2066  parm("StringZ:rFactC", 1.00 ); // kept fixed
2067  parm("StringZ:rFactB", 0.67 ); // kept fixed
2068  parm("StringPT:sigma", 0.286 );
2069  parm("StringPT:enhancedFraction", 0.01 ); // kept fixed
2070  parm("StringPT:enhancedWidth", 2.0 ); // kept fixed
2071  parm("TimeShower:alphaSvalue", 0.139 );
2072  mode("TimeShower:alphaSorder", 1 ); // kept fixed
2073  flag("TimeShower:alphaSuseCMW", false ); // kept fixed
2074  parm("TimeShower:pTmin", 0.409 );
2075  parm("TimeShower:pTminChgQ", 0.409 ); // kept same as pTmin
2076  }
2077 
2078  // Second e+e- tune by Nadine Fischer, using eeTune = 3 for flavour
2079  // composition (September 2013).
2080  else if (eeTune == 6) {
2081  parm("StringFlav:probStoUD", 0.19 ); // kept fixed
2082  parm("StringFlav:probQQtoQ", 0.09 ); // kept fixed
2083  parm("StringFlav:probSQtoQQ", 1.00 ); // kept fixed
2084  parm("StringFlav:probQQ1toQQ0", 0.027 ); // kept fixed
2085  parm("StringFlav:mesonUDvector", 0.62 ); // kept fixed
2086  parm("StringFlav:mesonSvector", 0.725 ); // kept fixed
2087  parm("StringFlav:mesonCvector", 1.06 ); // kept fixed
2088  parm("StringFlav:mesonBvector", 3.0 ); // kept fixed
2089  parm("StringFlav:etaSup", 0.63 ); // kept fixed
2090  parm("StringFlav:etaPrimeSup", 0.12 ); // kept fixed
2091  parm("StringFlav:popcornSpair", 0.5 ); // kept fixed
2092  parm("StringFlav:popcornSmeson", 0.5 ); // kept fixed
2093  flag("StringFlav:suppressLeadingB", false ); // kept fixed
2094  parm("StringZ:aLund", 0.351 );
2095  parm("StringZ:bLund", 0.942 );
2096  parm("StringZ:aExtraSquark", 0.00 ); // kept fixed
2097  parm("StringZ:aExtraDiquark", 0.547 );
2098  parm("StringZ:rFactC", 1.00 ); // kept fixed
2099  parm("StringZ:rFactB", 0.67 ); // kept fixed
2100  parm("StringPT:sigma", 0.283 );
2101  parm("StringPT:enhancedFraction", 0.01 ); // kept fixed
2102  parm("StringPT:enhancedWidth", 2.0 ); // kept fixed
2103  parm("TimeShower:alphaSvalue", 0.139);
2104  mode("TimeShower:alphaSorder", 1 ); // kept fixed
2105  flag("TimeShower:alphaSuseCMW", false ); // kept fixed
2106  parm("TimeShower:pTmin", 0.406 );
2107  parm("TimeShower:pTminChgQ", 0.406 ); // kept same as pTmin
2108  }
2109 
2110  // The Monash 2013 tune by Peter Skands, the e+e- part (January 2014).
2111  else if (eeTune == 7) {
2112  parm("StringFlav:probStoUD", 0.217 );
2113  parm("StringFlav:probQQtoQ", 0.081 );
2114  parm("StringFlav:probSQtoQQ", 0.915 );
2115  parm("StringFlav:probQQ1toQQ0", 0.0275);
2116  parm("StringFlav:mesonUDvector", 0.50 );
2117  parm("StringFlav:mesonSvector", 0.55 );
2118  parm("StringFlav:mesonCvector", 0.88 );
2119  parm("StringFlav:mesonBvector", 2.20 );
2120  parm("StringFlav:etaSup", 0.60 );
2121  parm("StringFlav:etaPrimeSup", 0.12 );
2122  parm("StringFlav:popcornSpair", 0.90 );
2123  parm("StringFlav:popcornSmeson", 0.50 );
2124  flag("StringFlav:suppressLeadingB", false ); // kept fixed
2125  parm("StringZ:aLund", 0.68 );
2126  parm("StringZ:bLund", 0.98 );
2127  parm("StringZ:aExtraSquark", 0.00 ); // kept fixed
2128  parm("StringZ:aExtraDiquark", 0.97 );
2129  parm("StringZ:rFactC", 1.32 );
2130  parm("StringZ:rFactB", 0.855 );
2131  parm("StringPT:sigma", 0.335 );
2132  parm("StringPT:enhancedFraction", 0.01 ); // kept fixed
2133  parm("StringPT:enhancedWidth", 2.0 ); // kept fixed
2134  parm("TimeShower:alphaSvalue", 0.1365);
2135  mode("TimeShower:alphaSorder", 1 ); // kept fixed
2136  flag("TimeShower:alphaSuseCMW", false ); // kept fixed
2137  parm("TimeShower:pTmin", 0.5 ); // kept fixed
2138  parm("TimeShower:pTminChgQ", 0.5 ); // kept fixed
2139  }
2140 
2141 }
2142 
2143 //--------------------------------------------------------------------------
2144 
2145 // Set the values related to a tune of pp/ppbar data,
2146 // i.e. mainly for initial-state radiation and multiparton interactions.
2147 
2148 void Settings::initTunePP( int ppTune) {
2149 
2150  // Restore all pp/ppbar settings to their original values.
2151  // Is first step for setting up a specific tune.
2152  if (ppTune != 0) resetTunePP();
2153 
2154  // Set up e+e- tune that goes with the corresponding pp tune.
2155  if (ppTune > 0) {
2156  int eeTune = 3;
2157  if (ppTune == 14 || ppTune >= 18) eeTune = 7;
2158  // The mode setting is for documentation, the real action is by initTuneEE.
2159  mode("Tune:ee", eeTune );
2160  initTuneEE( eeTune);
2161  }
2162 
2163  // Decide whether to use LHAPFD where possible.
2164  int preferLHAPDF = mode("Tune:preferLHAPDF");
2165 
2166  // Old ISR and MPI defaults from early and primitive comparisons with data.
2167  if (ppTune == 1) {
2168  word("PDF:pSet", "2" );
2169  parm("SigmaProcess:alphaSvalue", 0.1265);
2170  flag("SigmaTotal:zeroAXB", true );
2171  flag("SigmaDiffractive:dampen", false );
2172  parm("Diffraction:largeMassSuppress", 2.0 );
2173  flag("TimeShower:dampenBeamRecoil", false );
2174  flag("TimeShower:phiPolAsym", false );
2175  parm("SpaceShower:alphaSvalue", 0.127 );
2176  mode("SpaceShower:alphaSorder", 1 );
2177  flag("SpaceShower:alphaSuseCMW", false );
2178  flag("SpaceShower:samePTasMPI", true );
2179  parm("SpaceShower:pT0Ref", 2.2 );
2180  parm("SpaceShower:ecmRef", 1800.0);
2181  parm("SpaceShower:ecmPow", 0.16 );
2182  parm("SpaceShower:pTmaxFudge", 1.0 );
2183  parm("SpaceShower:pTdampFudge", 1.0 );
2184  flag("SpaceShower:rapidityOrder", false );
2185  flag("SpaceShower:rapidityOrderMPI", false );
2186  flag("SpaceShower:phiPolAsym", false );
2187  flag("SpaceShower:phiIntAsym", false );
2188  parm("MultipartonInteractions:alphaSvalue", 0.127 );
2189  parm("MultipartonInteractions:pT0Ref", 2.15 );
2190  parm("MultipartonInteractions:ecmRef", 1800. );
2191  parm("MultipartonInteractions:ecmPow", 0.16 );
2192  mode("MultipartonInteractions:bProfile", 2 );
2193  parm("MultipartonInteractions:expPow", 1.0 );
2194  parm("MultipartonInteractions:a1", 0.15 );
2195  parm("BeamRemnants:primordialKTsoft", 0.4 );
2196  parm("BeamRemnants:primordialKThard", 2.1 );
2197  parm("BeamRemnants:halfScaleForKT", 7.0 );
2198  parm("BeamRemnants:halfMassForKT", 2.0 );
2199  mode("ColourReconnection:mode", 0 );
2200  parm("ColourReconnection:range", 2.5 );
2201  }
2202 
2203  // "Tune 1" simple first tune by Peter Skands to ISR and MPI, July 2009.
2204  else if (ppTune == 2) {
2205  word("PDF:pSet", "2" );
2206  parm("SigmaProcess:alphaSvalue", 0.1265);
2207  flag("SigmaTotal:zeroAXB", true );
2208  flag("SigmaDiffractive:dampen", false );
2209  parm("Diffraction:largeMassSuppress", 2.0 );
2210  flag("TimeShower:dampenBeamRecoil", false );
2211  flag("TimeShower:phiPolAsym", false );
2212  parm("SpaceShower:alphaSvalue", 0.137 );
2213  mode("SpaceShower:alphaSorder", 1 );
2214  flag("SpaceShower:alphaSuseCMW", false );
2215  flag("SpaceShower:samePTasMPI", false );
2216  parm("SpaceShower:pT0Ref", 2.0 );
2217  parm("SpaceShower:ecmRef", 1800.0);
2218  parm("SpaceShower:ecmPow", 0.0 );
2219  parm("SpaceShower:pTmaxFudge", 1.0 );
2220  parm("SpaceShower:pTdampFudge", 1.0 );
2221  flag("SpaceShower:rapidityOrder", false );
2222  flag("SpaceShower:rapidityOrderMPI", false );
2223  flag("SpaceShower:phiPolAsym", false );
2224  flag("SpaceShower:phiIntAsym", false );
2225  parm("MultipartonInteractions:alphaSvalue", 0.127 );
2226  parm("MultipartonInteractions:pT0Ref", 2.25 );
2227  parm("MultipartonInteractions:ecmRef", 1800. );
2228  parm("MultipartonInteractions:ecmPow", 0.24 );
2229  mode("MultipartonInteractions:bProfile", 1 );
2230  parm("MultipartonInteractions:expPow", 1.0 );
2231  parm("MultipartonInteractions:a1", 0.15 );
2232  parm("BeamRemnants:primordialKTsoft", 0.5 );
2233  parm("BeamRemnants:primordialKThard", 2.0 );
2234  parm("BeamRemnants:halfScaleForKT", 1.0 );
2235  parm("BeamRemnants:halfMassForKT", 1.0 );
2236  mode("ColourReconnection:mode", 0 );
2237  parm("ColourReconnection:range", 10.0 );
2238  }
2239 
2240  // Tune 2C, July 2010.
2241  else if (ppTune == 3) {
2242  word("PDF:pSet", "8" );
2243  parm("SigmaProcess:alphaSvalue", 0.135 );
2244  flag("SigmaTotal:zeroAXB", true );
2245  flag("SigmaDiffractive:dampen", false );
2246  parm("Diffraction:largeMassSuppress", 2.0 );
2247  flag("TimeShower:dampenBeamRecoil", true );
2248  flag("TimeShower:phiPolAsym", true );
2249  parm("SpaceShower:alphaSvalue", 0.137 );
2250  mode("SpaceShower:alphaSorder", 1 );
2251  flag("SpaceShower:alphaSuseCMW", false );
2252  flag("SpaceShower:samePTasMPI", false );
2253  parm("SpaceShower:pT0Ref", 2.0 );
2254  parm("SpaceShower:ecmRef", 1800.0);
2255  parm("SpaceShower:ecmPow", 0.0 );
2256  parm("SpaceShower:pTmaxFudge", 1.0 );
2257  parm("SpaceShower:pTdampFudge", 1.0 );
2258  flag("SpaceShower:rapidityOrder", true );
2259  flag("SpaceShower:rapidityOrderMPI", true );
2260  flag("SpaceShower:phiPolAsym", true );
2261  flag("SpaceShower:phiIntAsym", true );
2262  parm("MultipartonInteractions:alphaSvalue", 0.135 );
2263  parm("MultipartonInteractions:pT0Ref", 2.32 );
2264  parm("MultipartonInteractions:ecmRef", 1800. );
2265  parm("MultipartonInteractions:ecmPow", 0.21 );
2266  mode("MultipartonInteractions:bProfile", 3 );
2267  parm("MultipartonInteractions:expPow", 1.6 );
2268  parm("MultipartonInteractions:a1", 0.15 );
2269  parm("BeamRemnants:primordialKTsoft", 0.5 );
2270  parm("BeamRemnants:primordialKThard", 2.0 );
2271  parm("BeamRemnants:halfScaleForKT", 1.0 );
2272  parm("BeamRemnants:halfMassForKT", 1.0 );
2273  mode("ColourReconnection:mode", 0 );
2274  parm("ColourReconnection:range", 3.0 );
2275  }
2276 
2277  // Tune 2M, July 2010.
2278  else if (ppTune == 4) {
2279  word("PDF:pSet", "4" );
2280  parm("SigmaProcess:alphaSvalue", 0.1265);
2281  flag("SigmaTotal:zeroAXB", true );
2282  flag("SigmaDiffractive:dampen", false );
2283  parm("Diffraction:largeMassSuppress", 2.0 );
2284  flag("TimeShower:dampenBeamRecoil", true );
2285  flag("TimeShower:phiPolAsym", true );
2286  parm("SpaceShower:alphaSvalue", 0.130 );
2287  mode("SpaceShower:alphaSorder", 1 );
2288  flag("SpaceShower:alphaSuseCMW", false );
2289  flag("SpaceShower:samePTasMPI", false );
2290  parm("SpaceShower:pT0Ref", 2.0 );
2291  parm("SpaceShower:ecmRef", 1800.0);
2292  parm("SpaceShower:ecmPow", 0.0 );
2293  parm("SpaceShower:pTmaxFudge", 1.0 );
2294  parm("SpaceShower:pTdampFudge", 1.0 );
2295  flag("SpaceShower:rapidityOrder", true );
2296  flag("SpaceShower:rapidityOrderMPI", true );
2297  flag("SpaceShower:phiPolAsym", true );
2298  flag("SpaceShower:phiIntAsym", true );
2299  parm("MultipartonInteractions:alphaSvalue", 0.127 );
2300  parm("MultipartonInteractions:pT0Ref", 2.455 );
2301  parm("MultipartonInteractions:ecmRef", 1800. );
2302  parm("MultipartonInteractions:ecmPow", 0.26 );
2303  mode("MultipartonInteractions:bProfile", 3 );
2304  parm("MultipartonInteractions:expPow", 1.15 );
2305  parm("MultipartonInteractions:a1", 0.15 );
2306  parm("BeamRemnants:primordialKTsoft", 0.5 );
2307  parm("BeamRemnants:primordialKThard", 2.0 );
2308  parm("BeamRemnants:halfScaleForKT", 1.0 );
2309  parm("BeamRemnants:halfMassForKT", 1.0 );
2310  mode("ColourReconnection:mode", 0 );
2311  parm("ColourReconnection:range", 3.0 );
2312  }
2313 
2314  // Tune 4C, October 2010.
2315  else if (ppTune == 5) {
2316  word("PDF:pSet", "8" );
2317  parm("SigmaProcess:alphaSvalue", 0.135 );
2318  flag("SigmaTotal:zeroAXB", true );
2319  flag("SigmaDiffractive:dampen", true );
2320  parm("SigmaDiffractive:maxXB", 65.0 );
2321  parm("SigmaDiffractive:maxAX", 65.0 );
2322  parm("SigmaDiffractive:maxXX", 65.0 );
2323  parm("Diffraction:largeMassSuppress", 2.0 );
2324  flag("TimeShower:dampenBeamRecoil", true );
2325  flag("TimeShower:phiPolAsym", true );
2326  parm("SpaceShower:alphaSvalue", 0.137 );
2327  mode("SpaceShower:alphaSorder", 1 );
2328  flag("SpaceShower:alphaSuseCMW", false );
2329  flag("SpaceShower:samePTasMPI", false );
2330  parm("SpaceShower:pT0Ref", 2.0 );
2331  parm("SpaceShower:ecmRef", 1800.0);
2332  parm("SpaceShower:ecmPow", 0.0 );
2333  parm("SpaceShower:pTmaxFudge", 1.0 );
2334  parm("SpaceShower:pTdampFudge", 1.0 );
2335  flag("SpaceShower:rapidityOrder", true );
2336  flag("SpaceShower:rapidityOrderMPI", true );
2337  flag("SpaceShower:phiPolAsym", true );
2338  flag("SpaceShower:phiIntAsym", true );
2339  parm("MultipartonInteractions:alphaSvalue", 0.135 );
2340  parm("MultipartonInteractions:pT0Ref", 2.085 );
2341  parm("MultipartonInteractions:ecmRef", 1800. );
2342  parm("MultipartonInteractions:ecmPow", 0.19 );
2343  mode("MultipartonInteractions:bProfile", 3 );
2344  parm("MultipartonInteractions:expPow", 2.0 );
2345  parm("MultipartonInteractions:a1", 0.15 );
2346  parm("BeamRemnants:primordialKTsoft", 0.5 );
2347  parm("BeamRemnants:primordialKThard", 2.0 );
2348  parm("BeamRemnants:halfScaleForKT", 1.0 );
2349  parm("BeamRemnants:halfMassForKT", 1.0 );
2350  mode("ColourReconnection:mode", 0 );
2351  parm("ColourReconnection:range", 1.5 );
2352  }
2353 
2354  // Tune 4Cx, January 2011.
2355  else if (ppTune == 6) {
2356  word("PDF:pSet", "8" );
2357  parm("SigmaProcess:alphaSvalue", 0.135 );
2358  flag("SigmaTotal:zeroAXB", true );
2359  flag("SigmaDiffractive:dampen", true );
2360  parm("SigmaDiffractive:maxXB", 65.0 );
2361  parm("SigmaDiffractive:maxAX", 65.0 );
2362  parm("SigmaDiffractive:maxXX", 65.0 );
2363  parm("Diffraction:largeMassSuppress", 2.0 );
2364  flag("TimeShower:dampenBeamRecoil", true );
2365  flag("TimeShower:phiPolAsym", true );
2366  parm("SpaceShower:alphaSvalue", 0.137 );
2367  mode("SpaceShower:alphaSorder", 1 );
2368  flag("SpaceShower:alphaSuseCMW", false );
2369  flag("SpaceShower:samePTasMPI", false );
2370  parm("SpaceShower:pT0Ref", 2.0 );
2371  parm("SpaceShower:ecmRef", 1800.0);
2372  parm("SpaceShower:ecmPow", 0.0 );
2373  parm("SpaceShower:pTmaxFudge", 1.0 );
2374  parm("SpaceShower:pTdampFudge", 1.0 );
2375  flag("SpaceShower:rapidityOrder", true );
2376  flag("SpaceShower:rapidityOrderMPI", true );
2377  flag("SpaceShower:phiPolAsym", true );
2378  flag("SpaceShower:phiIntAsym", true );
2379  parm("MultipartonInteractions:alphaSvalue", 0.135 );
2380  parm("MultipartonInteractions:pT0Ref", 2.15 );
2381  parm("MultipartonInteractions:ecmRef", 1800. );
2382  parm("MultipartonInteractions:ecmPow", 0.19 );
2383  mode("MultipartonInteractions:bProfile", 4 );
2384  parm("MultipartonInteractions:expPow", 1.0 );
2385  parm("MultipartonInteractions:a1", 0.15 );
2386  parm("BeamRemnants:primordialKTsoft", 0.5 );
2387  parm("BeamRemnants:primordialKThard", 2.0 );
2388  parm("BeamRemnants:halfScaleForKT", 1.0 );
2389  parm("BeamRemnants:halfMassForKT", 1.0 );
2390  mode("ColourReconnection:mode", 0 );
2391  parm("ColourReconnection:range", 1.5 );
2392  }
2393 
2394  // The Monash 2013 tune by Peter Skands, the pp part (January 2014).
2395  else if (ppTune == 14) {
2396  word("PDF:pSet", "13" ); // NNPDF
2397  parm("SigmaProcess:alphaSvalue", 0.130 ); // same as PDF
2398  flag("SigmaTotal:zeroAXB", true );
2399  flag("SigmaDiffractive:dampen", true );
2400  parm("SigmaDiffractive:maxXB", 65.0 );
2401  parm("SigmaDiffractive:maxAX", 65.0 );
2402  parm("SigmaDiffractive:maxXX", 65.0 );
2403  parm("Diffraction:largeMassSuppress", 4.0 );
2404  flag("TimeShower:dampenBeamRecoil", true );
2405  flag("TimeShower:phiPolAsym", true );
2406  parm("SpaceShower:alphaSvalue", 0.1365); // same as FSR
2407  mode("SpaceShower:alphaSorder", 1 );
2408  flag("SpaceShower:alphaSuseCMW", false );
2409  flag("SpaceShower:samePTasMPI", false );
2410  parm("SpaceShower:pT0Ref", 2.0 );
2411  parm("SpaceShower:ecmRef", 7000.0);
2412  parm("SpaceShower:ecmPow", 0.0 );
2413  parm("SpaceShower:pTmaxFudge", 1.0 );
2414  parm("SpaceShower:pTdampFudge", 1.0 );
2415  flag("SpaceShower:rapidityOrder", true );
2416  flag("SpaceShower:rapidityOrderMPI", true );
2417  flag("SpaceShower:phiPolAsym", true );
2418  flag("SpaceShower:phiIntAsym", true );
2419  parm("MultipartonInteractions:alphaSvalue", 0.130 ); // same as PDF
2420  parm("MultipartonInteractions:pT0Ref", 2.28 );
2421  parm("MultipartonInteractions:ecmRef", 7000. );
2422  parm("MultipartonInteractions:ecmPow", 0.215 );
2423  mode("MultipartonInteractions:bProfile", 3 );
2424  parm("MultipartonInteractions:expPow", 1.85 );
2425  parm("MultipartonInteractions:a1", 0.15 );
2426  parm("BeamRemnants:primordialKTsoft", 0.9 );
2427  parm("BeamRemnants:primordialKThard", 1.8 );
2428  parm("BeamRemnants:halfScaleForKT", 1.5 );
2429  parm("BeamRemnants:halfMassForKT", 1.0 );
2430  mode("ColourReconnection:mode", 0 );
2431  parm("ColourReconnection:range", 1.80 );
2432  }
2433 
2434  // Several ATLAS and CMS tunes start out from Tune 4C.
2435  else if (ppTune < 18) {
2436  parm("SigmaProcess:alphaSvalue", 0.135 );
2437  flag("SigmaTotal:zeroAXB", true );
2438  flag("SigmaDiffractive:dampen", true );
2439  parm("SigmaDiffractive:maxXB", 65.0 );
2440  parm("SigmaDiffractive:maxAX", 65.0 );
2441  parm("SigmaDiffractive:maxXX", 65.0 );
2442  parm("Diffraction:largeMassSuppress", 2.0 );
2443  flag("TimeShower:dampenBeamRecoil", true );
2444  flag("TimeShower:phiPolAsym", true );
2445  parm("SpaceShower:alphaSvalue", 0.137 );
2446  mode("SpaceShower:alphaSorder", 1 );
2447  flag("SpaceShower:alphaSuseCMW", false );
2448  flag("SpaceShower:samePTasMPI", false );
2449  parm("SpaceShower:pT0Ref", 2.0 );
2450  parm("SpaceShower:ecmRef", 1800.0);
2451  parm("SpaceShower:ecmPow", 0.0 );
2452  parm("SpaceShower:pTmaxFudge", 1.0 );
2453  parm("SpaceShower:pTdampFudge", 1.0 );
2454  flag("SpaceShower:rapidityOrder", true );
2455  flag("SpaceShower:rapidityOrderMPI", true );
2456  flag("SpaceShower:phiPolAsym", true );
2457  flag("SpaceShower:phiIntAsym", true );
2458  parm("MultipartonInteractions:alphaSvalue", 0.135 );
2459  parm("MultipartonInteractions:pT0Ref", 2.085 );
2460  parm("MultipartonInteractions:ecmRef", 1800. );
2461  parm("MultipartonInteractions:ecmPow", 0.19 );
2462  mode("MultipartonInteractions:bProfile", 3 );
2463  parm("MultipartonInteractions:expPow", 2.0 );
2464  parm("MultipartonInteractions:a1", 0.15 );
2465  parm("BeamRemnants:primordialKTsoft", 0.5 );
2466  parm("BeamRemnants:primordialKThard", 2.0 );
2467  parm("BeamRemnants:halfScaleForKT", 1.0 );
2468  parm("BeamRemnants:halfMassForKT", 1.0 );
2469  mode("ColourReconnection:mode", 0 );
2470  parm("ColourReconnection:range", 1.5 );
2471 
2472  // Several ATLAS tunes in the A2 and AU2 series, see
2473  // ATLAS note ATL-PHYS-PUB-2012-003 (August 2012).
2474  // ATLAS MB tune A2-CTEQ6L1.
2475  if (ppTune == 7) {
2476  if (preferLHAPDF == 1)
2477  word("PDF:pSet", "LHAPDF5:cteq6ll.LHpdf");
2478  else if (preferLHAPDF == 2)
2479  word("PDF:pSet", "LHAPDF6:cteq6l1");
2480  else word("PDF:pSet", "8" );
2481  flag("SpaceShower:rapidityOrder", false );
2482  flag("SpaceShower:rapidityOrderMPI", false );
2483  parm("MultipartonInteractions:pT0Ref", 2.18 );
2484  parm("MultipartonInteractions:ecmPow", 0.22 );
2485  mode("MultipartonInteractions:bProfile", 4 );
2486  parm("MultipartonInteractions:expPow", 1.0 );
2487  parm("MultipartonInteractions:a1", 0.06 );
2488  parm("ColourReconnection:range", 1.55 );
2489  }
2490 
2491  // ATLAS MB tune A2-MSTW2008LO.
2492  else if (ppTune == 8) {
2493  if (preferLHAPDF == 1)
2494  word("PDF:pSet", "LHAPDF5:MSTW2008lo68cl.LHgrid");
2495  else if (preferLHAPDF == 2)
2496  word("PDF:pSet", "LHAPDF6:MSTW2008lo68cl");
2497  else word("PDF:pSet", "5" );
2498  flag("SpaceShower:rapidityOrder", false );
2499  flag("SpaceShower:rapidityOrderMPI", false );
2500  parm("MultipartonInteractions:pT0Ref", 1.90 );
2501  parm("MultipartonInteractions:ecmPow", 0.30 );
2502  mode("MultipartonInteractions:bProfile", 4 );
2503  parm("MultipartonInteractions:expPow", 1.0 );
2504  parm("MultipartonInteractions:a1", 0.03 );
2505  parm("ColourReconnection:range", 2.28 );
2506  }
2507 
2508  // ATLAS UE tune AU2-CTEQ6L1.
2509  if (ppTune == 9) {
2510  if (preferLHAPDF == 1)
2511  word("PDF:pSet", "LHAPDF5:cteq6ll.LHpdf");
2512  else if (preferLHAPDF == 2)
2513  word("PDF:pSet", "LHAPDF6:cteq6l1");
2514  else word("PDF:pSet", "8" );
2515  flag("SpaceShower:rapidityOrder", false );
2516  flag("SpaceShower:rapidityOrderMPI", false );
2517  parm("MultipartonInteractions:pT0Ref", 2.13 );
2518  parm("MultipartonInteractions:ecmPow", 0.21 );
2519  mode("MultipartonInteractions:bProfile", 4 );
2520  parm("MultipartonInteractions:expPow", 1.0 );
2521  parm("MultipartonInteractions:a1", 0.00 );
2522  parm("ColourReconnection:range", 2.21 );
2523  }
2524 
2525  // ATLAS UE tune AU2-MSTW2008LO.
2526  else if (ppTune == 10) {
2527  if (preferLHAPDF == 1)
2528  word("PDF:pSet", "LHAPDF5:MSTW2008lo68cl.LHgrid");
2529  else if (preferLHAPDF == 2)
2530  word("PDF:pSet", "LHAPDF6:MSTW2008lo68cl");
2531  else word("PDF:pSet", "5" );
2532  flag("SpaceShower:rapidityOrder", false );
2533  flag("SpaceShower:rapidityOrderMPI", false );
2534  parm("MultipartonInteractions:pT0Ref", 1.87 );
2535  parm("MultipartonInteractions:ecmPow", 0.28 );
2536  mode("MultipartonInteractions:bProfile", 4 );
2537  parm("MultipartonInteractions:expPow", 1.0 );
2538  parm("MultipartonInteractions:a1", 0.01 );
2539  parm("ColourReconnection:range", 5.32 );
2540  }
2541 
2542  // ATLAS UE tune AU2-CT10.
2543  else if (ppTune == 11) {
2544  if (preferLHAPDF == 2)
2545  word("PDF:pSet", "LHAPDF6:CT10");
2546  else
2547  word("PDF:pSet", "LHAPDF5:CT10.LHgrid");
2548  flag("SpaceShower:rapidityOrder", false );
2549  flag("SpaceShower:rapidityOrderMPI", false );
2550  parm("MultipartonInteractions:pT0Ref", 1.70 );
2551  parm("MultipartonInteractions:ecmPow", 0.16 );
2552  mode("MultipartonInteractions:bProfile", 4 );
2553  parm("MultipartonInteractions:expPow", 1.0 );
2554  parm("MultipartonInteractions:a1", 0.10 );
2555  parm("ColourReconnection:range", 4.67 );
2556  }
2557 
2558  // ATLAS UE tune AU2-MRST2007LO*.
2559  else if (ppTune == 12) {
2560  if (preferLHAPDF == 1)
2561  word("PDF:pSet", "LHAPDF5:MRST2007lomod.LHgrid");
2562  else if (preferLHAPDF == 2)
2563  word("PDF:pSet", "LHAPDF6:MRST2007lomod");
2564  else word("PDF:pSet", "3" );
2565  flag("SpaceShower:rapidityOrder", false );
2566  flag("SpaceShower:rapidityOrderMPI", false );
2567  parm("MultipartonInteractions:pT0Ref", 2.39 );
2568  parm("MultipartonInteractions:ecmPow", 0.24 );
2569  mode("MultipartonInteractions:bProfile", 4 );
2570  parm("MultipartonInteractions:expPow", 1.0 );
2571  parm("MultipartonInteractions:a1", 0.01 );
2572  parm("ColourReconnection:range", 1.76 );
2573  }
2574 
2575  // ATLAS UE tune AU2-MRST2007LO**.
2576  else if (ppTune == 13) {
2577  if (preferLHAPDF == 1)
2578  word("PDF:pSet", "LHAPDF5:MRSTMCal.LHgrid");
2579  else if (preferLHAPDF == 2)
2580  word("PDF:pSet", "LHAPDF6:MRSTMCal");
2581  else word("PDF:pSet", "4" );
2582  flag("SpaceShower:rapidityOrder", false );
2583  flag("SpaceShower:rapidityOrderMPI", false );
2584  parm("MultipartonInteractions:pT0Ref", 2.57 );
2585  parm("MultipartonInteractions:ecmPow", 0.23 );
2586  mode("MultipartonInteractions:bProfile", 4 );
2587  parm("MultipartonInteractions:expPow", 1.0 );
2588  parm("MultipartonInteractions:a1", 0.01 );
2589  parm("ColourReconnection:range", 1.47 );
2590  }
2591 
2592  // The CMS UE tunes CUETP8S1-CTEQ6L1 and CUETP8S1-HERAPDF1.5LO,
2593  // see the note CMS PAS GEN-14-001 (April 2014).
2594  // CMS UE tune CUETP8S1-CTEQ6L1.
2595  else if (ppTune == 15) {
2596  if (preferLHAPDF == 1)
2597  word("PDF:pSet", "LHAPDF5:cteq6ll.LHpdf");
2598  else if (preferLHAPDF == 2)
2599  word("PDF:pSet", "LHAPDF6:cteq6l1");
2600  else word("PDF:pSet", "8" );
2601  parm("MultipartonInteractions:pT0Ref", 2.1006);
2602  parm("MultipartonInteractions:ecmPow", 0.2106);
2603  parm("MultipartonInteractions:expPow", 1.6089);
2604  parm("MultipartonInteractions:a1", 0.00 );
2605  parm("ColourReconnection:range", 3.3126);
2606  }
2607 
2608  // CMS UE tune CUETP8S1-HERAPDF1.5LO.
2609  else if (ppTune == 16) {
2610  if (preferLHAPDF == 2)
2611  word("PDF:pSet", "LHAPDF6:HERAPDF15LO_EIG");
2612  else
2613  word("PDF:pSet", "LHAPDF5:HERAPDF1.5LO_EIG.LHgrid");
2614  parm("MultipartonInteractions:pT0Ref", 2.0001);
2615  parm("MultipartonInteractions:ecmPow", 0.2499);
2616  parm("MultipartonInteractions:expPow", 1.6905);
2617  parm("MultipartonInteractions:a1", 0.00 );
2618  parm("ColourReconnection:range", 6.0964);
2619  }
2620 
2621  // ATLAS tune AZ to the Z0/gamma* pTspectrum, see the note
2622  // CERN-PH-EP-2014-075 [arXiv:1406.3660 [hep-ex]] (June 2014).
2623  else if (ppTune == 17) {
2624  parm("SpaceShower:alphaSvalue", 0.1237);
2625  parm("SpaceShower:pT0Ref", 0.59 );
2626  parm("MultipartonInteractions:pT0Ref", 2.18 );
2627  parm("BeamRemnants:primordialKThard", 1.71 );
2628  }
2629  }
2630 
2631  // Several ATLAS and CMS tunes and tunes close-packing of strings
2632  // and hadron rescattering with start out from Monash 2013 tune.
2633  else if (ppTune >= 18) {
2634  word("PDF:pSet", "13" ); // NNPDF
2635  parm("SigmaProcess:alphaSvalue", 0.130 ); // same as PDF
2636  flag("SigmaTotal:zeroAXB", true );
2637  flag("SigmaDiffractive:dampen", true );
2638  parm("SigmaDiffractive:maxXB", 65.0 );
2639  parm("SigmaDiffractive:maxAX", 65.0 );
2640  parm("SigmaDiffractive:maxXX", 65.0 );
2641  parm("Diffraction:largeMassSuppress", 4.0 );
2642  flag("TimeShower:dampenBeamRecoil", true );
2643  flag("TimeShower:phiPolAsym", true );
2644  parm("SpaceShower:alphaSvalue", 0.1365); // same as FSR
2645  mode("SpaceShower:alphaSorder", 1 );
2646  flag("SpaceShower:alphaSuseCMW", false );
2647  flag("SpaceShower:samePTasMPI", false );
2648  parm("SpaceShower:pT0Ref", 2.0 );
2649  parm("SpaceShower:ecmRef", 7000.0);
2650  parm("SpaceShower:ecmPow", 0.0 );
2651  parm("SpaceShower:pTmaxFudge", 1.0 );
2652  parm("SpaceShower:pTdampFudge", 1.0 );
2653  flag("SpaceShower:rapidityOrder", true );
2654  flag("SpaceShower:rapidityOrderMPI", true );
2655  flag("SpaceShower:phiPolAsym", true );
2656  flag("SpaceShower:phiIntAsym", true );
2657  parm("MultipartonInteractions:alphaSvalue", 0.130 ); // same as PDF
2658  parm("MultipartonInteractions:pT0Ref", 2.28 );
2659  parm("MultipartonInteractions:ecmRef", 7000. );
2660  parm("MultipartonInteractions:ecmPow", 0.215 );
2661  mode("MultipartonInteractions:bProfile", 3 );
2662  parm("MultipartonInteractions:expPow", 1.85 );
2663  parm("MultipartonInteractions:a1", 0.15 );
2664  parm("BeamRemnants:primordialKTsoft", 0.9 );
2665  parm("BeamRemnants:primordialKThard", 1.8 );
2666  parm("BeamRemnants:halfScaleForKT", 1.5 );
2667  parm("BeamRemnants:halfMassForKT", 1.0 );
2668  mode("ColourReconnection:mode", 0 );
2669  parm("ColourReconnection:range", 1.80 );
2670 
2671  // CMS tune MonashStar = CUETP8M1-NNPDF2.3LO.
2672  // See R.D. Field, presentation at MPI@LHC 2014, Krakow, Poland.
2673  if (ppTune == 18) {
2674  parm("MultipartonInteractions:pT0Ref", 2.4024);
2675  parm("MultipartonInteractions:ecmPow", 0.25208);
2676  parm("MultipartonInteractions:expPow", 1.60 );
2677  }
2678 
2679  // The ATLAS A14 tunes, central tune with CTEQL1.
2680  // See ATL-PHYS-PUB-2014-021 (November 2014).
2681  // Warning: note that TimeShower:alphaSvalue is set here, although
2682  // normally it would be in the domain of ee tunes. This makes the
2683  // order of Tune:ee and Tune:pp commands relevant.
2684  else if (ppTune == 19) {
2685  if (preferLHAPDF == 1)
2686  word("PDF:pSet", "LHAPDF5:cteq6ll.LHpdf");
2687  else if (preferLHAPDF == 2)
2688  word("PDF:pSet", "LHAPDF6:cteq6l1");
2689  else word("PDF:pSet", "8" );
2690  parm("SigmaProcess:alphaSvalue", 0.144 );
2691  parm("TimeShower:alphaSvalue", 0.126 );
2692  parm("SpaceShower:alphaSvalue", 0.125 );
2693  parm("SpaceShower:pT0Ref", 1.3 );
2694  parm("SpaceShower:pTmaxFudge", 0.95 );
2695  parm("SpaceShower:pTdampFudge", 1.21 );
2696  parm("MultipartonInteractions:alphaSvalue",0.118);
2697  parm("MultipartonInteractions:pT0Ref", 1.98 );
2698  parm("BeamRemnants:primordialKThard", 1.72 );
2699  parm("ColourReconnection:range", 2.08 );
2700  }
2701 
2702  // The ATLAS A14 tunes, central tune with MSTW2008LO.
2703  else if (ppTune == 20) {
2704  if (preferLHAPDF == 1)
2705  word("PDF:pSet", "LHAPDF5:MSTW2008lo68cl.LHgrid");
2706  else if (preferLHAPDF == 2)
2707  word("PDF:pSet", "LHAPDF6:MSTW2008lo68cl");
2708  else word("PDF:pSet", "5" );
2709  parm("SigmaProcess:alphaSvalue", 0.140 );
2710  parm("TimeShower:alphaSvalue", 0.129 );
2711  parm("SpaceShower:alphaSvalue", 0.129 );
2712  parm("SpaceShower:pT0Ref", 1.62 );
2713  parm("SpaceShower:pTmaxFudge", 0.92 );
2714  parm("SpaceShower:pTdampFudge", 1.14 );
2715  parm("MultipartonInteractions:alphaSvalue",0.130);
2716  parm("MultipartonInteractions:pT0Ref", 2.28 );
2717  parm("BeamRemnants:primordialKThard", 1.82 );
2718  parm("ColourReconnection:range", 1.87 );
2719  }
2720 
2721  // The ATLAS A14 tunes, central tune with NNPDF2.3LO.
2722  else if (ppTune == 21) {
2723  word("PDF:pSet", "13" );
2724  parm("SigmaProcess:alphaSvalue", 0.140 );
2725  parm("TimeShower:alphaSvalue", 0.127 );
2726  parm("SpaceShower:alphaSvalue", 0.127 );
2727  parm("SpaceShower:pT0Ref", 1.56 );
2728  parm("SpaceShower:pTmaxFudge", 0.91 );
2729  parm("SpaceShower:pTdampFudge", 1.05 );
2730  parm("MultipartonInteractions:alphaSvalue",0.126);
2731  parm("MultipartonInteractions:pT0Ref", 2.09 );
2732  parm("BeamRemnants:primordialKThard", 1.88 );
2733  parm("ColourReconnection:range", 1.71 );
2734  }
2735 
2736  // The ATLAS A14 tunes, central tune with HERAPDF1.5LO.
2737  else if (ppTune == 22) {
2738  if (preferLHAPDF == 2)
2739  word("PDF:pSet", "LHAPDF6:HERAPDF15LO_EIG");
2740  else
2741  word("PDF:pSet", "LHAPDF5:HERAPDF1.5LO_EIG.LHgrid");
2742  parm("SigmaProcess:alphaSvalue", 0.141 );
2743  parm("TimeShower:alphaSvalue", 0.130 );
2744  parm("SpaceShower:alphaSvalue", 0.128);
2745  parm("SpaceShower:pT0Ref", 1.61 );
2746  parm("SpaceShower:pTmaxFudge", 0.95 );
2747  parm("SpaceShower:pTdampFudge", 1.10 );
2748  parm("MultipartonInteractions:alphaSvalue",0.123);
2749  parm("MultipartonInteractions:pT0Ref", 2.14 );
2750  parm("BeamRemnants:primordialKThard", 1.83 );
2751  parm("ColourReconnection:range", 1.78 );
2752  }
2753 
2754  // The ATLAS A14 tunes, variation 1+.
2755  else if (ppTune == 23) {
2756  word("PDF:pSet", "13" );
2757  parm("SigmaProcess:alphaSvalue", 0.140 );
2758  parm("TimeShower:alphaSvalue", 0.127 );
2759  parm("SpaceShower:alphaSvalue", 0.127 );
2760  parm("SpaceShower:pT0Ref", 1.56 );
2761  parm("SpaceShower:pTmaxFudge", 0.91 );
2762  parm("SpaceShower:pTdampFudge", 1.05 );
2763  parm("MultipartonInteractions:alphaSvalue",0.131);
2764  parm("MultipartonInteractions:pT0Ref", 2.09 );
2765  parm("BeamRemnants:primordialKThard", 1.88 );
2766  parm("ColourReconnection:range", 1.73 );
2767  }
2768 
2769  // The ATLAS A14 tunes, variation 1-.
2770  else if (ppTune == 24) {
2771  word("PDF:pSet", "13" );
2772  parm("SigmaProcess:alphaSvalue", 0.140 );
2773  parm("TimeShower:alphaSvalue", 0.127 );
2774  parm("SpaceShower:alphaSvalue", 0.127 );
2775  parm("SpaceShower:pT0Ref", 1.56 );
2776  parm("SpaceShower:pTmaxFudge", 0.91 );
2777  parm("SpaceShower:pTdampFudge", 1.05 );
2778  parm("MultipartonInteractions:alphaSvalue",0.121);
2779  parm("MultipartonInteractions:pT0Ref", 2.09 );
2780  parm("BeamRemnants:primordialKThard", 1.88 );
2781  parm("ColourReconnection:range", 1.69 );
2782  }
2783 
2784  // The ATLAS A14 tunes, variation 2+.
2785  else if (ppTune == 25) {
2786  word("PDF:pSet", "13" );
2787  parm("SigmaProcess:alphaSvalue", 0.140 );
2788  parm("TimeShower:alphaSvalue", 0.139 );
2789  parm("SpaceShower:alphaSvalue", 0.127 );
2790  parm("SpaceShower:pT0Ref", 1.60 );
2791  parm("SpaceShower:pTmaxFudge", 0.91 );
2792  parm("SpaceShower:pTdampFudge", 1.04 );
2793  parm("MultipartonInteractions:alphaSvalue",0.126);
2794  parm("MultipartonInteractions:pT0Ref", 2.09 );
2795  parm("BeamRemnants:primordialKThard", 1.88 );
2796  parm("ColourReconnection:range", 1.71 );
2797  }
2798 
2799  // The ATLAS A14 tunes, variation 2-.
2800  else if (ppTune == 26) {
2801  word("PDF:pSet", "13" );
2802  parm("SigmaProcess:alphaSvalue", 0.140 );
2803  parm("TimeShower:alphaSvalue", 0.111 );
2804  parm("SpaceShower:alphaSvalue", 0.127 );
2805  parm("SpaceShower:pT0Ref", 1.50 );
2806  parm("SpaceShower:pTmaxFudge", 0.91 );
2807  parm("SpaceShower:pTdampFudge", 1.08 );
2808  parm("MultipartonInteractions:alphaSvalue",0.126);
2809  parm("MultipartonInteractions:pT0Ref", 2.09 );
2810  parm("BeamRemnants:primordialKThard", 1.88 );
2811  parm("ColourReconnection:range", 1.71 );
2812  }
2813 
2814  // The ATLAS A14 tunes, variation 3a+.
2815  else if (ppTune == 27) {
2816  word("PDF:pSet", "13" );
2817  parm("SigmaProcess:alphaSvalue", 0.140 );
2818  parm("TimeShower:alphaSvalue", 0.136 );
2819  parm("SpaceShower:alphaSvalue", 0.127 );
2820  parm("SpaceShower:pT0Ref", 1.67 );
2821  parm("SpaceShower:pTmaxFudge", 0.98 );
2822  parm("SpaceShower:pTdampFudge", 1.36 );
2823  parm("MultipartonInteractions:alphaSvalue",0.125);
2824  parm("MultipartonInteractions:pT0Ref", 2.09 );
2825  parm("BeamRemnants:primordialKThard", 1.88 );
2826  parm("ColourReconnection:range", 1.71 );
2827  }
2828 
2829  // The ATLAS A14 tunes, variation 3a-.
2830  else if (ppTune == 28) {
2831  word("PDF:pSet", "13" );
2832  parm("SigmaProcess:alphaSvalue", 0.140 );
2833  parm("TimeShower:alphaSvalue", 0.124 );
2834  parm("SpaceShower:alphaSvalue", 0.127 );
2835  parm("SpaceShower:pT0Ref", 1.51 );
2836  parm("SpaceShower:pTmaxFudge", 0.88 );
2837  parm("SpaceShower:pTdampFudge", 0.93 );
2838  parm("MultipartonInteractions:alphaSvalue",0.127);
2839  parm("MultipartonInteractions:pT0Ref", 2.09 );
2840  parm("BeamRemnants:primordialKThard", 1.88 );
2841  parm("ColourReconnection:range", 1.71 );
2842  }
2843 
2844  // The ATLAS A14 tunes, variation 3b+.
2845  else if (ppTune == 29) {
2846  word("PDF:pSet", "13" );
2847  parm("SigmaProcess:alphaSvalue", 0.140 );
2848  parm("TimeShower:alphaSvalue", 0.114 );
2849  parm("SpaceShower:alphaSvalue", 0.129 );
2850  parm("SpaceShower:pT0Ref", 1.56 );
2851  parm("SpaceShower:pTmaxFudge", 1.00 );
2852  parm("SpaceShower:pTdampFudge", 1.04 );
2853  parm("MultipartonInteractions:alphaSvalue",0.126);
2854  parm("MultipartonInteractions:pT0Ref", 2.09 );
2855  parm("BeamRemnants:primordialKThard", 1.88 );
2856  parm("ColourReconnection:range", 1.71 );
2857  }
2858 
2859  // The ATLAS A14 tunes, variation 3b-.
2860  else if (ppTune == 30) {
2861  word("PDF:pSet", "13" );
2862  parm("SigmaProcess:alphaSvalue", 0.140 );
2863  parm("TimeShower:alphaSvalue", 0.138 );
2864  parm("SpaceShower:alphaSvalue", 0.126 );
2865  parm("SpaceShower:pT0Ref", 1.56 );
2866  parm("SpaceShower:pTmaxFudge", 0.83 );
2867  parm("SpaceShower:pTdampFudge", 1.07 );
2868  parm("MultipartonInteractions:alphaSvalue",0.126);
2869  parm("MultipartonInteractions:pT0Ref", 2.09 );
2870  parm("BeamRemnants:primordialKThard", 1.88 );
2871  parm("ColourReconnection:range", 1.71 );
2872  }
2873 
2874  // The ATLAS A14 tunes, variation 3c+.
2875  else if (ppTune == 31) {
2876  word("PDF:pSet", "13" );
2877  parm("SigmaProcess:alphaSvalue", 0.140 );
2878  parm("TimeShower:alphaSvalue", 0.127 );
2879  parm("SpaceShower:alphaSvalue", 0.140 );
2880  parm("SpaceShower:pT0Ref", 1.56 );
2881  parm("SpaceShower:pTmaxFudge", 0.91 );
2882  parm("SpaceShower:pTdampFudge", 1.05 );
2883  parm("MultipartonInteractions:alphaSvalue",0.126);
2884  parm("MultipartonInteractions:pT0Ref", 2.09 );
2885  parm("BeamRemnants:primordialKThard", 1.88 );
2886  parm("ColourReconnection:range", 1.71 );
2887  }
2888 
2889  // The ATLAS A14 tunes, variation 3c-.
2890  else if (ppTune == 32) {
2891  word("PDF:pSet", "13" );
2892  parm("SigmaProcess:alphaSvalue", 0.140 );
2893  parm("TimeShower:alphaSvalue", 0.127 );
2894  parm("SpaceShower:alphaSvalue", 0.115 );
2895  parm("SpaceShower:pT0Ref", 1.56 );
2896  parm("SpaceShower:pTmaxFudge", 0.91 );
2897  parm("SpaceShower:pTdampFudge", 1.05 );
2898  parm("MultipartonInteractions:alphaSvalue",0.126);
2899  parm("MultipartonInteractions:pT0Ref", 2.09 );
2900  parm("BeamRemnants:primordialKThard", 1.88 );
2901  parm("ColourReconnection:range", 1.71 );
2902  }
2903 
2904  // Tune with close-packing of strings and rescattering (November 2016).
2905  // Gaussian pT.
2906  else if (ppTune == 33) {
2907  parm("MultipartonInteractions:pT0Ref", 2.34 );
2908  parm("ColourReconnection:range", 1.8 );
2909  flag("StringPT:thermalModel", false);
2910  parm("StringPT:sigma", 0.33 );
2911  parm("StringPT:widthPreStrange", 1.2 );
2912  parm("StringPT:widthPreDiquark", 1.2 );
2913  parm("StringPT:enhancedFraction", 0.0 );
2914  flag("StringPT:closePacking", true );
2915  parm("StringPT:expNSP", 0.01 );
2916  parm("StringPT:expMPI", 0.0 );
2917  flag("HadronLevel:HadronScatter", true );
2918  mode("HadronScatter:mode", 0 );
2919  parm("HadronScatter:maxProbDS", 0.25 );
2920  }
2921 
2922  // Tune with close-packing of strings and rescattering (November 2016).
2923  // Thermodynamical pT.
2924  else if (ppTune == 34) {
2925  parm("MultipartonInteractions:pT0Ref", 2.5 );
2926  parm("ColourReconnection:range", 1.1 );
2927  flag("StringPT:thermalModel", true );
2928  parm("StringPT:temperature", 0.21 );
2929  parm("StringFlav:BtoMratio", 0.357);
2930  parm("StringFlav:StrangeSuppression", 0.5 );
2931  flag("StringPT:closePacking", true );
2932  parm("StringPT:expNSP", 0.13 );
2933  parm("StringPT:expMPI", 0.0 );
2934  flag("HadronLevel:HadronScatter", true );
2935  mode("HadronScatter:mode", 0 );
2936  parm("HadronScatter:maxProbDS", 0.5 );
2937  }
2938 
2939  }
2940 
2941 }
2942 
2943 //--------------------------------------------------------------------------
2944 
2945 // Allow several alternative inputs for true/false.
2946 
2947 bool Settings::boolString(string tag) {
2948 
2949  string tagLow = toLower(tag);
2950  return ( tagLow == "true" || tagLow == "1" || tagLow == "on"
2951  || tagLow == "yes" || tagLow == "ok" );
2952 
2953 }
2954 
2955 //--------------------------------------------------------------------------
2956 
2957 // Extract XML value string following XML attribute.
2958 
2959 string Settings::attributeValue(string line, string attribute) {
2960 
2961  if (line.find(attribute) == string::npos) return "";
2962  int iBegAttri = line.find(attribute);
2963  int iBegQuote = line.find("\"", iBegAttri + 1);
2964  int iEndQuote = line.find("\"", iBegQuote + 1);
2965  return line.substr(iBegQuote + 1, iEndQuote - iBegQuote - 1);
2966 
2967 }
2968 
2969 //--------------------------------------------------------------------------
2970 
2971 // Extract XML bool value following XML attribute.
2972 
2973 bool Settings::boolAttributeValue(string line, string attribute) {
2974 
2975  string valString = attributeValue(line, attribute);
2976  if (valString == "") return false;
2977  return boolString(valString);
2978 
2979 }
2980 
2981 //--------------------------------------------------------------------------
2982 
2983 // Extract XML int value following XML attribute.
2984 
2985 int Settings::intAttributeValue(string line, string attribute) {
2986  string valString = attributeValue(line, attribute);
2987  if (valString == "") return 0;
2988  istringstream valStream(valString);
2989  int intVal;
2990  valStream >> intVal;
2991  return intVal;
2992 
2993 }
2994 
2995 //--------------------------------------------------------------------------
2996 
2997 // Extract XML double value following XML attribute.
2998 
2999 double Settings::doubleAttributeValue(string line, string attribute) {
3000  string valString = attributeValue(line, attribute);
3001  if (valString == "") return 0.;
3002  istringstream valStream(valString);
3003  double doubleVal;
3004  valStream >> doubleVal;
3005  return doubleVal;
3006 
3007 }
3008 
3009 //--------------------------------------------------------------------------
3010 
3011 // Extract XML bool vector value following XML attribute.
3012 
3013 vector<bool> Settings::boolVectorAttributeValue(string line,
3014  string attribute) {
3015  string valString = attributeValue(line, attribute);
3016  if (valString == "") return vector<bool>(1, false);
3017  vector<bool> vectorVal;
3018  size_t stringPos(0);
3019  while (stringPos != string::npos) {
3020  stringPos = valString.find(",");
3021  istringstream valStream(valString.substr(0, stringPos));
3022  valString = valString.substr(stringPos + 1);
3023  vectorVal.push_back(boolString(valStream.str()));
3024  }
3025  return vectorVal;
3026 
3027 }
3028 
3029 //--------------------------------------------------------------------------
3030 
3031 // Extract XML int vector value following XML attribute.
3032 
3033 vector<int> Settings::intVectorAttributeValue(string line,
3034  string attribute) {
3035  string valString = attributeValue(line, attribute);
3036  if (valString == "") return vector<int>(1, 0);
3037  int intVal;
3038  vector<int> vectorVal;
3039  size_t stringPos(0);
3040  while (stringPos != string::npos) {
3041  stringPos = valString.find(",");
3042  istringstream valStream(valString.substr(0, stringPos));
3043  valString = valString.substr(stringPos + 1);
3044  valStream >> intVal;
3045  vectorVal.push_back(intVal);
3046  }
3047  return vectorVal;
3048 
3049 }
3050 
3051 //--------------------------------------------------------------------------
3052 
3053 // Extract XML double vector value following XML attribute.
3054 
3055 vector<double> Settings::doubleVectorAttributeValue(string line,
3056  string attribute) {
3057  string valString = attributeValue(line, attribute);
3058  if (valString == "") return vector<double>(1, 0.);
3059  double doubleVal;
3060  vector<double> vectorVal;
3061  size_t stringPos(0);
3062  while (stringPos != string::npos) {
3063  stringPos = valString.find(",");
3064  istringstream valStream(valString.substr(0, stringPos));
3065  valString = valString.substr(stringPos + 1);
3066  valStream >> doubleVal;
3067  vectorVal.push_back(doubleVal);
3068  }
3069  return vectorVal;
3070 
3071 }
3072 
3073 //--------------------------------------------------------------------------
3074 
3075 // Extract XML string vector value following XML attribute.
3076 
3077 vector<string> Settings::stringVectorAttributeValue(string line,
3078  string attribute) {
3079  string valString = attributeValue(line, attribute);
3080  if (valString == "") return vector<string>(1, " ");
3081  string stringVal;
3082  vector<string> vectorVal;
3083  size_t stringPos(0);
3084  while (stringPos != string::npos) {
3085  stringPos = valString.find(",");
3086  if (stringPos != string::npos) {
3087  vectorVal.push_back(valString.substr(0, stringPos));
3088  valString = valString.substr(stringPos + 1);
3089  } else vectorVal.push_back(valString);
3090  }
3091  return vectorVal;
3092 
3093 }
3094 
3095 //==========================================================================
3096 
3097 } // end namespace Pythia8