dynamics#
import ampform.dynamics
Lineshape functions that describe the dynamics of an interaction.
See also
- class BlattWeisskopfSquared(z, angular_momentum, *args, evaluate: bool = False, **kwargs)[source]#
Bases:
ExprBlatt-Weisskopf function \(B_L^2(z)\), up to \(L \leq 8\).
- Parameters:
z – Argument of the Blatt-Weisskopf function \(B_L^2(z)\). A usual choice is \(z = (d q)^2\) with \(d\) the impact parameter and \(q\) the breakup-momentum (see
BreakupMomentumSquared).angular_momentum – Angular momentum \(L\) of the decaying particle.
Note that equal powers of \(z\) appear in the nominator and the denominator, while some sources have nominator \(1\), instead of \(z^L\). Compare for instance Equation (50.27) in PDG2021, §Resonances, p.9.
Each of these cases for \(L\) has been taken from [12], p.59, [2], p.415, and [13]. For a good overview of where to use these Blatt-Weisskopf functions, see [14].
See also Form factor.
(1)#\[\begin{split} \begin{array}{rcl} B_{L}^2\left(z\right) &=& \begin{cases} 1 & \text{for}\: L = 0 \\\frac{2 z}{z + 1} & \text{for}\: L = 1 \\\frac{13 z^{2}}{9 z + \left(z - 3\right)^{2}} & \text{for}\: L = 2 \\\frac{277 z^{3}}{z \left(z - 15\right)^{2} + \left(2 z - 5\right) \left(18 z - 45\right)} & \text{for}\: L = 3 \\\frac{12746 z^{4}}{25 z \left(2 z - 21\right)^{2} + \left(z^{2} - 45 z + 105\right)^{2}} & \text{for}\: L = 4 \\\frac{998881 z^{5}}{z^{5} + 15 z^{4} + 315 z^{3} + 6300 z^{2} + 99225 z + 893025} & \text{for}\: L = 5 \\\frac{118394977 z^{6}}{z^{6} + 21 z^{5} + 630 z^{4} + 18900 z^{3} + 496125 z^{2} + 9823275 z + 108056025} & \text{for}\: L = 6 \\\frac{19727003738 z^{7}}{z^{7} + 28 z^{6} + 1134 z^{5} + 47250 z^{4} + 1819125 z^{3} + 58939650 z^{2} + 1404728325 z + 18261468225} & \text{for}\: L = 7 \\\frac{4392846440677 z^{8}}{z^{8} + 36 z^{7} + 1890 z^{6} + 103950 z^{5} + 5457375 z^{4} + 255405150 z^{3} + 9833098275 z^{2} + 273922023375 z + 4108830350625} & \text{for}\: L = 8 \end{cases} \\ \end{array}\end{split}\]
- class EnergyDependentWidth(s, mass0, gamma0, m_a, m_b, angular_momentum, meson_radius, *args, evaluate: bool = False, **kwargs)[source]#
Bases:
ExprMass-dependent width, coupled to the pole position of the resonance.
See Equation (50.28) in PDG2021, §Resonances, p.9 and [14], equation (6). Default value for
phsp_factorisPhaseSpaceFactor.Note that the
BlattWeisskopfSquaredof AmpForm is normalized in the sense that equal powers of \(z\) appear in the nominator and the denominator, while the definition in the PDG (as well as some other sources), always have \(1\) in the nominator of the Blatt-Weisskopf. In that case, one needs an additional factor \(\left(q/q_0\right)^{2L}\) in the definition for \(\Gamma(m)\).With that in mind, the “mass-dependent” width in a
relativistic_breit_wigner_with_ffbecomes:(2)#\[\begin{split} \begin{array}{rcl} \Gamma_{0}\left(s\right) &=& \frac{\Gamma_{0} B_{L}^2\left(q^2\left(s\right)\right) \rho\left(s\right)}{B_{L}^2\left(q^2_{0}\left(m_{0}^{2}\right)\right) \rho_{0}\left(m_{0}^{2}\right)} \\ \end{array}\end{split}\]where \(B_L^2\) is defined by (1), \(q\) is defined by (1), and \(\rho\) is (by default) defined by (2).
- phsp_factor: PhaseSpaceFactorProtocol[source]#
- relativistic_breit_wigner(s, mass0, gamma0) Expr[source]#
Relativistic Breit-Wigner lineshape.
See Without form factor and [14].
(3)#\[\frac{\Gamma_{0} m_{0}}{- i \Gamma_{0} m_{0} + m_{0}^{2} - s}\]
- relativistic_breit_wigner_with_ff(s, mass0, gamma0, m_a, m_b, angular_momentum, meson_radius, phsp_factor: ~ampform.dynamics.phasespace.PhaseSpaceFactorProtocol = <class 'ampform.dynamics.phasespace.PhaseSpaceFactor'>) Expr[source]#
Relativistic Breit-Wigner with
BlattWeisskopfSquaredfactor.See With form factor and PDG2021, §Resonances, p.9.
The general form of a relativistic Breit-Wigner with Blatt-Weisskopf form factor is:
(4)#\[\frac{\Gamma_{0} m_{0} \sqrt{B_{L}^2\left(d^{2} q^2\left(s\right)\right)}}{m_{0}^{2} - i m_{0} \Gamma_{0}\left(s\right) - s}\]where \(\Gamma(s)\) is defined by (2), \(B_L^2\) is defined by (1), and \(q^2\) is defined by (1).
- formulate_form_factor(s, m_a, m_b, angular_momentum, meson_radius) Expr[source]#
Formulate a Blatt-Weisskopf form factor.
Returns the production process factor \(n_a\) from Equation (50.26) in PDG2021, §Resonances, p.9, which features the
sqrtof aBlattWeisskopfSquared.(5)#\[\sqrt{B_{L}^2\left(d^{2} q^2\left(s\right)\right)}\]
Submodules and Subpackages
- builder
- kmatrix
- phasespace