StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
HadronWidths.h
1 // HadronWidths.h 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 // Header file for computing mass-dependent widths and branching ratios
7 
8 #ifndef Pythia8_HadronWidths_H
9 #define Pythia8_HadronWidths_H
10 
11 #include "Pythia8/MathTools.h"
12 #include "Pythia8/ParticleData.h"
13 #include "Pythia8/PhysicsBase.h"
14 
15 namespace Pythia8 {
16 
17 //==========================================================================
18 
19 // The HadronWidths class is used to compute mass-dependent widths
20 // and branching ratios of hadrons.
21 
22 class HadronWidths : public PhysicsBase {
23 
24 public:
25 
26  HadronWidths() = default;
27  HadronWidths(const HadronWidths&) = delete;
28  HadronWidths(HadronWidths&&) = delete;
29 
30  // Load hadron widths data from an xml file.
31  bool init(string path);
32  bool init(istream& stream);
33 
34  // Check whether input data is valid and matches particle data.
35  bool check();
36 
37  // Get a list of all implemented resonances.
38  vector<int> getResonances() const;
39 
40  // Get whether the specified incoming particles can form a resonance.
41  bool hasResonances(int idA, int idB) const;
42 
43  // Get resonances that can be formed by the specified incoming particles.
44  vector<int> possibleResonances(int idA, int idB) const;
45 
46  // Returns whether the specified particle is handled by HadronWidths.
47  bool hasData(int id) const {
48  auto iter = entries.find(abs(id));
49  return iter != entries.end();
50  }
51 
52  // Get whether the resonance can decay into the specified products.
53  bool canDecay(int id, int prodA, int prodB) const;
54 
55  // Get the total width of the specified particle at the specified mass.
56  double width(int id, double m) const;
57 
58  // Get the partial width for the specified decay channel of the particle.
59  double partialWidth(int idR, int prodA, int prodB, double m) const;
60 
61  // Get the branching ratio for the specified decay channel of the particle.
62  double br(int idR, int prodA, int prodB, double m) const;
63 
64  // Get the mass distribution density for the particle at the specified mass.
65  double mDistr(int id, double m) const;
66 
67  // Sample masses for the outgoing system with a given eCM.
68  bool pickMasses(int idA, int idB, double eCM,
69  double& mAOut, double& mBOut, int lType = 1);
70 
71  // Pick a decay channel for the specified particle, together with phase
72  // space configuration. Returns whether successful.
73  bool pickDecay(int idDec, double m, int& idAOut, int& idBOut,
74  double& mAOut, double& mBOut);
75 
76  // Calculate the total width of the particle without using interpolation.
77  double widthCalc(int id, double m) const;
78 
79  // Calculate partial width of the particle without using interpolation.
80  double widthCalc(int id, int prodA, int prodB, double m) const;
81 
82  // Regenerate parameterization for the specified particle.
83  bool parameterize(int id, int precision);
84 
85  // Regenerate parameterization for all particles.
86  void parameterizeAll(int precision);
87 
88  // Write all width data to an xml file.
89  bool save(ostream& stream) const;
90  bool save(string file = "HadronWidths.dat") const {
91  ofstream stream(file); return save(stream); }
92 
93 private:
94 
95  // Struct for mass dependent partial width and other decay channel info.
96  struct ResonanceDecayChannel {
97  LinearInterpolator partialWidth;
98  int prodA, prodB;
99 
100  // 2l+1, where l is the angular momentum of the outgoing two-body system.
101  int lType;
102 
103  // Minimum mass for this channel.
104  double mThreshold;
105  };
106 
107  // Structure for total width parameterization and map to decay channels.
108  struct HadronWidthEntry {
109  LinearInterpolator width;
110  map<pair<int, int>, ResonanceDecayChannel> decayChannels;
111  };
112 
113  // Map from particle id to corresponding HadronWidthEntry.
114  map<int, HadronWidthEntry> entries;
115 
116  // Gets key for the decay and flips idR if necessary
117  pair<int, int> getKey(int& idR, int idA, int idB) const;
118 
119  // Map from signatures to candidate resonances. Used for optimization.
120  map<int, vector<int>> signatureToParticles;
121 
122  // Get signature of system based on total baryon number and electric charge.
123  int getSignature(int baryonNumber, int charge) const;
124 
125  // Get total available phase space.
126  double psSize(double eCM, ParticleDataEntry& prodA,
127  ParticleDataEntry& prodB, double lType) const;
128 
129  // Calculate partial width of the particle without using interpolation.
130  double widthCalc(int id, DecayChannel& channel, double m) const;
131 
132  // Generate parameterization for particle and its decay products if needed.
133  bool parameterizeRecursive(int id, int precision);
134 
135 };
136 
137 //==========================================================================
138 
139 } // end namespace Pythia8
140 
141 #endif // Pythia8_HadronWidths_H