cal_angle

This module provides functions which are useful when calculating the angular variables.

The full data structure provided is

{
    "particle": {
        A: {"p": ..., "m": ...},
        (C, D): {"p": ..., "m": ...},
        ...
    },

    "decay":{
        [A->(C, D)+B, (C, D)->C+D]:
        {
            A->(C, D)+B: {
                (C, D): {
                    "ang": {
                        "alpha":[...],
                        "beta": [...],
                        "gamma": [...]
                    },
                    "z": [[x1,y1,z1],...],
                    "x": [[x2,y2,z2],...]
                },
                B: {...}
            },

            (C, D)->C+D: {
                C: {
                    ...,
                    "aligned_angle": {
                    "alpha": [...],
                    "beta": [...],
                    "gamma": [...]
                }
            },
                D: {...}
            },
            A->(B, D)+C: {...},
            (B, D)->B+D: {...}
        },
    ...
    }
}

Inner nodes are named as tuple of particles.

class CalAngleData[source]

Bases: dict

get_angle(decay, p)[source]

get hilicity angle of decay which product particle p

get_decay()[source]
get_mass(name)[source]
get_momentum(name)[source]
get_weight()[source]
mass_hist(name, bins='sqrt', **kwargs)[source]
savetxt(file_name, order=None, cp_trans=False, save_charge=False)[source]
Getp(M_0, M_1, M_2)[source]

Consider a two-body decay \(M_0\rightarrow M_1M_2\). In the rest frame of \(M_0\), the momentum of \(M_1\) and \(M_2\) are definite.

Parameters:
  • M_0 – The invariant mass of \(M_0\)

  • M_1 – The invariant mass of \(M_1\)

  • M_2 – The invariant mass of \(M_2\)

Returns:

the momentum of \(M_1\) (or \(M_2\))

Getp2(M_0, M_1, M_2)[source]

Consider a two-body decay \(M_0\rightarrow M_1M_2\). In the rest frame of \(M_0\), the momentum of \(M_1\) and \(M_2\) are definite.

Parameters:
  • M_0 – The invariant mass of \(M_0\)

  • M_1 – The invariant mass of \(M_1\)

  • M_2 – The invariant mass of \(M_2\)

Returns:

the momentum of \(M_1\) (or \(M_2\))

add_mass(data: dict, _decay_chain: DecayChain | None = None) dict[source]
add particles mass array for data momentum.

{top:{p:momentum},inner:{p:..},outs:{p:..}} => {top:{p:momentum,m:mass},…}

add_relative_momentum(data: dict)[source]

add add rest frame momentum scalar from data momentum.

