StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
ParticleData.cc
1 // ParticleData.cc is a part of the PYTHIA event generator.
2 // Copyright (C) 2020 Torbjorn Sjostrand.
3 // PYTHIA is licenced under the GNU GPL v2 or later, see COPYING for details.
4 // Please respect the MCnet Guidelines, see GUIDELINES for details.
5 
6 // Function definitions (not found in the header) for the
7 // DecayChannel, ParticleDataEntry and ParticleData classes.
8 
9 #include "Pythia8/ParticleData.h"
10 #include "Pythia8/ResonanceWidths.h"
11 #include "Pythia8/StandardModel.h"
12 #include "Pythia8/SusyResonanceWidths.h"
13 #include "Pythia8/ResonanceWidthsDM.h"
14 
15 // Allow string and character manipulation.
16 #include <cctype>
17 
18 namespace Pythia8 {
19 
20 //==========================================================================
21 
22 // DecayChannel class.
23 // This class holds info on a single decay channel.
24 
25 //--------------------------------------------------------------------------
26 
27 // Check whether id1 occurs anywhere in product list.
28 
29 bool DecayChannel::contains(int id1) const {
30 
31  bool found1 = false;
32  for (int i = 0; i < nProd; ++i) if (prod[i] == id1) found1 = true;
33  return found1;
34 
35 }
36 
37 //--------------------------------------------------------------------------
38 
39 // Check whether id1 and id2 occur anywhere in product list.
40 // iF ID1 == ID2 then two copies of this particle must be present.
41 
42 bool DecayChannel::contains(int id1, int id2) const {
43 
44  bool found1 = false;
45  bool found2 = false;
46  for (int i = 0; i < nProd; ++i) {
47  if (!found1 && prod[i] == id1) {found1 = true; continue;}
48  if (!found2 && prod[i] == id2) {found2 = true; continue;}
49  }
50  return found1 && found2;
51 
52 }
53 
54 //--------------------------------------------------------------------------
55 
56 // Check whether id1, id2 and id3 occur anywhere in product list.
57 // iF ID1 == ID2 then two copies of this particle must be present, etc.
58 
59 bool DecayChannel::contains(int id1, int id2, int id3) const {
60 
61  bool found1 = false;
62  bool found2 = false;
63  bool found3 = false;
64  for (int i = 0; i < nProd; ++i) {
65  if (!found1 && prod[i] == id1) {found1 = true; continue;}
66  if (!found2 && prod[i] == id2) {found2 = true; continue;}
67  if (!found3 && prod[i] == id3) {found3 = true; continue;}
68  }
69  return found1 && found2 && found3;
70 
71 }
72 
73 //==========================================================================
74 
75 // ParticleDataEntry class.
76 // This class holds info on a single particle species.
77 
78 //--------------------------------------------------------------------------
79 
80 // Constants: could be changed here if desired, but normally should not.
81 // These are of technical nature, as described for each.
82 
83 // A particle is invisible if it has neither strong nor electric charge,
84 // and is not made up of constituents that have it. Only relevant for
85 // long-lived particles. This list may need to be extended.
86 const int ParticleDataEntry::INVISIBLENUMBER = 62;
87 const int ParticleDataEntry::INVISIBLETABLE[80] = {
88  12, 14, 16, 18, 23, 25, 32, 33,
89  35, 36, 39, 41, 45, 46, 51, 52,
90  53, 54, 55, 56, 57, 58, 59, 60,
91  1000012, 1000014, 1000016, 1000018, 1000022, 1000023, 1000025, 1000035,
92  1000045, 1000039, 2000012, 2000014, 2000016, 2000018, 4900012, 4900014,
93  4900016, 4900021, 4900022, 4900101, 4900102, 4900103, 4900104, 4900105,
94  4900106, 4900107, 4900108, 4900111, 4900113, 4900211, 4900213, 4900991,
95  5000039, 5100039, 9900012, 9900014, 9900016, 9900023, 0, 0,
96  0, 0, 0, 0, 0, 0, 0, 0,
97  0, 0, 0, 0, 0, 0, 0, 0
98 };
99 
100 // For some particles we know it is necessary to switch off width,
101 // although they do have one, so do not warn.
102 const int ParticleDataEntry::KNOWNNOWIDTH[3] = {10313, 10323, 10333};
103 
104 // Particles with a read-in tau0 (in mm/c) below this mayDecay by default.
105 const double ParticleDataEntry::MAXTAU0FORDECAY = 1000.;
106 
107 // Particles with a read-in m0 above this isResonance by default.
108 const double ParticleDataEntry::MINMASSRESONANCE = 20.;
109 
110 // Narrow states are assigned nominal mass.
111 const double ParticleDataEntry::NARROWMASS = 1e-6;
112 
113 // Constituent quark masses (d, u, s, c, b, -, -, -, g).
114 const double ParticleDataEntry::CONSTITUENTMASSTABLE[10]
115  = {0., 0.325, 0.325, 0.50, 1.60, 5.00, 0., 0., 0., 0.7};
116 
117 //--------------------------------------------------------------------------
118 
119 // Destructor: delete any ResonanceWidths object.
120 
121 ParticleDataEntry::~ParticleDataEntry() {
122  if (resonancePtr != 0) delete resonancePtr;
123 }
124 
125 //--------------------------------------------------------------------------
126 
127 // Set initial default values for some quantities.
128 
129 void ParticleDataEntry::setDefaults() {
130 
131  // A particle is a resonance if it is heavy enough.
132  isResonanceSave = (m0Save > MINMASSRESONANCE);
133 
134  // A particle may decay if it is shortlived enough.
135  mayDecaySave = (tau0Save < MAXTAU0FORDECAY);
136 
137  // A particle's lifetime is calculated from its decay width.
138  tauCalcSave = true;
139 
140  // A particle by default has no external decays.
141  doExternalDecaySave = false;
142 
143  // A particle is invisible if in current table of such.
144  isVisibleSave = true;
145  for (int i = 0; i < INVISIBLENUMBER; ++i)
146  if (idSave == INVISIBLETABLE[i]) isVisibleSave = false;
147 
148  // Normally a resonance should not have width forced to fixed value.
149  doForceWidthSave = false;
150 
151  // Set up constituent masses.
152  setConstituentMass();
153 
154  // No Breit-Wigner mass selection before initialized. Status tau0.
155  modeBWnow = 0;
156  modeTau0now = 0;
157 
158 }
159 
160 //--------------------------------------------------------------------------
161 
162 // Find out if a particle is a hadron.
163 // Only covers normal hadrons, not e.g. R-hadrons.
164 
165 bool ParticleDataEntry::isHadron() const {
166 
167  if (idSave <= 100 || (idSave >= 1000000 && idSave <= 9000000)
168  || idSave >= 9900000) return false;
169  if (idSave == 130 || idSave == 310) return true;
170  if (idSave%10 == 0 || (idSave/10)%10 == 0 || (idSave/100)%10 == 0)
171  return false;
172  return true;
173 
174 }
175 
176 //--------------------------------------------------------------------------
177 
178 // Find out if a particle is a meson.
179 // Only covers normal hadrons, not e.g. R-hadrons.
180 
181 bool ParticleDataEntry::isMeson() const {
182 
183  if (idSave <= 100 || (idSave >= 1000000 && idSave <= 9000000)
184  || idSave >= 9900000) return false;
185  if (idSave == 130 || idSave == 310) return true;
186  if (idSave%10 == 0 || (idSave/10)%10 == 0 || (idSave/100)%10 == 0
187  || (idSave/1000)%10 != 0) return false;
188  return true;
189 
190 }
191 
192 //--------------------------------------------------------------------------
193 
194 // Find out if a particle is a baryon.
195 // Only covers normal hadrons, not e.g. R-hadrons.
196 
197 bool ParticleDataEntry::isBaryon() const {
198 
199  if (idSave <= 1000 || (idSave >= 1000000 && idSave <= 9000000)
200  || idSave >= 9900000) return false;
201  if (idSave%10 == 0 || (idSave/10)%10 == 0 || (idSave/100)%10 == 0
202  || (idSave/1000)%10 == 0) return false;
203  return true;
204 
205 }
206 
207 //--------------------------------------------------------------------------
208 
209 // Find out if a particle is quarkonia.
210 
211 bool ParticleDataEntry::isOnium() const {
212  if (idSave%2 != 1) return false;
213  if (idSave > 1000000) return false;
214  if ((idSave/10)%10 < 4) return false;
215  if ((idSave/10)%10 > 6) return false;
216  if ((idSave/10)%10 != (idSave/100)%10) return false;
217  if ((idSave/1000)%10 != 0) return false;
218  return true;
219 
220 }
221 
222 //--------------------------------------------------------------------------
223 
224 // Extract the heaviest (= largest id) quark in a hadron.
225 
226 int ParticleDataEntry::heaviestQuark(int idIn) const {
227 
228  if (!isHadron()) return 0;
229  int hQ = 0;
230 
231  // Meson.
232  if ( (idSave/1000)%10 == 0 ) {
233  hQ = (idSave/100)%10;
234  if (idSave == 130) hQ = 3;
235  if (hQ%2 == 1) hQ = -hQ;
236 
237  // Baryon.
238  } else hQ = (idSave/1000)%10;
239 
240  // Done.
241  return (idIn > 0) ? hQ : -hQ;
242 
243 }
244 
245 //--------------------------------------------------------------------------
246 
247 // Calculate three times baryon number, i.e. net quark - antiquark number.
248 
249 int ParticleDataEntry::baryonNumberType(int idIn) const {
250 
251  // Quarks.
252  if (isQuark()) return (idIn > 0) ? 1 : -1;
253 
254  // Diquarks
255  if (isDiquark()) return (idIn > 0) ? 2 : -2;
256 
257  // Baryons.
258  if (isBaryon()) return (idIn > 0) ? 3 : -3;
259 
260  // Done.
261  return 0;
262 
263 }
264 
265 //--------------------------------------------------------------------------
266 
267 // Find number of quarks of given kind inside quark, diquark or hadron.
268 // Note: naive answer for flavour-diagonal meson mixing.
269 
270 int ParticleDataEntry::nQuarksInCode(int idQIn) const {
271 
272  // Do not keep track of sign.
273  int idQ = abs(idQIn);
274  int idNow = abs(idSave);
275  int nQ = 0;
276 
277  // Quarks.
278  if (isQuark()) return (idQ == idNow) ? 1 : 0;
279 
280  // Diquarks.
281  if (isDiquark()) {
282  if ( (idNow/1000) % 10 == idQ) ++nQ;
283  if ( (idNow/100) % 10 == idQ) ++nQ;
284  return nQ;
285  }
286 
287  // Mesons.
288  if (isMeson()) {
289  if ( (idNow/100) % 10 == idQ) ++nQ;
290  if ( (idNow/10) % 10 == idQ) ++nQ;
291  return nQ;
292  }
293 
294  // Baryons.
295  if (isBaryon()) {
296  if ( (idNow/1000) % 10 == idQ) ++nQ;
297  if ( (idNow/100) % 10 == idQ) ++nQ;
298  if ( (idNow/10) % 10 == idQ) ++nQ;
299  return nQ;
300  }
301 
302  // Done. Room for improvements e.g. w.r.t. R-hadrons.
303  return 0;
304 
305 }
306 
307 //--------------------------------------------------------------------------
308 
309 // Prepare the Breit-Wigner mass selection by precalculating
310 // frequently-used expressions.
311 
312 void ParticleDataEntry::initBWmass() {
313 
314  // Optionally set decay vertices also for short-lived particles.
315  // (Lifetimes are explicitly tabulated for long-lived ones.)
316  if (modeTau0now == 0) modeTau0now = (particleDataPtr->setRapidDecayVertex
317  && tau0Save == 0. && channels.size() > 0) ? 2 : 1;
318  if (modeTau0now == 2) tau0Save = (mWidthSave > NARROWMASS)
319  ? HBARC * FM2MM / mWidthSave : particleDataPtr->intermediateTau0;
320 
321  // Find Breit-Wigner mode for current particle.
322  modeBWnow = particleDataPtr->modeBreitWigner;
323  if ( m0Save < NARROWMASS ) mWidthSave = 0.;
324  if ( mWidthSave < NARROWMASS || (mMaxSave > mMinSave
325  && mMaxSave - mMinSave < NARROWMASS) ) modeBWnow = 0;
326  if (modeBWnow == 0) {
327  mMinSave = mMaxSave = m0Save;
328  return;
329  }
330 
331  // Find atan expressions to be used in random mass selection.
332  if (modeBWnow < 3) {
333  atanLow = atan( 2. * (mMinSave - m0Save) / mWidthSave );
334  double atanHigh = (mMaxSave > mMinSave)
335  ? atan( 2. * (mMaxSave - m0Save) / mWidthSave ) : 0.5 * M_PI;
336  atanDif = atanHigh - atanLow;
337  } else {
338  atanLow = atan( (pow2(mMinSave) - pow2(m0Save))
339  / (m0Save * mWidthSave) );
340  double atanHigh = (mMaxSave > mMinSave)
341  ? atan( (pow2(mMaxSave) - pow2(m0Save)) / (m0Save * mWidthSave) )
342  : 0.5 * M_PI;
343  atanDif = atanHigh - atanLow;
344  }
345 
346  // Done if no threshold factor.
347  if (modeBWnow%2 == 1) return;
348 
349  // Find average mass threshold for threshold-factor correction.
350  double bRatSum = 0.;
351  double mThrSum = 0;
352  for (int i = 0; i < int(channels.size()); ++i)
353  if (channels[i].onMode() > 0) {
354  bRatSum += channels[i].bRatio();
355  double mChannelSum = 0.;
356  for (int j = 0; j < channels[i].multiplicity(); ++j)
357  mChannelSum += particleDataPtr->m0( channels[i].product(j) );
358  mThrSum += channels[i].bRatio() * mChannelSum;
359  }
360  mThr = (bRatSum == 0.) ? 0. : mThrSum / bRatSum;
361 
362  // Switch off Breit-Wigner if very close to threshold.
363  if (mThr + NARROWMASS > m0Save && !isResonanceSave) {
364  modeBWnow = 0;
365  bool knownProblem = false;
366  for (int i = 0; i < 3; ++i) if (idSave == KNOWNNOWIDTH[i])
367  knownProblem = true;
368  if (!knownProblem) {
369  ostringstream osWarn;
370  osWarn << "for id = " << idSave;
371  particleDataPtr->infoPtr->errorMsg("Warning in ParticleDataEntry::"
372  "initBWmass: switching off width", osWarn.str(), true);
373  }
374  }
375 
376 }
377 
378 //--------------------------------------------------------------------------
379 
380 // Function to give mass of a particle, either at the nominal value
381 // or picked according to a (linear or quadratic) Breit-Wigner.
382 
383 double ParticleDataEntry::mSel() const {
384 
385  // Nominal value. (Width check should not be needed, but just in case.)
386  if (modeBWnow == 0 || mWidthSave < NARROWMASS) return m0Save;
387  double mNow, m2Now;
388 
389  // Mass according to a Breit-Wigner linear in m.
390  if (modeBWnow == 1) {
391  mNow = m0Save + 0.5 * mWidthSave
392  * tan( atanLow + atanDif * particleDataPtr->rndmPtr->flat() );
393 
394  // Ditto, but make Gamma proportional to sqrt(m^2 - m_threshold^2).
395  } else if (modeBWnow == 2) {
396  double mWidthNow, fixBW, runBW;
397  double m0ThrS = m0Save*m0Save - mThr*mThr;
398  do {
399  mNow = m0Save + 0.5 * mWidthSave
400  * tan( atanLow + atanDif * particleDataPtr->rndmPtr->flat() );
401  mWidthNow = mWidthSave * sqrtpos( (mNow*mNow - mThr*mThr) / m0ThrS );
402  fixBW = mWidthSave / (pow2(mNow - m0Save) + pow2(0.5 * mWidthSave));
403  runBW = mWidthNow / (pow2(mNow - m0Save) + pow2(0.5 * mWidthNow));
404  } while (runBW < particleDataPtr->rndmPtr->flat()
405  * particleDataPtr->maxEnhanceBW * fixBW);
406 
407  // Mass according to a Breit-Wigner quadratic in m.
408  } else if (modeBWnow == 3) {
409  m2Now = m0Save*m0Save + m0Save * mWidthSave
410  * tan( atanLow + atanDif * particleDataPtr->rndmPtr->flat() );
411  mNow = sqrtpos( m2Now);
412 
413  // Ditto, but m_0 Gamma_0 -> m Gamma(m) with threshold factor as above.
414  } else {
415  double mwNow, fixBW, runBW;
416  double m2Ref = m0Save * m0Save;
417  double mwRef = m0Save * mWidthSave;
418  double m2Thr = mThr * mThr;
419  do {
420  m2Now = m2Ref + mwRef * tan( atanLow + atanDif
421  * particleDataPtr->rndmPtr->flat() );
422  mNow = sqrtpos( m2Now);
423  mwNow = mNow * mWidthSave
424  * sqrtpos( (m2Now - m2Thr) / (m2Ref - m2Thr) );
425  fixBW = mwRef / (pow2(m2Now - m2Ref) + pow2(mwRef));
426  runBW = mwNow / (pow2(m2Now - m2Ref) + pow2(mwNow));
427  } while (runBW < particleDataPtr->rndmPtr->flat()
428  * particleDataPtr->maxEnhanceBW * fixBW);
429  }
430 
431  // Done.
432  return mNow;
433 }
434 
435 //--------------------------------------------------------------------------
436 
437 // Function to calculate running mass at given mass scale.
438 
439 double ParticleDataEntry::mRun(double mHat) const {
440 
441  // Except for six quarks return nominal mass.
442  if (idSave > 6) return m0Save;
443  double mQRun = particleDataPtr->mQRun[idSave];
444  double Lam5 = particleDataPtr->Lambda5Run;
445 
446  // For d, u, s quarks start running at 2 GeV (RPP 2006 p. 505).
447  if (idSave < 4) return mQRun * pow ( log(2. / Lam5)
448  / log(max(2., mHat) / Lam5), 12./23.);
449 
450  // For c, b and t quarks start running at respective mass.
451  return mQRun * pow ( log(mQRun / Lam5)
452  / log(max(mQRun, mHat) / Lam5), 12./23.);
453 
454 }
455 
456 //--------------------------------------------------------------------------
457 
458 // Rescale all branching ratios to assure normalization to unity.
459 
460 void ParticleDataEntry::rescaleBR(double newSumBR) {
461 
462  // Sum up branching ratios. Find rescaling factor. Rescale.
463  double oldSumBR = 0.;
464  for ( int i = 0; i < int(channels.size()); ++ i)
465  oldSumBR += channels[i].bRatio();
466  double rescaleFactor = newSumBR / oldSumBR;
467  for ( int i = 0; i < int(channels.size()); ++ i)
468  channels[i].rescaleBR(rescaleFactor);
469 
470 }
471 
472 //--------------------------------------------------------------------------
473 
474 // Prepare to pick a decay channel.
475 
476 bool ParticleDataEntry::preparePick(int idSgn, double mHat, int idInFlav) {
477 
478  // Reset sum of allowed widths/branching ratios.
479  currentBRSum = 0.;
480 
481  // For resonances the widths are calculated dynamically.
482  if (isResonanceSave && resonancePtr != 0) {
483  resonancePtr->widthStore(idSgn, mHat, idInFlav);
484  for (int i = 0; i < int(channels.size()); ++i)
485  currentBRSum += channels[i].currentBR();
486 
487  // Else use normal fixed branching ratios.
488  } else {
489  int onMode;
490  double currentBRNow;
491  for (int i = 0; i < int(channels.size()); ++i) {
492  onMode = channels[i].onMode();
493  currentBRNow = 0.;
494  if ( idSgn > 0 && (onMode == 1 || onMode == 2) )
495  currentBRNow = channels[i].bRatio();
496  else if ( idSgn < 0 && (onMode == 1 || onMode == 3) )
497  currentBRNow = channels[i].bRatio();
498  channels[i].currentBR(currentBRNow);
499  currentBRSum += currentBRNow;
500  }
501  }
502 
503  // Failure if no channels found with positive branching ratios.
504  return (currentBRSum > 0.);
505 
506 }
507 
508 //--------------------------------------------------------------------------
509 
510 // Pick a decay channel according to branching ratios from preparePick.
511 
512 DecayChannel& ParticleDataEntry::pickChannel() {
513 
514  // Find channel in table.
515  int size = channels.size();
516  double rndmBR = currentBRSum * particleDataPtr->rndmPtr->flat();
517  int i = -1;
518  do rndmBR -= channels[++i].currentBR();
519  while (rndmBR > 0. && i < size);
520 
521  // Emergency if no channel found. Done.
522  if (i == size) i = 0;
523  return channels[i];
524 
525 }
526 
527 //--------------------------------------------------------------------------
528 
529 // Access methods stored in ResonanceWidths. Could have been
530 // inline in .h, except for problems with forward declarations.
531 
532 void ParticleDataEntry::setResonancePtr(
533  ResonanceWidths* resonancePtrIn) {
534  if (resonancePtr == resonancePtrIn) return;
535  if (resonancePtr != 0) delete resonancePtr;
536  resonancePtr = resonancePtrIn;
537 }
538 
539 void ParticleDataEntry::resInit(Info* infoPtr) {
540  if (resonancePtr != 0) resonancePtr->init(infoPtr);
541 }
542 
543 double ParticleDataEntry::resWidth(int idSgn, double mHat, int idIn,
544  bool openOnly, bool setBR) {
545  return (resonancePtr != 0) ? resonancePtr->width( idSgn, mHat,
546  idIn, openOnly, setBR) : 0.;
547 }
548 
549 double ParticleDataEntry::resWidthOpen(int idSgn, double mHat, int idIn) {
550  return (resonancePtr != 0) ? resonancePtr->widthOpen( idSgn, mHat, idIn)
551  : 0.;
552 }
553 
554 double ParticleDataEntry::resWidthStore(int idSgn, double mHat, int idIn) {
555  return (resonancePtr != 0) ? resonancePtr->widthStore( idSgn, mHat, idIn)
556  : 0.;
557 }
558 
559 double ParticleDataEntry::resOpenFrac(int idSgn) {
560  return (resonancePtr != 0) ? resonancePtr->openFrac(idSgn) : 1.;
561 }
562 
563 double ParticleDataEntry::resWidthRescaleFactor() {
564  return (resonancePtr != 0) ? resonancePtr->widthRescaleFactor() : 1.;
565 }
566 
567 double ParticleDataEntry::resWidthChan(double mHat, int idAbs1,
568  int idAbs2) {
569  return (resonancePtr != 0) ? resonancePtr->widthChan( mHat, idAbs1,
570  idAbs2) : 0.;
571 }
572 
573 //--------------------------------------------------------------------------
574 
575 // Constituent masses for (d, u, s, c, b) quarks and diquarks.
576 // Hardcoded in CONSTITUENTMASSTABLE so that they are not overwritten
577 // by mistake, and separated from the "normal" masses.
578 // Called both by setDefaults and setM0 so kept as separate method.
579 
580 void ParticleDataEntry::setConstituentMass() {
581 
582  // Equate with the normal masses as default guess.
583  constituentMassSave = m0Save;
584 
585  // Quark masses trivial. Also gluon mass.
586  if (idSave < 6) constituentMassSave = CONSTITUENTMASSTABLE[idSave];
587  if (idSave == 21) constituentMassSave = CONSTITUENTMASSTABLE[9];
588 
589  // Diquarks as simple sum of constituent quarks.
590  if (idSave > 1000 && idSave < 10000 && (idSave/10)%10 == 0) {
591  int id1 = idSave/1000;
592  int id2 = (idSave/100)%10;
593  if (id1 <6 && id2 < 6) constituentMassSave
594  = CONSTITUENTMASSTABLE[id1] + CONSTITUENTMASSTABLE[id2];
595  }
596 
597 }
598 
599 //==========================================================================
600 
601 // ParticleData class.
602 // This class holds a map of all ParticleDataEntries,
603 // each entry containing info on a particle species.
604 
605 //--------------------------------------------------------------------------
606 
607 // Get data to be distributed among particles during setup.
608 // Note: this routine is called twice. Firstly from init(...), but
609 // the data should not be used at that point, so is likely overkill.
610 // Secondly, from initWidths, after user has had time to change.
611 
612 void ParticleData::initCommon() {
613 
614  // Mass generation: fixed mass or linear/quadratic Breit-Wigner.
615  modeBreitWigner = settingsPtr->mode("ParticleData:modeBreitWigner");
616 
617  // Maximum tail enhancement when adding threshold factor to Breit-Wigner.
618  maxEnhanceBW = settingsPtr->parm("ParticleData:maxEnhanceBW");
619 
620  // Find initial MSbar masses for five light flavours.
621  mQRun[1] = settingsPtr->parm("ParticleData:mdRun");
622  mQRun[2] = settingsPtr->parm("ParticleData:muRun");
623  mQRun[3] = settingsPtr->parm("ParticleData:msRun");
624  mQRun[4] = settingsPtr->parm("ParticleData:mcRun");
625  mQRun[5] = settingsPtr->parm("ParticleData:mbRun");
626  mQRun[6] = settingsPtr->parm("ParticleData:mtRun");
627 
628  // Find Lambda5 value to use in running of MSbar masses.
629  double alphaSvalue = settingsPtr->parm("ParticleData:alphaSvalueMRun");
630  AlphaStrong alphaS;
631  alphaS.init( alphaSvalue, 1, 5, false);
632  Lambda5Run = alphaS.Lambda5();
633 
634  // Set secondary vertices also for rapidly decaying particles.
635  setRapidDecayVertex = settingsPtr->flag("HadronLevel:Rescatter")
636  || ( settingsPtr->flag("Fragmentation:setVertices")
637  && settingsPtr->flag("HadronVertex:rapidDecays") );
638  intermediateTau0 = settingsPtr->parm("HadronVertex:intermediateTau0");
639 
640 }
641 
642 //--------------------------------------------------------------------------
643 
644 // Initialize pointer for particles to the full database, the Breit-Wigners
645 // of normal hadrons and the ResonanceWidths of resonances. For the latter
646 // the order of initialization is essential to get secondary widths right.
647 
648 void ParticleData::initWidths( vector<ResonanceWidths*> resonancePtrs) {
649 
650  // Initialize some common data (but preserve history of read statements).
651  initCommon();
652 
653  // Pointer to database and Breit-Wigner mass initialization for each
654  // particle.
655  ResonanceWidths* resonancePtr = 0;
656  for (map<int, ParticleDataEntry>::iterator pdtEntry = pdt.begin();
657  pdtEntry != pdt.end(); ++pdtEntry) {
658  ParticleDataEntry& pdtNow = pdtEntry->second;
659  pdtNow.initBWmass();
660 
661  // Remove any existing resonances.
662  resonancePtr = pdtNow.getResonancePtr();
663  if (resonancePtr != 0) pdtNow.setResonancePtr(0);
664  }
665 
666  // Begin set up new resonance objects.
667  // Z0, W+- and top are almost always needed.
668  resonancePtr = new ResonanceGmZ(23);
669  setResonancePtr( 23, resonancePtr);
670  resonancePtr = new ResonanceW(24);
671  setResonancePtr( 24, resonancePtr);
672  resonancePtr = new ResonanceTop(6);
673  setResonancePtr( 6, resonancePtr);
674 
675  // Higgs in SM.
676  if (!settingsPtr->flag("Higgs:useBSM")) {
677  resonancePtr = new ResonanceH(0, 25);
678  setResonancePtr( 25, resonancePtr);
679 
680  // Higgses in BSM.
681  } else {
682  resonancePtr = new ResonanceH(1, 25);
683  setResonancePtr( 25, resonancePtr);
684  resonancePtr = new ResonanceH(2, 35);
685  setResonancePtr( 35, resonancePtr);
686  resonancePtr = new ResonanceH(3, 36);
687  setResonancePtr( 36, resonancePtr);
688  resonancePtr = new ResonanceHchg(37);
689  setResonancePtr( 37, resonancePtr);
690  resonancePtr = new ResonanceH(4, 45);
691  setResonancePtr( 45, resonancePtr);
692  resonancePtr = new ResonanceH(5, 46);
693  setResonancePtr( 46, resonancePtr);
694  }
695 
696  // A fourth generation: b', t', tau', nu'_tau.
697  resonancePtr = new ResonanceFour(7);
698  setResonancePtr( 7, resonancePtr);
699  resonancePtr = new ResonanceFour(8);
700  setResonancePtr( 8, resonancePtr);
701  resonancePtr = new ResonanceFour(17);
702  setResonancePtr( 17, resonancePtr);
703  resonancePtr = new ResonanceFour(18);
704  setResonancePtr( 18, resonancePtr);
705 
706  // New gauge bosons: Z', W', R.
707  resonancePtr = new ResonanceZprime(32);
708  setResonancePtr( 32, resonancePtr);
709  resonancePtr = new ResonanceWprime(34);
710  setResonancePtr( 34, resonancePtr);
711  resonancePtr = new ResonanceRhorizontal(41);
712  setResonancePtr( 41, resonancePtr);
713 
714  // A leptoquark.
715  resonancePtr = new ResonanceLeptoquark(42);
716  setResonancePtr( 42, resonancePtr);
717 
718  // Mediators for Dark Matter.
719  resonancePtr = new ResonanceS(54);
720  setResonancePtr( 54, resonancePtr);
721  resonancePtr = new ResonanceZp(55);
722  setResonancePtr( 55, resonancePtr);
723  resonancePtr = new ResonanceSl(56);
724  setResonancePtr( 56, resonancePtr);
725  resonancePtr = new ResonanceCha(57);
726  setResonancePtr( 57, resonancePtr);
727  resonancePtr = new ResonanceDM2(58);
728  setResonancePtr( 58, resonancePtr);
729  resonancePtr = new ResonanceChaD(59);
730  setResonancePtr( 59, resonancePtr);
731 
732  // 93 = Z0copy and 94 = W+-copy used to pick decay channels
733  // for W/Z production in parton showers.
734  resonancePtr = new ResonanceGmZ(93);
735  setResonancePtr( 93, resonancePtr);
736  resonancePtr = new ResonanceW(94);
737  setResonancePtr( 94, resonancePtr);
738 
739  // Supersymmetry:
740  // - Squarks;
741  for(int i = 1; i < 7; i++){
742  resonancePtr = new ResonanceSquark(1000000 + i);
743  setResonancePtr( 1000000 + i, resonancePtr);
744  resonancePtr = new ResonanceSquark(2000000 + i);
745  setResonancePtr( 2000000 + i, resonancePtr);
746  }
747 
748  // - Sleptons and sneutrinos;
749  for(int i = 1; i < 7; i++){
750  resonancePtr = new ResonanceSlepton(1000010 + i);
751  setResonancePtr( 1000010 + i, resonancePtr);
752  resonancePtr = new ResonanceSlepton(2000010 + i);
753  setResonancePtr( 2000010 + i, resonancePtr);
754  }
755 
756  // - Gluino;
757  resonancePtr = new ResonanceGluino(1000021);
758  setResonancePtr( 1000021, resonancePtr);
759 
760  // - Charginos;
761  resonancePtr = new ResonanceChar(1000024);
762  setResonancePtr( 1000024, resonancePtr);
763  resonancePtr = new ResonanceChar(1000037);
764  setResonancePtr( 1000037, resonancePtr);
765 
766  // - Neutralinos.
767  if (isResonance(1000022)) {
768  resonancePtr = new ResonanceNeut(1000022);
769  setResonancePtr( 1000022, resonancePtr);
770  }
771  resonancePtr = new ResonanceNeut(1000023);
772  setResonancePtr( 1000023, resonancePtr);
773  resonancePtr = new ResonanceNeut(1000025);
774  setResonancePtr( 1000025, resonancePtr);
775  resonancePtr = new ResonanceNeut(1000035);
776  setResonancePtr( 1000035, resonancePtr);
777  resonancePtr = new ResonanceNeut(1000045);
778  setResonancePtr( 1000045, resonancePtr);
779 
780  // Excited quarks and leptons.
781  for (int i = 1; i < 7; ++i) {
782  resonancePtr = new ResonanceExcited(4000000 + i);
783  setResonancePtr( 4000000 + i, resonancePtr);
784  }
785  for (int i = 11; i < 17; ++i) {
786  resonancePtr = new ResonanceExcited(4000000 + i);
787  setResonancePtr( 4000000 + i, resonancePtr);
788  }
789 
790  // An excited graviton/gluon in extra-dimensional scenarios.
791  resonancePtr = new ResonanceGraviton(5100039);
792  setResonancePtr( 5100039, resonancePtr);
793  resonancePtr = new ResonanceKKgluon(5100021);
794  setResonancePtr( 5100021, resonancePtr);
795 
796  // A left-right-symmetric scenario with new righthanded neutrinos,
797  // righthanded gauge bosons and doubly charged Higgses.
798  resonancePtr = new ResonanceNuRight(9900012);
799  setResonancePtr( 9900012, resonancePtr);
800  resonancePtr = new ResonanceNuRight(9900014);
801  setResonancePtr( 9900014, resonancePtr);
802  resonancePtr = new ResonanceNuRight(9900016);
803  setResonancePtr( 9900016, resonancePtr);
804  resonancePtr = new ResonanceZRight(9900023);
805  setResonancePtr( 9900023, resonancePtr);
806  resonancePtr = new ResonanceWRight(9900024);
807  setResonancePtr( 9900024, resonancePtr);
808  resonancePtr = new ResonanceHchgchgLeft(9900041);
809  setResonancePtr( 9900041, resonancePtr);
810  resonancePtr = new ResonanceHchgchgRight(9900042);
811  setResonancePtr( 9900042, resonancePtr);
812 
813  // Attach user-defined external resonances and do basic initialization.
814  for (int i = 0; i < int(resonancePtrs.size()); ++i) {
815  int idNow = resonancePtrs[i]->id();
816  resonancePtrs[i]->initBasic(idNow);
817  setResonancePtr( idNow, resonancePtrs[i]);
818  }
819 
820  // Set up lists to order resonances in ascending mass.
821  vector<int> idOrdered;
822  vector<double> m0Ordered;
823 
824  // Put Z0 and W+- first, since known to be SM and often off-shell.
825  idOrdered.push_back(23);
826  m0Ordered.push_back(m0(23));
827  idOrdered.push_back(24);
828  m0Ordered.push_back(m0(24));
829 
830  // Loop through particle table to find resonances.
831  for (map<int, ParticleDataEntry>::iterator pdtEntry = pdt.begin();
832  pdtEntry != pdt.end(); ++pdtEntry) {
833  ParticleDataEntry& pdtNow = pdtEntry->second;
834  int idNow = pdtNow.id();
835 
836  // Set up a simple default object for uninitialized resonances.
837  if (pdtNow.isResonance() && pdtNow.getResonancePtr() == 0) {
838  resonancePtr = new ResonanceGeneric(idNow);
839  setResonancePtr( idNow, resonancePtr);
840  }
841 
842  // Insert resonances in ascending mass, to respect decay hierarchies.
843  if (pdtNow.getResonancePtr() != 0 && idNow != 23 && idNow != 24) {
844  double m0Now = pdtNow.m0();
845  idOrdered.push_back(idNow);
846  m0Ordered.push_back(m0Now);
847  for (int i = int(idOrdered.size()) - 2; i > 1; --i) {
848  if (m0Ordered[i] < m0Now) break;
849  swap( idOrdered[i], idOrdered[i + 1]);
850  swap( m0Ordered[i], m0Ordered[i + 1]);
851  }
852  }
853  }
854 
855  // Initialize the resonances in ascending mass order. Reset mass generation.
856  for (int i = 0; i < int(idOrdered.size()); ++i) {
857  resInit( idOrdered[i]);
858  ParticleDataEntry* pdtPtrNow = particleDataEntryPtr( idOrdered[i]);
859  pdtPtrNow->initBWmass();
860  }
861 
862 }
863 
864 //--------------------------------------------------------------------------
865 
866 // Read in database from specific XML file (which may refer to others).
867 
868 bool ParticleData::readXML(string inFile, bool reset) {
869 
870  // Load XML file into memory
871  if (!loadXML(inFile,reset)) return false;
872 
873  // Process XML file (now stored in memory)
874  if (!processXML(reset)) return false;
875 
876  // Done.
877  return true;
878 }
879 
880 //--------------------------------------------------------------------------
881 
882 // Read in database from specific XML stream (which may refer to others).
883 
884 bool ParticleData::readXML(istream &inStr, bool reset) {
885 
886  // Load XML file into memory
887  if (!loadXML(inStr,reset)) return false;
888 
889  // Process XML file (now stored in memory)
890  if (!processXML(reset)) return false;
891 
892  // Done.
893  return true;
894 }
895 
896 //--------------------------------------------------------------------------
897 
898 // Read in database from pre-initialised particleData object.
899 
900 bool ParticleData::copyXML(const ParticleData &particleDataIn) {
901 
902  // First Reset everything.
903  pdt.clear();
904  xmlFileSav.clear();
905  readStringHistory.resize(0);
906  readStringSubrun.clear();
907  isInit = false;
908  xmlFileSav=particleDataIn.xmlFileSav;
909 
910  // Process XML file (now stored in memory)
911  if (!processXML(true)) return false;
912 
913  // Done.
914  return true;
915 }
916 //--------------------------------------------------------------------------
917 
918 // Load a specific XML file into memory (which may refer to others).
919 
920 bool ParticleData::loadXML(istream& is, bool reset) {
921 
922  // Normally reset whole database before beginning.
923  if (reset) {
924  pdt.clear();
925  xmlFileSav.clear();
926  readStringHistory.resize(0);
927  readStringSubrun.clear();
928  isInit = false;
929  }
930 
931  // Check that instream is OK.
932  if (!is.good()) {
933  infoPtr->errorMsg("Error in ParticleData::readXML:"
934  " did not find data");
935  return false;
936  }
937 
938  // Read in one line at a time.
939  particlePtr = 0;
940  string line;
941  while ( getline(is, line) ) {
942 
943  // Get first word of a line.
944  istringstream getfirst(line);
945  string word1;
946  getfirst >> word1;
947 
948  // Check for occurence of a file also to be read.
949  if (word1 == "<file") {
950  string file = attributeValue(line, "name");
951  }
952 
953  // Else save line to memory.
954  else {
955  xmlFileSav.push_back(line);
956  }
957  }
958 
959  //Done.
960  return true;
961 
962 }
963 
964 
965 //--------------------------------------------------------------------------
966 
967 // Load a specific XML file into memory (which may refer to others).
968 
969 bool ParticleData::loadXML(string inFile, bool reset) {
970 
971  const char* cstring = inFile.c_str();
972  ifstream is(cstring);
973 
974  return loadXML(is, reset);
975 }
976 
977 //--------------------------------------------------------------------------
978 
979 // Process XML contents stored in memory.
980 
981 bool ParticleData::processXML(bool reset) {
982 
983  // Number of lines saved.
984  int nLines = xmlFileSav.size();
985 
986  // Process each line sequentially.
987  particlePtr = 0;
988  int i=-1;
989  while (++i < nLines) {
990 
991  // Retrieve line.
992  string line = xmlFileSav[i];
993 
994  // Get first word of a line.
995  istringstream getfirst(line);
996  string word1;
997  getfirst >> word1;
998 
999  // Check for occurence of a particle. Add any continuation lines.
1000  if (word1 == "<particle") {
1001  while (line.find(">") == string::npos) {
1002  if (++i >= nLines) break;
1003  string addLine = xmlFileSav[i];
1004  line += addLine;
1005  }
1006 
1007  // Read in particle properties.
1008  int idTmp = intAttributeValue( line, "id");
1009  string nameTmp = attributeValue( line, "name");
1010  string antiNameTmp = attributeValue( line, "antiName");
1011  if (antiNameTmp == "") antiNameTmp = "void";
1012  int spinTypeTmp = intAttributeValue( line, "spinType");
1013  int chargeTypeTmp = intAttributeValue( line, "chargeType");
1014  int colTypeTmp = intAttributeValue( line, "colType");
1015  double m0Tmp = doubleAttributeValue( line, "m0");
1016  double mWidthTmp = doubleAttributeValue( line, "mWidth");
1017  double mMinTmp = doubleAttributeValue( line, "mMin");
1018  double mMaxTmp = doubleAttributeValue( line, "mMax");
1019  double tau0Tmp = doubleAttributeValue( line, "tau0");
1020  bool varWidthTmp = boolAttributeValue( line, "varWidth");
1021 
1022  // Erase if particle already exists.
1023  if (isParticle(idTmp)) pdt.erase(idTmp);
1024 
1025  // Store new particle. Save pointer, to be used for decay channels.
1026  addParticle( idTmp, nameTmp, antiNameTmp, spinTypeTmp, chargeTypeTmp,
1027  colTypeTmp, m0Tmp, mWidthTmp, mMinTmp, mMaxTmp, tau0Tmp,
1028  varWidthTmp);
1029  particlePtr = particleDataEntryPtr(idTmp);
1030 
1031  // Check for occurence of a decay channel. Add any continuation lines.
1032  } else if (word1 == "<channel") {
1033  while (line.find(">") == string::npos) {
1034  if (++i >= nLines) break;
1035  string addLine = xmlFileSav[i];
1036  line += addLine;
1037  }
1038 
1039  // Read in channel properties - products so far only as a string.
1040  int onMode = intAttributeValue( line, "onMode");
1041  double bRatio = doubleAttributeValue( line, "bRatio");
1042  int meMode = intAttributeValue( line, "meMode");
1043  string products = attributeValue( line, "products");
1044 
1045  // Read in decay products from stream. Must have at least one.
1046  istringstream prodStream(products);
1047  int prod0 = 0; int prod1 = 0; int prod2 = 0; int prod3 = 0;
1048  int prod4 = 0; int prod5 = 0; int prod6 = 0; int prod7 = 0;
1049  prodStream >> prod0 >> prod1 >> prod2 >> prod3 >> prod4 >> prod5
1050  >> prod6 >> prod7;
1051  if (prod0 == 0) {
1052  infoPtr->errorMsg("Error in ParticleData::readXML:"
1053  " incomplete decay channel", line);
1054  return false;
1055  }
1056 
1057  // Store new channel (if particle already known).
1058  if (particlePtr == 0) {
1059  infoPtr->errorMsg("Error in ParticleData::readXML:"
1060  " orphan decay channel", line);
1061  return false;
1062  }
1063  particlePtr->addChannel(onMode, bRatio, meMode, prod0, prod1,
1064  prod2, prod3, prod4, prod5, prod6, prod7);
1065 
1066  // End of loop over lines in input file and loop over files.
1067  };
1068  };
1069 
1070  // All particle data at this stage defines baseline original.
1071  if (reset) for (map<int, ParticleDataEntry>::iterator pdtEntry
1072  = pdt.begin(); pdtEntry != pdt.end(); ++pdtEntry) {
1073  particlePtr = &pdtEntry->second;
1074  particlePtr->setHasChanged(false);
1075  }
1076 
1077  // Done.
1078  isInit = true;
1079  return true;
1080 
1081 }
1082 
1083 //--------------------------------------------------------------------------
1084 
1085 // Print out complete database in numerical order as an XML file.
1086 
1087 void ParticleData::listXML(string outFile) {
1088 
1089  // Convert file name to ofstream.
1090  const char* cstring = outFile.c_str();
1091  ofstream os(cstring);
1092 
1093  // Iterate through the particle data table.
1094  for (map<int, ParticleDataEntry>::iterator pdtEntry
1095  = pdt.begin(); pdtEntry != pdt.end(); ++pdtEntry) {
1096  particlePtr = &pdtEntry->second;
1097 
1098  // Print particle properties.
1099  os << "<particle id=\"" << particlePtr->id() << "\""
1100  << " name=\"" << particlePtr->name() << "\"";
1101  if (particlePtr->hasAnti())
1102  os << " antiName=\"" << particlePtr->name(-1) << "\"";
1103  os << " spinType=\"" << particlePtr->spinType() << "\""
1104  << " chargeType=\"" << particlePtr->chargeType() << "\""
1105  << " colType=\"" << particlePtr->colType() << "\"\n";
1106  // Pick format for mass and width based on mass value.
1107  double m0Now = particlePtr->m0();
1108  if (m0Now == 0 || (m0Now > 0.1 && m0Now < 1000.))
1109  os << fixed << setprecision(5);
1110  else os << scientific << setprecision(3);
1111  os << " m0=\"" << m0Now << "\"";
1112  if (particlePtr->mWidth() > 0.)
1113  os << " mWidth=\"" << particlePtr->mWidth() << "\""
1114  << " mMin=\"" << particlePtr->mMin() << "\""
1115  << " mMax=\"" << particlePtr->mMax() << "\"";
1116  if (particlePtr->tau0() > 0.) os << scientific << setprecision(5)
1117  << " tau0=\"" << particlePtr->tau0() << "\"";
1118  os << ">\n";
1119 
1120  // Loop through the decay channel table for each particle.
1121  if (particlePtr->sizeChannels() > 0) {
1122  for (int i = 0; i < particlePtr->sizeChannels(); ++i) {
1123  const DecayChannel& channel = particlePtr->channel(i);
1124  int mult = channel.multiplicity();
1125 
1126  // Print decay channel properties.
1127  os << " <channel onMode=\"" << channel.onMode() << "\""
1128  << fixed << setprecision(7)
1129  << " bRatio=\"" << channel.bRatio() << "\"";
1130  if (channel.meMode() > 0)
1131  os << " meMode=\"" << channel.meMode() << "\"";
1132  os << " products=\"";
1133  for (int j = 0; j < mult; ++j) {
1134  os << channel.product(j);
1135  if (j < mult - 1) os << " ";
1136  }
1137 
1138  // Finish off line and loop over allowed decay channels.
1139  os << "\"/>\n";
1140  }
1141  }
1142 
1143  // Finish off existing particle.
1144  os << "</particle>\n\n";
1145 
1146  }
1147 
1148 }
1149 
1150 //--------------------------------------------------------------------------
1151 
1152 // Read in database from specific free format file.
1153 
1154 bool ParticleData::readFF(istream& is, bool reset) {
1155 
1156  // Normally reset whole database before beginning.
1157  if (reset) {
1158  pdt.clear();
1159  readStringHistory.resize(0);
1160  readStringSubrun.clear();
1161  isInit = false;
1162  }
1163 
1164  if (!is.good()) {
1165  infoPtr->errorMsg("Error in ParticleData::readFF:"
1166  " did not find stream");
1167  return false;
1168  }
1169 
1170  // Read in one line at a time.
1171  particlePtr = 0;
1172  string line;
1173  bool readParticle = false;
1174  while ( getline(is, line) ) {
1175 
1176  // Empty lines begins new particle.
1177  if (line.find_first_not_of(" \n\t\v\b\r\f\a") == string::npos) {
1178  readParticle = true;
1179  continue;
1180  }
1181 
1182  // Prepare to use standard read from line.
1183  istringstream readLine(line);
1184 
1185  // Read in a line with particle information.
1186  if (readParticle) {
1187 
1188  // Properties to be read.
1189  int idTmp;
1190  string nameTmp, antiNameTmp;
1191  int spinTypeTmp, chargeTypeTmp, colTypeTmp;
1192  double m0Tmp, mWidthTmp, mMinTmp, mMaxTmp, tau0Tmp;
1193  bool varWidthTmp;
1194  string mayTmp;
1195 
1196  // Do the reading.
1197  readLine >> idTmp >> nameTmp >> antiNameTmp >> spinTypeTmp
1198  >> chargeTypeTmp >> colTypeTmp >> m0Tmp >> mWidthTmp
1199  >> mMinTmp >> mMaxTmp >> tau0Tmp >> varWidthTmp;
1200 
1201  // Error printout if something went wrong.
1202  if (!readLine) {
1203  infoPtr->errorMsg("Error in ParticleData::readFF:"
1204  " incomplete particle", line);
1205  return false;
1206  }
1207 
1208  // Erase if particle already exists.
1209  if (isParticle(idTmp)) pdt.erase(idTmp);
1210 
1211  // Store new particle. Save pointer, to be used for decay channels.
1212  addParticle( idTmp, nameTmp, antiNameTmp, spinTypeTmp, chargeTypeTmp,
1213  colTypeTmp, m0Tmp, mWidthTmp, mMinTmp, mMaxTmp, tau0Tmp, varWidthTmp);
1214  particlePtr = particleDataEntryPtr(idTmp);
1215  readParticle = false;
1216 
1217  // Read in a line with decay channel information.
1218  } else {
1219 
1220  // Properties to be read.
1221  int onMode = 0;
1222  double bRatio = 0.;
1223  int meMode = 0;
1224  int prod0 = 0; int prod1 = 0; int prod2 = 0; int prod3 = 0;
1225  int prod4 = 0; int prod5 = 0; int prod6 = 0; int prod7 = 0;
1226 
1227  // Read in data from stream. Need at least one decay product.
1228  readLine >> onMode >> bRatio >> meMode >> prod0;
1229  if (!readLine) {
1230  infoPtr->errorMsg("Error in ParticleData::readFF:"
1231  " incomplete decay channel", line);
1232  return false;
1233  }
1234  readLine >> prod1 >> prod2 >> prod3 >> prod4 >> prod5
1235  >> prod6 >> prod7;
1236 
1237  // Store new channel.
1238  if (particlePtr == 0) {
1239  infoPtr->errorMsg("Error in ParticleData::readFF:"
1240  " orphan decay channel", line);
1241  return false;
1242  }
1243  particlePtr->addChannel(onMode, bRatio, meMode, prod0, prod1,
1244  prod2, prod3, prod4, prod5, prod6, prod7);
1245 
1246  }
1247 
1248  // End of while loop over lines in the file.
1249  }
1250 
1251 
1252  // Done.
1253  isInit = true;
1254  return true;
1255 
1256 }
1257 
1258 
1259 //--------------------------------------------------------------------------
1260 
1261 // Read in database from specific free format file.
1262 
1263 bool ParticleData::readFF(string inFile, bool reset) {
1264 
1265  const char* cstring = inFile.c_str();
1266  ifstream is(cstring);
1267 
1268  return readFF(is,reset);
1269 }
1270 
1271 //--------------------------------------------------------------------------
1272 
1273 // Print out complete database in numerical order as a free format file.
1274 
1275 void ParticleData::listFF(string outFile) {
1276 
1277  // Convert file name to ofstream.
1278  const char* cstring = outFile.c_str();
1279  ofstream os(cstring);
1280 
1281  // Iterate through the particle data table.
1282  for (map<int, ParticleDataEntry>::iterator pdtEntry
1283  = pdt.begin(); pdtEntry != pdt.end(); ++pdtEntry) {
1284  particlePtr = &pdtEntry->second;
1285 
1286  // Pick format for mass and width based on mass value.
1287  double m0Now = particlePtr->m0();
1288  if (m0Now == 0 || (m0Now > 0.1 && m0Now < 1000.))
1289  os << fixed << setprecision(5);
1290  else os << scientific << setprecision(3);
1291 
1292  // Print particle properties.
1293  os << "\n" << setw(8) << particlePtr->id() << " "
1294  << left << setw(16) << particlePtr->name() << " "
1295  << setw(16) << particlePtr->name(-1) << " "
1296  << right << setw(2) << particlePtr->spinType() << " "
1297  << setw(2) << particlePtr->chargeType() << " "
1298  << setw(2) << particlePtr->colType() << " "
1299  << setw(10) << particlePtr->m0() << " "
1300  << setw(10) << particlePtr->mWidth() << " "
1301  << setw(10) << particlePtr->mMin() << " "
1302  << setw(10) << particlePtr->mMax() << " "
1303  << scientific << setprecision(5)
1304  << setw(12) << particlePtr->tau0() << "\n";
1305 
1306  // Loop through the decay channel table for each particle.
1307  if (particlePtr->sizeChannels() > 0) {
1308  for (int i = 0; i < particlePtr->sizeChannels(); ++i) {
1309  const DecayChannel& channel = particlePtr->channel(i);
1310  os << " " << setw(6) << channel.onMode()
1311  << " " << fixed << setprecision(7) << setw(10)
1312  << channel.bRatio() << " "
1313  << setw(3) << channel.meMode() << " ";
1314  for (int j = 0; j < channel.multiplicity(); ++j)
1315  os << setw(8) << channel.product(j) << " ";
1316  os << "\n";
1317  }
1318  }
1319 
1320  }
1321 
1322 }
1323 
1324 //--------------------------------------------------------------------------
1325 
1326 // Read in updates from a character string, like a line of a file.
1327 // Is used by readString (and readFile) in Pythia.
1328 
1329  bool ParticleData::readString(string lineIn, bool warn) {
1330 
1331  // If empty line then done.
1332  if (lineIn.find_first_not_of(" \n\t\v\b\r\f\a") == string::npos) return true;
1333 
1334  // Take copy that will be modified.
1335  string line = lineIn;
1336 
1337  // If first character is not a digit then taken to be a comment.
1338  int firstChar = line.find_first_not_of(" \n\t\v\b\r\f\a");
1339  if (!isdigit(line[firstChar])) return true;
1340 
1341  // Replace colons and equal signs by blanks to make parsing simpler.
1342  for ( int j = 0; j < int(line.size()); ++ j)
1343  if (line[j] == ':' || line[j] == '=') line[j] = ' ';
1344 
1345  // Get particle id and property name.
1346  int idTmp;
1347  string property;
1348  istringstream getWord(line);
1349  getWord >> idTmp >> property;
1350  toLowerRep(property);
1351 
1352  // Check that valid particle.
1353  if ( (!isParticle(idTmp) && property != "all" && property != "new")
1354  || idTmp <= 0) {
1355  if (warn) cout << "\n PYTHIA Error: input particle not found in Particle"
1356  << " Data Table:\n " << lineIn << "\n";
1357  readingFailedSave = true;
1358  return false;
1359  }
1360 
1361  // Store history of readString statements.
1362  readStringHistory.push_back(line);
1363  int subrun = max( -1, settingsPtr->mode("Main:subrun"));
1364  if (readStringSubrun.find(subrun) == readStringSubrun.end())
1365  readStringSubrun[subrun] = vector<string>();
1366  readStringSubrun[subrun].push_back(line);
1367 
1368  // Identify particle property and read + set its value, case by case.
1369  if (property == "name") {
1370  string nameTmp;
1371  getWord >> nameTmp;
1372  pdt[idTmp].setName(nameTmp);
1373  return true;
1374  }
1375  if (property == "antiname") {
1376  string antiNameTmp;
1377  getWord >> antiNameTmp;
1378  pdt[idTmp].setAntiName(antiNameTmp);
1379  return true;
1380  }
1381  if (property == "names") {
1382  string nameTmp, antiNameTmp;
1383  getWord >> nameTmp >> antiNameTmp;
1384  pdt[idTmp].setNames(nameTmp, antiNameTmp);
1385  return true;
1386  }
1387  if (property == "spintype") {
1388  int spinTypeTmp;
1389  getWord >> spinTypeTmp;
1390  pdt[idTmp].setSpinType(spinTypeTmp);
1391  return true;
1392  }
1393  if (property == "chargetype") {
1394  int chargeTypeTmp;
1395  getWord >> chargeTypeTmp;
1396  pdt[idTmp].setChargeType(chargeTypeTmp);
1397  return true;
1398  }
1399  if (property == "coltype") {
1400  int colTypeTmp;
1401  getWord >> colTypeTmp;
1402  pdt[idTmp].setColType(colTypeTmp);
1403  return true;
1404  }
1405  if (property == "m0") {
1406  double m0Tmp;
1407  getWord >> m0Tmp;
1408  pdt[idTmp].setM0(m0Tmp);
1409  return true;
1410  }
1411  if (property == "mwidth") {
1412  double mWidthTmp;
1413  getWord >> mWidthTmp;
1414  pdt[idTmp].setMWidth(mWidthTmp);
1415  return true;
1416  }
1417  if (property == "mmin") {
1418  double mMinTmp;
1419  getWord >> mMinTmp;
1420  pdt[idTmp].setMMin(mMinTmp);
1421  return true;
1422  }
1423  if (property == "mmax") {
1424  double mMaxTmp;
1425  getWord >> mMaxTmp;
1426  pdt[idTmp].setMMax(mMaxTmp);
1427  return true;
1428  }
1429  if (property == "tau0") {
1430  double tau0Tmp;
1431  getWord >> tau0Tmp;
1432  pdt[idTmp].setTau0(tau0Tmp);
1433  return true;
1434  }
1435  if (property == "isresonance") {
1436  string isresTmp;
1437  getWord >> isresTmp;
1438  bool isResonanceTmp = boolString(isresTmp);
1439  pdt[idTmp].setIsResonance(isResonanceTmp);
1440  return true;
1441  }
1442  if (property == "maydecay") {
1443  string mayTmp;
1444  getWord >> mayTmp;
1445  bool mayDecayTmp = boolString(mayTmp);
1446  pdt[idTmp].setMayDecay(mayDecayTmp);
1447  return true;
1448  }
1449  if (property == "taucalc") {
1450  string tauTmp;
1451  getWord >> tauTmp;
1452  bool tauCalcTmp = boolString(tauTmp);
1453  pdt[idTmp].setTauCalc(tauCalcTmp);
1454  return true;
1455  }
1456  if (property == "doexternaldecay") {
1457  string extdecTmp;
1458  getWord >> extdecTmp;
1459  bool doExternalDecayTmp = boolString(extdecTmp);
1460  pdt[idTmp].setDoExternalDecay(doExternalDecayTmp);
1461  return true;
1462  }
1463  if (property == "isvisible") {
1464  string isvisTmp;
1465  getWord >> isvisTmp;
1466  bool isVisibleTmp = boolString(isvisTmp);
1467  pdt[idTmp].setIsVisible(isVisibleTmp);
1468  return true;
1469  }
1470  if (property == "doforcewidth") {
1471  string doforceTmp;
1472  getWord >> doforceTmp;
1473  bool doForceWidthTmp = boolString(doforceTmp);
1474  pdt[idTmp].setDoForceWidth(doForceWidthTmp);
1475  return true;
1476  }
1477  if (property == "varwidth") {
1478  string varwidthTmp;
1479  getWord >> varwidthTmp;
1480  bool varWidthTmp = boolString(varwidthTmp);
1481  pdt[idTmp].setVarWidth(varWidthTmp);
1482  return true;
1483  }
1484 
1485  // Addition or complete replacement of a particle.
1486  if (property == "all" || property == "new") {
1487 
1488  // Default values for properties to be read.
1489  string nameTmp = "void";
1490  string antiNameTmp = "void";
1491  int spinTypeTmp = 0;
1492  int chargeTypeTmp = 0;
1493  int colTypeTmp = 0;
1494  double m0Tmp = 0.;
1495  double mWidthTmp = 0.;
1496  double mMinTmp = 0.;
1497  double mMaxTmp = 0.;
1498  double tau0Tmp = 0.;
1499  bool varWidthTmp = false;
1500 
1501  // Read in data from stream.
1502  getWord >> nameTmp >> antiNameTmp >> spinTypeTmp >> chargeTypeTmp
1503  >> colTypeTmp >> m0Tmp >> mWidthTmp >> mMinTmp >> mMaxTmp
1504  >> tau0Tmp >> varWidthTmp;
1505 
1506  // To keep existing decay channels, only overwrite particle data.
1507  if (property == "all" && isParticle(idTmp)) {
1508  setAll( idTmp, nameTmp, antiNameTmp, spinTypeTmp, chargeTypeTmp,
1509  colTypeTmp, m0Tmp, mWidthTmp, mMinTmp, mMaxTmp, tau0Tmp, varWidthTmp);
1510 
1511  // Else start over completely from scratch.
1512  } else {
1513  if (isParticle(idTmp)) pdt.erase(idTmp);
1514  addParticle( idTmp, nameTmp, antiNameTmp, spinTypeTmp, chargeTypeTmp,
1515  colTypeTmp, m0Tmp, mWidthTmp, mMinTmp, mMaxTmp, tau0Tmp, varWidthTmp);
1516  }
1517  return true;
1518  }
1519 
1520  // Set onMode of all decay channels in one go.
1521  if (property == "onmode") {
1522  int onMode = 0;
1523  string onModeIn;
1524  getWord >> onModeIn;
1525  // For onMode allow the optional possibility of Bool input.
1526  if (isdigit(onModeIn[0])) {
1527  istringstream getOnMode(onModeIn);
1528  getOnMode >> onMode;
1529  } else onMode = (boolString(onModeIn)) ? 1 : 0;
1530  for (int i = 0; i < pdt[idTmp].sizeChannels(); ++i)
1531  pdt[idTmp].channel(i).onMode(onMode);
1532  return true;
1533  }
1534 
1535  // Selective search for matching decay channels.
1536  int matchKind = 0;
1537  if (property == "offifany" || property == "onifany" ||
1538  property == "onposifany" || property == "onnegifany") matchKind = 1;
1539  if (property == "offifall" || property == "onifall" ||
1540  property == "onposifall" || property == "onnegifall") matchKind = 2;
1541  if (property == "offifmatch" || property == "onifmatch" ||
1542  property == "onposifmatch" || property == "onnegifmatch") matchKind = 3;
1543  if (matchKind > 0) {
1544  int onMode = 0;
1545  if (property == "onifany" || property == "onifall"
1546  || property == "onifmatch") onMode = 1;
1547  if (property == "onposifany" || property == "onposifall"
1548  || property == "onposifmatch") onMode = 2;
1549  if (property == "onnegifany" || property == "onnegifall"
1550  || property == "onnegifmatch") onMode = 3;
1551 
1552  // Read in particles to match.
1553  vector<int> idToMatch;
1554  int idRead;
1555  for ( ; ; ) {
1556  getWord >> idRead;
1557  if (!getWord) break;
1558  idToMatch.push_back(abs(idRead));
1559  }
1560  int nToMatch = idToMatch.size();
1561 
1562  // Loop over all decay channels.
1563  for (int i = 0; i < pdt[idTmp].sizeChannels(); ++i) {
1564  int multi = pdt[idTmp].channel(i).multiplicity();
1565 
1566  // Look for any match at all.
1567  if (matchKind == 1) {
1568  bool foundMatch = false;
1569  for (int j = 0; j < multi; ++j) {
1570  int idNow = abs(pdt[idTmp].channel(i).product(j));
1571  for (int k = 0; k < nToMatch; ++k)
1572  if (idNow == idToMatch[k]) {foundMatch = true; break;}
1573  if (foundMatch) break;
1574  }
1575  if (foundMatch) pdt[idTmp].channel(i).onMode(onMode);
1576 
1577  // Look for match of all products provided.
1578  } else {
1579  int nUnmatched = nToMatch;
1580  if (multi < nToMatch);
1581  else if (multi > nToMatch && matchKind == 3);
1582  else {
1583  vector<int> idUnmatched;
1584  for (int k = 0; k < nToMatch; ++k)
1585  idUnmatched.push_back(idToMatch[k]);
1586  for (int j = 0; j < multi; ++j) {
1587  int idNow = abs(pdt[idTmp].channel(i).product(j));
1588  for (int k = 0; k < nUnmatched; ++k)
1589  if (idNow == idUnmatched[k]) {
1590  idUnmatched[k] = idUnmatched[--nUnmatched];
1591  break;
1592  }
1593  if (nUnmatched == 0) break;
1594  }
1595  }
1596  if (nUnmatched == 0) pdt[idTmp].channel(i).onMode(onMode);
1597  }
1598  }
1599  return true;
1600  }
1601 
1602  // Rescale all branching ratios by common factor.
1603  if (property == "rescalebr") {
1604  double factor;
1605  getWord >> factor;
1606  pdt[idTmp].rescaleBR(factor);
1607  return true;
1608  }
1609 
1610  // Reset decay table in preparation for new input.
1611  if (property == "onechannel") pdt[idTmp].clearChannels();
1612 
1613  // Add or change a decay channel: get channel number and new property.
1614  if (property == "addchannel" || property == "onechannel"
1615  || isdigit(property[0])) {
1616  int channel;
1617  if (property == "addchannel" || property == "onechannel") {
1618  pdt[idTmp].addChannel();
1619  channel = pdt[idTmp].sizeChannels() - 1;
1620  property = "all";
1621  } else{
1622  istringstream getChannel(property);
1623  getChannel >> channel;
1624  getWord >> property;
1625  toLowerRep(property);
1626  }
1627 
1628  // Check that channel exists.
1629  if (channel < 0 || channel >= pdt[idTmp].sizeChannels()) return false;
1630 
1631  // Find decay channel property and value, case by case.
1632  // At same time also do case where all should be replaced.
1633  if (property == "onmode" || property == "all") {
1634  int onMode = 0;
1635  string onModeIn;
1636  getWord >> onModeIn;
1637  // For onMode allow the optional possibility of Bool input.
1638  if (isdigit(onModeIn[0])) {
1639  istringstream getOnMode(onModeIn);
1640  getOnMode >> onMode;
1641  } else onMode = (boolString(onModeIn)) ? 1 : 0;
1642  pdt[idTmp].channel(channel).onMode(onMode);
1643  if (property == "onmode") return true;
1644  }
1645  if (property == "bratio" || property == "all") {
1646  double bRatio;
1647  getWord >> bRatio;
1648  pdt[idTmp].channel(channel).bRatio(bRatio);
1649  if (property == "bratio") return true;
1650  }
1651  if (property == "memode" || property == "all") {
1652  int meMode;
1653  getWord >> meMode;
1654  pdt[idTmp].channel(channel).meMode(meMode);
1655  if (property == "memode") return true;
1656  }
1657 
1658  // Scan for products until end of line.
1659  if (property == "products" || property == "all") {
1660  int nProd = 0;
1661  for (int iProd = 0; iProd < 8; ++iProd) {
1662  int idProd;
1663  getWord >> idProd;
1664  if (!getWord) break;
1665  pdt[idTmp].channel(channel).product(iProd, idProd);
1666  ++nProd;
1667  }
1668  for (int iProd = nProd; iProd < 8; ++iProd)
1669  pdt[idTmp].channel(channel).product(iProd, 0);
1670  pdt[idTmp].channel(channel).multiplicity(nProd);
1671  return true;
1672  }
1673 
1674  // Rescale an existing branching ratio.
1675  if (property == "rescalebr") {
1676  double factor;
1677  getWord >> factor;
1678  pdt[idTmp].channel(channel).rescaleBR(factor);
1679  return true;
1680  }
1681  }
1682 
1683  // Return false if failed to recognize property.
1684  if (warn) cout << "\n PYTHIA Error: input property not found in Particle"
1685  << " Data Table:\n " << lineIn << "\n";
1686  readingFailedSave = true;
1687  return false;
1688 
1689 }
1690 
1691 //--------------------------------------------------------------------------
1692 
1693 // Print out complete or changed table of database in numerical order.
1694 
1695 void ParticleData::list(ostream& str, bool changedOnly, bool changedRes) {
1696 
1697  // Table header; output for bool as off/on.
1698  if (!changedOnly) {
1699  str << "\n -------- PYTHIA Particle Data Table (complete) --------"
1700  << "------------------------------------------------------------"
1701  << "--------------\n \n";
1702 
1703  } else {
1704  str << "\n -------- PYTHIA Particle Data Table (changed only) ----"
1705  << "------------------------------------------------------------"
1706  << "--------------\n \n";
1707  }
1708  str << " id name antiName spn chg col m0"
1709  << " mWidth mMin mMax tau0 res dec ext "
1710  << "vis wid\n no onMode bRatio meMode products \n";
1711 
1712  // Iterate through the particle data table. Option to skip unchanged.
1713  int nList = 0;
1714  for (map<int, ParticleDataEntry>::iterator pdtEntry
1715  = pdt.begin(); pdtEntry != pdt.end(); ++pdtEntry) {
1716  particlePtr = &pdtEntry->second;
1717  if ( !changedOnly || particlePtr->hasChanged() ||
1718  ( changedRes && particlePtr->getResonancePtr() != 0 ) ) {
1719 
1720  // Pick format for mass and width based on mass value.
1721  double m0Now = particlePtr->m0();
1722  if (m0Now == 0 || (m0Now > 0.1 && m0Now < 1000.))
1723  str << fixed << setprecision(5);
1724  else str << scientific << setprecision(3);
1725 
1726  // Print particle properties.
1727  ++nList;
1728  str << "\n" << setw(8) << particlePtr->id() << " " << left;
1729  if (particlePtr->name(-1) == "void")
1730  str << setw(33) << particlePtr->name() << " ";
1731  else str << setw(16) << particlePtr->name() << " "
1732  << setw(16) << particlePtr->name(-1) << " ";
1733  str << right << setw(2) << particlePtr->spinType() << " "
1734  << setw(2) << particlePtr->chargeType() << " "
1735  << setw(2) << particlePtr->colType() << " "
1736  << setw(10) << particlePtr->m0() << " "
1737  << setw(10) << particlePtr->mWidth() << " "
1738  << setw(10) << particlePtr->mMin() << " "
1739  << setw(10) << particlePtr->mMax() << " "
1740  << scientific << setprecision(5)
1741  << setw(12) << particlePtr->tau0() << " " << setw(2)
1742  << particlePtr->isResonance() << " " << setw(2)
1743  << (particlePtr->mayDecay() && particlePtr->canDecay())
1744  << " " << setw(2) << particlePtr->doExternalDecay() << " "
1745  << setw(2) << particlePtr->isVisible()<< " "
1746  << setw(2) << particlePtr->doForceWidth() << "\n";
1747 
1748  // Loop through the decay channel table for each particle.
1749  if (particlePtr->sizeChannels() > 0) {
1750  for (int i = 0; i < int(particlePtr->sizeChannels()); ++i) {
1751  const DecayChannel& channel = particlePtr->channel(i);
1752  str << " " << setprecision(7) << setw(5) << i
1753  << setw(6) << channel.onMode() << fixed<< setw(12)
1754  << channel.bRatio() << setw(5) << channel.meMode() << " ";
1755  for (int j = 0; j < channel.multiplicity(); ++j)
1756  str << setw(8) << channel.product(j) << " ";
1757  str << "\n";
1758  }
1759  }
1760  }
1761 
1762  }
1763 
1764  // End of loop over database contents.
1765  if (changedOnly && nList == 0) cout << "\n no particle data has been "
1766  << "changed from its default value \n";
1767  cout << "\n -------- End PYTHIA Particle Data Table -----------------"
1768  << "--------------------------------------------------------------"
1769  << "----------\n" << endl;
1770 
1771 }
1772 
1773 //--------------------------------------------------------------------------
1774 
1775 // Print out partial table of database in input order.
1776 
1777 void ParticleData::list(vector<int> idList) {
1778 
1779  // Table header; output for bool as off/on.
1780  cout << "\n -------- PYTHIA Particle Data Table (partial) ---------"
1781  << "------------------------------------------------------------"
1782  << "--------------\n \n";
1783  cout << " id name antiName spn chg col m0"
1784  << " mWidth mMin mMax tau0 res dec ext "
1785  << "vis wid\n no onMode bRatio meMode products \n";
1786 
1787  // Iterate through the given list of input particles.
1788  for (int i = 0; i < int(idList.size()); ++i) {
1789  particlePtr = particleDataEntryPtr(idList[i]);
1790 
1791  // Pick format for mass and width based on mass value.
1792  double m0Now = particlePtr->m0();
1793  if (m0Now == 0 || (m0Now > 0.1 && m0Now < 1000.))
1794  cout << fixed << setprecision(5);
1795  else cout << scientific << setprecision(3);
1796 
1797  // Print particle properties.
1798  cout << "\n" << setw(8) << particlePtr->id() << " " << left;
1799  if (particlePtr->name(-1) == "void")
1800  cout << setw(33) << particlePtr->name() << " ";
1801  else cout << setw(16) << particlePtr->name() << " "
1802  << setw(16) << particlePtr->name(-1) << " ";
1803  cout << right << setw(2) << particlePtr->spinType() << " "
1804  << setw(2) << particlePtr->chargeType() << " "
1805  << setw(2) << particlePtr->colType() << " "
1806  << setw(10) << particlePtr->m0() << " "
1807  << setw(10) << particlePtr->mWidth() << " "
1808  << setw(10) << particlePtr->mMin() << " "
1809  << setw(10) << particlePtr->mMax() << " "
1810  << scientific << setprecision(5)
1811  << setw(12) << particlePtr->tau0() << " " << setw(2)
1812  << particlePtr->isResonance() << " " << setw(2)
1813  << (particlePtr->mayDecay() && particlePtr->canDecay())
1814  << " " << setw(2) << particlePtr->doExternalDecay() << " "
1815  << setw(2) << particlePtr->isVisible() << " "
1816  << setw(2) << particlePtr->doForceWidth() << "\n";
1817 
1818  // Loop through the decay channel table for each particle.
1819  if (particlePtr->sizeChannels() > 0) {
1820  for (int j = 0; j < int(particlePtr->sizeChannels()); ++j) {
1821  const DecayChannel& channel = particlePtr->channel(j);
1822  cout << " " << setprecision(7) << setw(5) << j
1823  << setw(6) << channel.onMode() << fixed << setw(12)
1824  << channel.bRatio() << setw(5) << channel.meMode() << " ";
1825  for (int k = 0; k < channel.multiplicity(); ++k)
1826  cout << setw(8) << channel.product(k) << " ";
1827  cout << "\n";
1828  }
1829  }
1830 
1831  }
1832 
1833  // End of loop over database contents.
1834  cout << "\n -------- End PYTHIA Particle Data Table -----------------"
1835  << "--------------------------------------------------------------"
1836  << "----------\n" << endl;
1837 
1838 }
1839 
1840 //--------------------------------------------------------------------------
1841 
1842 // Check that table makes sense: e.g. consistent names and mass ranges,
1843 // that branching ratios sum to unity, that charge is conserved and
1844 // that phase space is open in each channel.
1845 // verbosity = 0: mimimal amount of checks, e.g. that no channels open.
1846 // = 1: further warning if individual channels closed
1847 // (except for resonances).
1848 // = 2: also print branching-ratio-averaged threshold mass.
1849 // = 11, 12: as 1, 2, but include resonances in detailed checks.
1850 
1851 void ParticleData::checkTable(int verbosity) {
1852 
1853  // Header.
1854  cout << "\n -------- PYTHIA Check of Particle Data Table ------------"
1855  <<"------\n\n";
1856  int nErr = 0;
1857 
1858  // Loop through all particles.
1859  for (map<int, ParticleDataEntry>::iterator pdtEntry
1860  = pdt.begin(); pdtEntry != pdt.end(); ++pdtEntry) {
1861  particlePtr = &pdtEntry->second;
1862 
1863  // Extract some particle properties. Set some flags.
1864  int idNow = particlePtr->id();
1865  bool hasAntiNow = particlePtr->hasAnti();
1866  int spinTypeNow = particlePtr->spinType();
1867  int chargeTypeNow = particlePtr->chargeType();
1868  int baryonTypeNow = particlePtr->baryonNumberType();
1869  double m0Now = particlePtr->m0();
1870  double mMinNow = particlePtr->mMin();
1871  double mMaxNow = particlePtr->mMax();
1872  double mWidthNow = particlePtr->mWidth();
1873  double tau0Now = particlePtr->tau0();
1874  bool isResonanceNow = particlePtr->isResonance();
1875  bool hasPrinted = false;
1876  bool studyCloser = verbosity > 10 || !isResonanceNow;
1877 
1878  // Check that particle name consistent with charge information.
1879  string particleName = particlePtr->name(1);
1880  if (particleName.size() > 16) {
1881  cout << " Warning: particle " << idNow << " has name " << particleName
1882  << " of length " << particleName.size() << "\n";
1883  hasPrinted = true;
1884  ++nErr;
1885  }
1886  int nPos = 0;
1887  int nNeg = 0;
1888  for (int i = 0; i < int(particleName.size()); ++i) {
1889  if (particleName[i] == '+') ++nPos;
1890  if (particleName[i] == '-') ++nNeg;
1891  }
1892  if ( (nPos > 0 && nNeg > 0) || ( nPos + nNeg > 0
1893  && 3 * (nPos - nNeg) != chargeTypeNow )) {
1894  cout << " Warning: particle " << idNow << " has name " << particleName
1895  << " inconsistent with charge type " << chargeTypeNow << "\n";
1896  hasPrinted = true;
1897  ++nErr;
1898  }
1899 
1900  // Check that antiparticle name consistent with charge information.
1901  if (hasAntiNow) {
1902  particleName = particlePtr->name(-1);
1903  if (particleName.size() > 16) {
1904  cout << " Warning: particle " << idNow << " has name " << particleName
1905  << " of length " << particleName.size() << "\n";
1906  hasPrinted = true;
1907  ++nErr;
1908  }
1909  nPos = 0;
1910  nNeg = 0;
1911  for (int i = 0; i < int(particleName.size()); ++i) {
1912  if (particleName[i] == '+') ++nPos;
1913  if (particleName[i] == '-') ++nNeg;
1914  }
1915  if ( (nPos > 0 && nNeg > 0) || ( nPos + nNeg > 0
1916  && 3 * (nPos - nNeg) != -chargeTypeNow )) {
1917  cout << " Warning: particle " << -idNow << " has name "
1918  << particleName << " inconsistent with charge type "
1919  << -chargeTypeNow << "\n";
1920  hasPrinted = true;
1921  ++nErr;
1922  }
1923  }
1924 
1925  // Check that mass, mass range and width are consistent.
1926  if (particlePtr->useBreitWigner()) {
1927  if (mMinNow > m0Now) {
1928  cout << " Error: particle " << idNow << " has mMin "
1929  << fixed << setprecision(5) << mMinNow
1930  << " larger than m0 " << m0Now << "\n";
1931  hasPrinted = true;
1932  ++nErr;
1933  }
1934  if (mMaxNow > mMinNow && mMaxNow < m0Now) {
1935  cout << " Error: particle " << idNow << " has mMax "
1936  << fixed << setprecision(5) << mMaxNow
1937  << " smaller than m0 " << m0Now << "\n";
1938  hasPrinted = true;
1939  ++nErr;
1940  }
1941  if (mMaxNow > mMinNow && mMaxNow - mMinNow < mWidthNow) {
1942  cout << " Warning: particle " << idNow << " has mMax - mMin "
1943  << fixed << setprecision(5) << mMaxNow - mMinNow
1944  << " smaller than mWidth " << mWidthNow << "\n";
1945  hasPrinted = true;
1946  ++nErr;
1947  }
1948  }
1949 
1950  // Check that particle does not both have width and lifetime.
1951  if (mWidthNow > 0. && tau0Now > 0.) {
1952  cout << " Warning: particle " << idNow << " has both nonvanishing width "
1953  << scientific << setprecision(5) << mWidthNow << " and lifetime "
1954  << tau0Now << "\n";
1955  hasPrinted = true;
1956  ++nErr;
1957  }
1958 
1959  // Begin study decay channels.
1960  if (particlePtr->sizeChannels() > 0) {
1961 
1962  // Loop through all decay channels.
1963  double bRatioSum = 0.;
1964  double bRatioPos = 0.;
1965  double bRatioNeg = 0.;
1966  bool hasPosBR = false;
1967  bool hasNegBR = false;
1968  double threshMass = 0.;
1969  bool openChannel = false;
1970  for (int i = 0; i < particlePtr->sizeChannels(); ++i) {
1971 
1972  // Extract channel properties.
1973  int onMode = particlePtr->channel(i).onMode();
1974  double bRatio = particlePtr->channel(i).bRatio();
1975  int meMode = particlePtr->channel(i).meMode();
1976  int mult = particlePtr->channel(i).multiplicity();
1977  int prod[8];
1978  for (int j = 0; j < 8; ++j)
1979  prod[j] = particlePtr->channel(i).product(j);
1980 
1981  // Sum branching ratios. Include off-channels.
1982  if (onMode == 0 || onMode == 1) bRatioSum += bRatio;
1983  else if (onMode == 2) {bRatioPos += bRatio; hasPosBR = true;}
1984  else if (onMode == 3) {bRatioNeg += bRatio; hasNegBR = true;}
1985 
1986  // Error printout when unknown decay product code.
1987  for (int j = 0; j < 8; ++j) {
1988  if ( prod[j] != 0 && !isParticle(prod[j]) ) {
1989  cout << " Error: unknown decay product for " << idNow
1990  << " -> " << prod[j] << "\n";
1991  hasPrinted = true;
1992  ++nErr;
1993  continue;
1994  }
1995  }
1996 
1997  // Error printout when no decay products or 0 interspersed.
1998  int nLast = 0;
1999  for (int j = 0; j < 8; ++j)
2000  if (prod[j] != 0) nLast = j + 1;
2001  if (mult == 0 || mult != nLast) {
2002  cout << " Error: corrupted decay product list for "
2003  << particlePtr->id() << " -> ";
2004  for (int j = 0; j < 8; ++j) cout << prod[j] << " ";
2005  cout << "\n";
2006  hasPrinted = true;
2007  ++nErr;
2008  continue;
2009  }
2010 
2011  // Check charge conservation and open phase space in decay channel.
2012  int chargeTypeSum = -chargeTypeNow;
2013  int baryonTypeSum = -baryonTypeNow;
2014  double avgFinalMass = 0.;
2015  double minFinalMass = 0.;
2016  bool canHandle = true;
2017  for (int j = 0; j < mult; ++j) {
2018  chargeTypeSum += chargeType( prod[j] );
2019  baryonTypeSum += baryonNumberType( prod[j] );
2020  avgFinalMass += m0( prod[j] );
2021  minFinalMass += m0Min( prod[j] );
2022  if (prod[j] == 81 || prod[j] == 82 || prod[j] == 83)
2023  canHandle = false;
2024  }
2025  threshMass += bRatio * avgFinalMass;
2026 
2027  // Error printout when charge or baryon number not conserved.
2028  if (chargeTypeSum != 0 && canHandle) {
2029  cout << " Error: 3*charge changed by " << chargeTypeSum
2030  << " in " << idNow << " -> ";
2031  for (int j = 0; j < mult; ++j) cout << prod[j] << " ";
2032  cout << "\n";
2033  hasPrinted = true;
2034  ++nErr;
2035  continue;
2036  }
2037  if ( baryonTypeSum != 0 && canHandle && particlePtr->isHadron() ) {
2038  cout << " Error: 3*baryon number changed by " << baryonTypeSum
2039  << " in " << idNow << " -> ";
2040  for (int j = 0; j < mult; ++j) cout << prod[j] << " ";
2041  cout << "\n";
2042  hasPrinted = true;
2043  ++nErr;
2044  continue;
2045  }
2046 
2047  // Begin check that some matrix elements are used correctly.
2048  bool correctME = true;
2049 
2050  // Check matrix element mode 0: recommended not into partons.
2051  if (meMode == 0 && !isResonanceNow) {
2052  bool hasPartons = false;
2053  for (int j = 0; j < mult; ++j) {
2054  int idAbs = abs(prod[j]);
2055  if ( idAbs < 10 || idAbs == 21 || idAbs == 81 || idAbs == 82
2056  || idAbs == 83 || (idAbs > 1000 && idAbs < 10000
2057  && (idAbs/10)%10 == 0) ) hasPartons = true;
2058  }
2059  if (hasPartons) correctME = false;
2060  }
2061 
2062  // Check matrix element mode 1: rho/omega -> pi+ pi- pi0.
2063  bool useME1 = ( mult == 3 && spinTypeNow == 3 && idNow > 100
2064  && idNow < 1000
2065  && particlePtr->channel(i).contains(211, -211, 111) );
2066  if ( meMode == 1 && !useME1 ) correctME = false;
2067  if ( meMode != 1 && useME1 ) correctME = false;
2068 
2069  // Check matrix element mode 2: polarization in V -> PS + PS.
2070  bool useME2 = ( mult == 2 && spinTypeNow == 3 && idNow > 100
2071  && idNow < 1000 && spinType(prod[0]) == 1
2072  && spinType(prod[1]) == 1 );
2073  if ( meMode == 2 && !useME2 ) correctME = false;
2074  if ( meMode != 2 && useME2 ) correctME = false;
2075 
2076  // Check matrix element mode 11, 12 and 13: Dalitz decay with
2077  // one or more particles in addition to the lepton pair,
2078  // or double Dalitz decay.
2079  bool useME11 = ( mult == 3 && !isResonanceNow
2080  && (prod[1] == 11 || prod[1] == 13 || prod[1] == 15)
2081  && prod[2] == -prod[1] );
2082  bool useME12 = ( mult > 3 && !isResonanceNow
2083  && (prod[mult - 2] == 11 || prod[mult - 2] == 13
2084  || prod[mult - 2] == 15) && prod[mult - 1] == -prod[mult - 2] );
2085  bool useME13 = ( mult == 4 && !isResonanceNow
2086  && (prod[0] == 11 || prod[0] == 13) && prod[1] == -prod[0]
2087  && (prod[2] == 11 || prod[2] == 13) && prod[3] == -prod[2] );
2088  if (useME13) useME12 = false;
2089  if ( meMode == 11 && !useME11 ) correctME = false;
2090  if ( meMode != 11 && useME11 ) correctME = false;
2091  if ( meMode == 12 && !useME12 ) correctME = false;
2092  if ( meMode != 12 && useME12 ) correctME = false;
2093  if ( meMode == 13 && !useME13 ) correctME = false;
2094  if ( meMode != 13 && useME13 ) correctME = false;
2095 
2096  // Check matrix element mode 21: tau -> nu_tau hadrons.
2097  bool useME21 = (idNow == 15 && mult > 2 && prod[0] == 16
2098  && abs(prod[1]) > 100);
2099  if ( meMode == 21 && !useME21 ) correctME = false;
2100  if ( meMode != 21 && useME21 ) correctME = false;
2101 
2102  // Check matrix element mode 22, but only for semileptonic decay.
2103  // For a -> b c d require types u = 2, ubar = -2, d = 1, dbar = -1.
2104  if ( isLepton(prod[0]) && isLepton(prod[1]) ) {
2105  bool useME22 = false;
2106  int typeA = 0;
2107  int typeB = 0;
2108  int typeC = 0;
2109  if (particlePtr->isLepton()) {
2110  typeA = (idNow > 0) ? 1 + (idNow-1)%2 : -1 - (1-idNow)%2;
2111  } else if (particlePtr->isHadron()) {
2112  int hQ = particlePtr->heaviestQuark();
2113  // Special case: for B_c either bbar or c decays.
2114  if (idNow == 541 && heaviestQuark(prod[2]) == -5) hQ = 4;
2115  typeA = (hQ > 0) ? 1 + (hQ-1)%2 : -1 - (1-hQ)%2;
2116  }
2117  typeB = (prod[0] > 0) ? 1 + (prod[0]-1)%2 : -1 - (1-prod[0])%2;
2118  typeC = (prod[1] > 0) ? 1 + (prod[1]-1)%2 : -1 - (1-prod[1])%2;
2119  // Special cases.
2120  if ( (idNow == 130 || idNow == 310) && typeC * typeA < 0)
2121  typeA = -typeA;
2122  if (mult == 3 && idNow == 2112 && prod[2] == 2212)
2123  typeA = 3 - typeA;
2124  if (mult == 3 && idNow == 3222 && prod[2] == 3122)
2125  typeA = 3 - typeA;
2126  if (mult > 2 && typeC == typeA && typeB * typeC < 0
2127  && (typeB + typeC)%2 != 0) useME22 = true;
2128  if ( meMode == 22 && !useME22 ) correctME = false;
2129  if ( meMode != 22 && useME22 ) correctME = false;
2130  }
2131 
2132  // Check for matrix element mode 31.
2133  if (meMode == 31) {
2134  int nGamma = 0;
2135  for (int j = 0; j < mult; ++j) if (prod[j] == 22) ++nGamma;
2136  if (nGamma != 1) correctME = false;
2137  }
2138 
2139  // Check for unknown mode, or resonance-only mode.
2140  if ( !isResonanceNow && (meMode < 0 || (meMode > 2 && meMode <= 10)
2141  || (meMode > 13 && meMode <= 20) || (meMode > 23 && meMode <= 30)
2142  || (meMode > 31 && meMode <= 41) || meMode == 51 || meMode == 61
2143  || meMode == 71 || (meMode > 80 && meMode <= 90)
2144  || (!particlePtr->isOctetHadron() && meMode > 92) ) )
2145  correctME = false;
2146 
2147  // Print if incorrect matrix element mode.
2148  if ( !correctME ) {
2149  cout << " Warning: meMode " << meMode << " used for "
2150  << idNow << " -> ";
2151  for (int j = 0; j < mult; ++j) cout << prod[j] << " ";
2152  cout << "\n";
2153  hasPrinted = true;
2154  ++nErr;
2155  }
2156 
2157  // Warning printout when no phase space for decay.
2158  if ( studyCloser && verbosity > 0 && canHandle && onMode > 0
2159  && particlePtr->m0Min() - minFinalMass < 0. ) {
2160  if (particlePtr->m0Max() - minFinalMass < 0.)
2161  cout << " Error: decay never possible for ";
2162  else cout << " Warning: decay sometimes not possible for ";
2163  cout << idNow << " -> ";
2164  for (int j = 0; j < mult; ++j) cout << prod[j] << " ";
2165  cout << "\n";
2166  hasPrinted = true;
2167  ++nErr;
2168  }
2169 
2170  // End loop over decay channels.
2171  if (onMode > 0 && bRatio > 0.) openChannel = true;
2172  }
2173 
2174  // Optional printout of threshold.
2175  if (verbosity%10 > 1 && particlePtr->useBreitWigner()) {
2176  threshMass /= bRatioSum;
2177  cout << " Info: particle " << idNow << fixed << setprecision(5)
2178  << " has average mass threshold " << threshMass
2179  << " while mMin is " << mMinNow << "\n";
2180  hasPrinted = true;
2181  }
2182 
2183  // Error printout when no acceptable decay channels found.
2184  if (studyCloser && !openChannel) {
2185  cout << " Error: no acceptable decay channel found for particle "
2186  << idNow << "\n";
2187  hasPrinted = true;
2188  ++nErr;
2189  }
2190 
2191  // Warning printout when branching ratios do not sum to unity.
2192  if (studyCloser && (!hasAntiNow || (!hasPosBR && !hasNegBR))
2193  && abs(bRatioSum + bRatioPos - 1.) > 1e-8) {
2194  cout << " Warning: particle " << idNow << fixed << setprecision(8)
2195  << " has branching ratio sum " << bRatioSum << "\n";
2196  hasPrinted = true;
2197  ++nErr;
2198  } else if (studyCloser && hasAntiNow
2199  && (abs(bRatioSum + bRatioPos - 1.) > 1e-8
2200  || abs(bRatioSum + bRatioNeg - 1.) > 1e-8)) {
2201  cout << " Warning: particle " << idNow << fixed << setprecision(8)
2202  << " has branching ratio sum " << bRatioSum + bRatioPos
2203  << " while its antiparticle has " << bRatioSum + bRatioNeg
2204  << "\n";
2205  hasPrinted = true;
2206  ++nErr;
2207  }
2208 
2209  // End study of decay channels and loop over particles.
2210  }
2211  if (hasPrinted) cout << "\n";
2212  }
2213 
2214  // Final output. Done.
2215  cout << " Total number of errors and warnings is " << nErr << "\n";
2216  cout << "\n -------- End PYTHIA Check of Particle Data Table --------"
2217  << "------\n" << endl;
2218 
2219 }
2220 
2221 //--------------------------------------------------------------------------
2222 
2223 // Return the id of the sequentially next particle stored in table.
2224 
2225 int ParticleData::nextId(int idIn) const {
2226 
2227  // Return 0 for negative or unknown codes. Return first for 0.
2228  if (idIn < 0 || (idIn > 0 && !isParticle(idIn))) return 0;
2229  if (idIn == 0) return pdt.begin()->first;
2230 
2231  // Find pointer to current particle and step up. Return 0 if impossible.
2232  map<int, ParticleDataEntry>::const_iterator pdtIn = pdt.find(idIn);
2233  if (pdtIn == pdt.end()) return 0;
2234  ++pdtIn;
2235  if (pdtIn == pdt.end()) return 0;
2236  return pdtIn->first;
2237 
2238 }
2239 
2240 //--------------------------------------------------------------------------
2241 
2242 // Fractional width associated with open channels of one or two resonances.
2243 
2244 double ParticleData::resOpenFrac(int id1In, int id2In, int id3In) {
2245 
2246  // Default value.
2247  double answer = 1.;
2248 
2249  // First resonance.
2250  if ( ParticleDataEntry* ptr = findParticle(id1In) )
2251  answer = ptr->resOpenFrac(id1In);
2252 
2253  // Possibly second resonance.
2254  if ( ParticleDataEntry* ptr = findParticle(id2In) )
2255  answer *= ptr->resOpenFrac(id2In);
2256 
2257  // Possibly third resonance.
2258  if ( ParticleDataEntry* ptr = findParticle(id3In) )
2259  answer *= ptr->resOpenFrac(id3In);
2260 
2261  // Done.
2262  return answer;
2263 
2264 }
2265 
2266 //--------------------------------------------------------------------------
2267 
2268 // Extract XML value string following XML attribute.
2269 
2270 string ParticleData::attributeValue(string line, string attribute) {
2271 
2272  if (line.find(attribute) == string::npos) return "";
2273  int iBegAttri = line.find(attribute);
2274  int iBegQuote = line.find("\"", iBegAttri + 1);
2275  int iEndQuote = line.find("\"", iBegQuote + 1);
2276  return line.substr(iBegQuote + 1, iEndQuote - iBegQuote - 1);
2277 
2278 }
2279 
2280 //--------------------------------------------------------------------------
2281 
2282 // Extract XML bool value following XML attribute.
2283 
2284 bool ParticleData::boolAttributeValue(string line, string attribute) {
2285 
2286  string valString = attributeValue(line, attribute);
2287  if (valString == "") return false;
2288  return boolString(valString);
2289 
2290 }
2291 
2292 //--------------------------------------------------------------------------
2293 
2294 // Extract XML int value following XML attribute.
2295 
2296 int ParticleData::intAttributeValue(string line, string attribute) {
2297  string valString = attributeValue(line, attribute);
2298  if (valString == "") return 0;
2299  istringstream valStream(valString);
2300  int intVal;
2301  valStream >> intVal;
2302  return intVal;
2303 
2304 }
2305 
2306 //--------------------------------------------------------------------------
2307 
2308 // Extract XML double value following XML attribute.
2309 
2310 double ParticleData::doubleAttributeValue(string line, string attribute) {
2311  string valString = attributeValue(line, attribute);
2312  if (valString == "") return 0.;
2313  istringstream valStream(valString);
2314  double doubleVal;
2315  valStream >> doubleVal;
2316  return doubleVal;
2317 
2318 }
2319 
2320 //==========================================================================
2321 
2322 } // end namespace Pythia8