|
||||
File indexing completed on 2025-01-18 10:06:24
0001 // HadronWidths.h is a part of the PYTHIA event generator. 0002 // Copyright (C) 2024 Torbjorn Sjostrand. 0003 // PYTHIA is licenced under the GNU GPL v2 or later, see COPYING for details. 0004 // Please respect the MCnet Guidelines, see GUIDELINES for details. 0005 0006 // Header file for computing mass-dependent widths and branching ratios 0007 0008 #ifndef Pythia8_HadronWidths_H 0009 #define Pythia8_HadronWidths_H 0010 0011 #include "Pythia8/MathTools.h" 0012 #include "Pythia8/ParticleData.h" 0013 #include "Pythia8/PhysicsBase.h" 0014 0015 namespace Pythia8 { 0016 0017 //========================================================================== 0018 0019 // The HadronWidths class is used to compute mass-dependent widths 0020 // and branching ratios of hadrons. 0021 0022 class HadronWidths : public PhysicsBase { 0023 0024 public: 0025 0026 // Load hadron widths data from an xml file. 0027 bool init(string path); 0028 bool init(istream& stream); 0029 0030 // Check whether input data is valid and matches particle data. 0031 bool check(); 0032 0033 // Get whether the specified incoming particles can form a resonance. 0034 bool hasResonances(int idA, int idB) const; 0035 0036 // Get all implemented resonances. 0037 set<int> getResonances() const; 0038 0039 // Get resonances that can be formed by the specified incoming particles. 0040 set<int> getResonances(int idA, int idB) const; 0041 0042 // Returns whether the specified particle is handled by HadronWidths. 0043 bool hasData(int id) const { 0044 auto iter = entries.find(abs(id)); 0045 return iter != entries.end(); 0046 } 0047 0048 // Get whether the resonance can decay into the specified products. 0049 bool canDecay(int id, int prodA, int prodB) const; 0050 0051 // Get the total width of the specified particle at the specified mass. 0052 double width(int id, double m) const; 0053 0054 // Get the partial width for the specified decay channel of the particle. 0055 double partialWidth(int idR, int prodA, int prodB, double m) const; 0056 0057 // Get the branching ratio for the specified decay channel of the particle. 0058 double br(int idR, int prodA, int prodB, double m) const; 0059 0060 // Get the mass distribution density for the particle at the specified mass. 0061 double mDistr(int id, double m) const; 0062 0063 // Sample masses for the outgoing system with a given eCM. 0064 bool pickMasses(int idA, int idB, double eCM, 0065 double& mAOut, double& mBOut, int lType = 1); 0066 0067 // Pick a decay channel for the specified particle, together with phase 0068 // space configuration. Returns whether successful. 0069 bool pickDecay(int idDec, double m, int& idAOut, int& idBOut, 0070 double& mAOut, double& mBOut); 0071 0072 // Calculate the total width of the particle without using interpolation. 0073 double widthCalc(int id, double m) const; 0074 0075 // Calculate partial width of the particle without using interpolation. 0076 double widthCalc(int id, int prodA, int prodB, double m) const; 0077 0078 // Regenerate parameterization for particle, using the specified number of 0079 // interpolation points. If needed, its decay products are automatically 0080 // parameterized as well. 0081 bool parameterize(int id, int precision = 50); 0082 0083 // Regenerate parameterization for all particles. 0084 void parameterizeAll(int precision = 50); 0085 0086 // Write all width data to an xml file. 0087 bool save(ostream& stream) const; 0088 bool save(string file = "HadronWidths.dat") const { 0089 ofstream stream(file); return save(stream); } 0090 0091 private: 0092 0093 // Struct for mass dependent partial width and other decay channel info. 0094 struct ResonanceDecayChannel { 0095 LinearInterpolator partialWidth; 0096 int prodA, prodB; 0097 0098 // 2l+1, where l is the angular momentum of the outgoing two-body system. 0099 int lType; 0100 0101 // Minimum mass for this channel. 0102 double mThreshold; 0103 }; 0104 0105 // Structure for total width parameterization and map to decay channels. 0106 struct HadronWidthEntry { 0107 LinearInterpolator width; 0108 map<pair<int, int>, ResonanceDecayChannel> decayChannels; 0109 bool isUserDefined; 0110 }; 0111 0112 // Map from particle id to corresponding HadronWidthEntry. 0113 map<int, HadronWidthEntry> entries; 0114 0115 // Gets key for the decay and flips idR if necessary 0116 pair<int, int> getKey(int& idR, int idA, int idB) const; 0117 0118 // Map from signatures to candidate resonances. Used for optimization. 0119 map<int, vector<int>> signatureToParticles; 0120 0121 // Get signature of system based on total baryon number and electric charge. 0122 int getSignature(int baryonNumber, int charge) const; 0123 0124 // Get total available phase space. 0125 double psSize(double eCM, ParticleDataEntryPtr prodA, 0126 ParticleDataEntryPtr prodB, double lType) const; 0127 0128 // Calculate partial width of the particle without using interpolation. 0129 double widthCalc(int id, DecayChannel& channel, double m) const; 0130 0131 // Recursive method for parameterization. 0132 bool parameterizeRecursive(int id, int precision); 0133 0134 }; 0135 0136 //========================================================================== 0137 0138 } // end namespace Pythia8 0139 0140 #endif // Pythia8_HadronWidths_H
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |