core

Basic Amplitude Calculations. A partial wave analysis process has following structure:

DecayGroup: addition (+)
DecayChain: multiplication (x)

Decay, Particle(Propagator)

class AmpBase[source]

Bases: object

Base class for amplitude

add_var(names, is_complex=False, shape=(), **kwargs)[source]

default add_var method

amp_shape()[source]
get_factor_variable()[source]
get_params_head()[source]
get_var(name)[source]
get_variable_name(name='')[source]
class AmpDecay(core, outs, name=None, disable=False, p_break=False, c_break=True, curve_style=None, **kwargs)[source]

Bases: Decay, AmpBase

base class for decay with amplitude

amp_index(base_map)[source]
amp_shape()[source]
get_params_head()[source]
list_helicity_inner()[source]
n_helicity_inner()[source]
class AmpDecayChain(*args, is_cp=False, aligned=True, **kwargs)[source]

Bases: DecayChain, AmpBase

get_params_head()[source]
class AngSam3Decay(core, outs, name=None, disable=False, p_break=False, c_break=True, curve_style=None, **kwargs)[source]

Bases: AmpDecay, AmpBase

get_amp(data, data_extra=None, **kwargs)[source]
init_params()[source]
model_name = 'default'
class DecayChain(*args, is_cp=False, aligned=True, **kwargs)[source]

Bases: AmpDecayChain

A list of Decay as a chain decay

amp_index(base_map=None)[source]
amp_shape()[source]
factor_iteration(deep=1)[source]
get_all_factor()[source]
get_amp(data_c, data_p, all_data=None, base_map=None)[source]
get_amp_particle(data_p, data_c, all_data=None)[source]
get_amp_total(charge=1)[source]
get_angle_amp(data_c, data_p, all_data=None, base_map=None)[source]
get_base_map(base_map=None)[source]
get_cp_amp_total(charge=1)[source]
get_factor()[source]
get_factor_angle_amp(data_c, data_p, all_data=None, base_map=None)[source]
get_factor_variable()[source]
get_m_dep(data_c, data_p, all_data=None, base_map=None)[source]
init_params(name='')[source]
model_name = 'default'
product_gls()[source]
class DecayGroup(chains)[source]

Bases: DecayGroup, AmpBase

A Group of Decay Chains with the same final particles.

add_used_chains(used_chains)[source]
amp_index(gen=None, base_map=None)[source]
chains_particle()[source]
factor_iteration(deep=2)[source]
generate_phasespace(num=100000)[source]
get_amp(data)[source]

calculate the amplitude as complex number

get_amp2(data)[source]
get_amp3(data)[source]
get_angle_amp(data)[source]
get_base_map(gen=None, base_map=None)[source]
get_density_matrix()[source]
get_factor()[source]
get_factor_angle_amp(data)[source]
get_factor_variable()[source]
get_id_swap_transpose(key, n)[source]
get_m_dep(data)[source]

get mass dependent items

get_res_map()[source]
get_swap_factor(key)[source]
get_swap_transpose(trans, n)[source]
init_params(name='')[source]
partial_weight(data, combine=None)[source]
partial_weight_interference(data)[source]
set_used_chains(used_chains)[source]
set_used_res(res, only=False)[source]
sum_amp(data, cached=True)[source]

calculat the amplitude modular square

sum_amp_polarization(data)[source]

sum amplitude suqare with density _get_cg_matrix