from {"particle": {A: {"m": ...}, ...}, "decay": {A->B+C: {...}, ...}

to {"particle": {A: {"m": ...}, ...},"decay": {[A->B+C,...]: {A->B+C:{...,"|q|": ...},...},...}

add_weight(data: dict, weight: float = 1.0) dict[source]
add inner data weights for data.

{…} => {…, “weight”: weights}

aligned_angle_ref_rule1(decay_group, decay_chain_struct, decay_data, data)[source]
aligned_angle_ref_rule2(decay_group, decay_chain_struct, decay_data, data)[source]
cal_angle(data, decay_group: DecayGroup, using_topology=True, random_z=True, r_boost=True, final_rest=True, align_ref=None, only_left_angle=False)

Calculate helicity angle for particle momentum, add aligned angle.

Params data:

dict as {particle: {“p”:…}}

Returns:

Dictionary of data

cal_angle_from_momentum(p, decs: DecayGroup, using_topology=True, center_mass=False, r_boost=True, random_z=False, batch=65000, align_ref=None, only_left_angle=False) CalAngleData[source]

Transform 4-momentum data in files for the amplitude model automatically via DecayGroup.

Parameters:
  • p – 4-momentum data

  • decs – DecayGroup

Returns:

Dictionary of data

cal_angle_from_momentum_base(p, decs: DecayGroup, using_topology=True, center_mass=False, r_boost=True, random_z=False, batch=65000, align_ref=None, only_left_angle=False) CalAngleData[source]

Transform 4-momentum data in files for the amplitude model automatically via DecayGroup.

Parameters:
  • p – 4-momentum data

  • decs – DecayGroup

Returns:

Dictionary of data

cal_angle_from_momentum_id_swap(p, decs: DecayGroup, using_topology=True, center_mass=False, r_boost=True, random_z=False, batch=65000, align_ref=None, only_left_angle=False) CalAngleData[source]
cal_angle_from_momentum_single(p, decs: DecayGroup, using_topology=True, center_mass=False, r_boost=True, random_z=True, align_ref=None, only_left_angle=False) CalAngleData[source]

Transform 4-momentum data in files for the amplitude model automatically via DecayGroup.

Parameters:
  • p – 4-momentum data

  • decs – DecayGroup

Returns:

Dictionary of data

cal_angle_from_particle(data, decay_group: DecayGroup, using_topology=True, random_z=True, r_boost=True, final_rest=True, align_ref=None, only_left_angle=False)[source]

Calculate helicity angle for particle momentum, add aligned angle.

Params data:

dict as {particle: {“p”:…}}

Returns:

Dictionary of data

cal_chain_boost(data, decay_chain: DecayChain) dict[source]

calculate chain boost for a decay chain

cal_helicity_angle(data: dict, decay_chain: DecayChain, base_z=array([0., 0., 1.]), base_x=array([1., 0., 0.])) dict[source]

Calculate helicity angle for A -> B + C: \(\theta_{B}^{A}, \phi_{B}^{A}\) from momentum.

from {A:{p:momentum},B:{p:...},C:{p:...}}

to {A->B+C:{B:{"ang":{"alpha":...,"beta":...,"gamma":...},"x":...,"z"},...}}

cal_single_boost(data, decay_chain: DecayChain) dict[source]
cp_swap_p(p4, finals, id_particles, cp_particles)[source]
get_chain_data(data, decay_chain=None)[source]

get all independent data for a decay chain

get_key_content(dic, key_path)[source]

get key content. E.g. get_key_content(data, ‘/particle/(B, C)/m’)

>>> data = {"particle": {"(B, C)": {"p": 0.1, "m": 1.0},"B":1.0}}
>>> get_key_content(data, '/particle/(B, C)/m')
1.0
get_keys(dic, key_path='')[source]

get_keys of nested dictionary

>>> a = {"a": 1, "b": {"c": 2}}
>>> get_keys(a)
['/a', '/b/c']
get_relative_momentum(data: dict, decay_chain: DecayChain)[source]

add add rest frame momentum scalar from data momentum.

from {"particle": {A: {"m": ...}, ...}, "decay": {A->B+C: {...}, ...}

to {"particle": {A: {"m": ...}, ...},"decay": {A->B+C:{...,"|q|": ...},...}

identical_particles_swap(id_particles)[source]
identical_particles_swap_p(p4, id_particles)[source]
infer_momentum(data, decay_chain: DecayChain) dict[source]
infer momentum of all particles in the decay chain from outer particles momentum.

{outs:{p:momentum}} => {top:{p:momentum},inner:{p:..},outs:{p:..}}

parity_trans(p, charges)[source]
prepare_data_from_dat_file(fnames)[source]

[deprecated] angle for amplitude.py

prepare_data_from_dat_file4(fnames)[source]

[deprecated] angle for amplitude4.py

prepare_data_from_decay(fnames, decs, particles=None, dtype=None, charges=None, **kwargs)[source]

Transform 4-momentum data in files for the amplitude model automatically via DecayGroup.

Parameters:
  • fnames – File name(s).

  • decs – DecayGroup

  • particles – List of Particle. The final particles.

  • dtype – Data type.

Returns:

Dictionary

struct_momentum(p, center_mass=True) dict[source]
restructure momentum as dict

{outs:momentum} => {outs:{p:momentum}}