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