\[P = \sum_{m, m', \cdots } A_{m, \cdots} \rho_{m, m'} A^{*}_{m', \cdots}\]
sum_with_polarization(amp)[source]
temp_used_res(res)[source]
class FloatParams(x=0, /)[source]

Bases: float

class HelicityDecay(*args, has_barrier_factor=True, l_list=None, barrier_factor_mass=False, has_ql=True, has_bprime=True, aligned=False, allow_cc=True, ls_list=None, barrier_factor_norm=False, params_polar=None, below_threshold=False, force_min_l=False, params_head=None, no_q0=False, helicity_inner_full=False, ls_selector=None, **kwargs)[source]

Bases: AmpDecay

default decay model

The total amplitude is

\[A = H_{\lambda_{B},\lambda_{C}}^{A \rightarrow B+C} D^{J_A*}_{\lambda_{A}, \lambda_{B}-\lambda_{C}} (\varphi,\theta,0)\]

The helicity coupling is

\[H_{\lambda_{B},\lambda_{C}}^{A \rightarrow B+C} = \sum_{ls} g_{ls} \sqrt{\frac{2l+1}{2 J_{A}+1}} \langle l 0; s \delta|J_{A} \delta\rangle \langle J_{B} \lambda_{B} ;J_{C} -\lambda_{C} | s \delta \rangle q^{l} B_{l}'(q, q_0, d)\]

The fit parameters is \(g_{ls}\)

There are some options

(1). has_bprime=False will remove the \(B_{l}'(q, q_0, d)\) part.

(2). has_barrier_factor=False will remove the \(q^{l} B_{l}'(q, q_0, d)\) part.

(3). barrier_factor_norm=True will replace \(q^l\) with \((q/q_{0})^l\)

(4). below_threshold=True will replace the mass used to calculate \(q_0\) with

\[m_0^{eff} = m^{min} + \frac{m^{max} - m^{min}}{2}(1+tanh \frac{m_0 - \frac{m^{max} + m^{min}}{2}}{m^{max} - m^{min}})\]

(5). l_list=[l1, l2] and ls_list=[[l1, s1], [l2, s2]] options give the list of all possible LS used in the decay.

(6). no_q0=True will set the \(q_0=1\).

add_algin(ret, data)[source]
build_ls2hel_eq()[source]
build_simple_data()[source]
check_valid_jp()[source]
factor_iter_names(deep=1, extra=[])[source]
get_amp(data, data_p, **kwargs)[source]
get_angle_amp(data, data_p, **kwargs)[source]
get_angle_g_ls()[source]
get_angle_helicity_amp(data, data_p, **kwargs)[source]
get_angle_ls_amp(data, data_p, **kwargs)[source]
get_barrier_factor(mass, q, q0, d)[source]
get_barrier_factor2(mass, q2, q02, d)[source]
get_barrier_factor_mass(mass)[source]
get_cg_matrix(out_sym=False)[source]

The matrix indexed by \([(l,s),(\lambda_b,\lambda_c)]\). The matrix element is

\[\sqrt{\frac{ 2 l + 1 }{ 2 j_a + 1 }} \langle j_b, j_c, \lambda_b, - \lambda_c | s, \lambda_b - \lambda_c \rangle \langle l, s, 0, \lambda_b - \lambda_c | j_a, \lambda_b - \lambda_c \rangle\]

This is actually the pre-factor of \(g_ls\) in the amplitude formula.

Returns:

2-d array of real numbers

get_factor()[source]
get_factor_H(data, data_p, **kwargs)[source]
get_factor_angle_amp(data, data_p, **kwargs)[source]
get_factor_angle_helicity_amp(data, data_p, **kwargs)[source]
get_factor_m_dep(data, data_p, **kwargs)[source]
get_factor_variable()[source]
get_g_ls()[source]
get_helicity_amp(data, data_p, **kwargs)[source]
get_ls_amp(data, data_p, **kwargs)[source]
get_ls_amp_org(data, data_p, **kwargs)[source]
get_ls_list()[source]

get possible ls for decay, with l_list filter possible l

get_m_dep(data, data_p, **kwargs)[source]
get_params_head()[source]
get_relative_momentum(data, from_data=False)[source]
get_relative_momentum2(data, from_data=False)[source]
get_total_ls_list()[source]
init_params()[source]
mask_factor_vars()[source]
model_name = 'default'
set_ls(ls)[source]
class Particle(*args, running_width=True, bw_l=None, width_norm=False, params_head=None, **kwargs)[source]

Bases: BaseParticle, AmpBase

\[R(m) = \frac{1}{m_0^2 - m^2 - i m_0 \Gamma(m)}\]

Argand diagram

(Source code, png, hires.png, pdf)

../_images/tf_pwa-amp-core-1.png

Pole position

(Source code, png, hires.png, pdf)

../_images/tf_pwa-amp-core-2.png
amp_shape()[source]
get_amp(data, data_c, **kwargs)[source]
get_factor()[source]
get_mass()[source]
get_num_var()[source]
get_subdecay_mass(idx=0)[source]
get_sympy_dom(m, m0, g0, m1=None, m2=None, sheet=0)[source]
get_sympy_var()[source]
get_width()[source]
init_params()[source]
is_fixed_shape()[source]
model_name = 'BWR'
pole_function(sheet=0, modules='numpy')[source]
solve_pole(init=None, sheet=0, return_complex=True)[source]
class ParticleX(*args, running_width=True, bw_l=None, width_norm=False, params_head=None, **kwargs)[source]

Bases: Particle

simple particle model for mass, (used in expr)

\[R(m) = m\]

(Source code, png, hires.png, pdf)

../_images/tf_pwa-amp-core-3.png
get_amp(data, *args, **kwargs)[source]
model_name = 'x'
class SimpleResonances(*args, **kwargs)[source]

Bases: Particle

get_amp(*args, **kwargs)[source]
data_device(data)[source]
get_decay(core, outs, **kwargs)[source]

method for getting decay of model

get_decay_chain(decays, **kwargs)[source]

method for getting decay of model

get_decay_model(model, num_outs=2)[source]
get_name(self, names)[source]
get_particle(*args, model='default', **kwargs)[source]

method for getting particle of model

get_particle_model(name)[source]
get_particle_model_name(p)[source]
get_relative_p(m_0, m_1, m_2)[source]

relative momentum for 0 -> 1 + 2

get_relative_p2(m_0, m_1, m_2)[source]

relative momentum for 0 -> 1 + 2

index_generator(base_map=None)[source]
load_decfile_particle(fname)[source]
ls_selector_qr(decay, ls_list)[source]
regist_decay(name=None, num_outs=2, f=None)

register a decay model

Params name:

model name used in configuration

Params f:

Model class

regist_particle(name=None, f=None)

register a particle model

Params name:

model name used in configuration

Params f:

Model class

register_decay(name=None, num_outs=2, f=None)[source]

register a decay model

Params name:

model name used in configuration

Params f:

Model class

register_decay_chain(name=None, f=None)[source]

register a decay model

Params name:

model name used in configuration

Params f:

Model class

register_particle(name=None, f=None)[source]

register a particle model

Params name:

model name used in configuration

Params f:

Model class

rename_data_dict(data, idx_map)[source]
simple_cache_fun(f)[source]
simple_deepcopy(dic)[source]
simple_resonance(name, fun=None, params=None)[source]

convert simple fun f(m) into a resonances model

Params name:

model name used in configuration

Params fun:

Model function

Params params:

arguments name list for parameters

trans_model(model)[source]
value_and_grad(f, var)[source]
variable_scope(vm=None)[source]

variabel name scope