Introduction

phasr is a python package able to calculate scattering phase shifts for arbitrary radial potentials and calculate (elastic) electron nucleus scattering crosssections based on it. In combination with a fitting routine it can be used to fit elastic electron nucleus scattering data to extract charge density parameterizations. It is also possible to increase the precision sufficiently to quantitatively resolve parity violating electron scattering (optimisation for this use case is planed for a future update).

Versions/Change-Log:

Version 0.2.0. : Basic features implemented.
Version 0.1.4. : Build dirac_solvers module.
Version 0.1.3. : Build nuclei module.
Version 0.1.2. : Build skeleton.
Version 0.1.0. : First init. Empty.

About this Package

This package was developed with the goal of extracting charge density parameterizations based on (elastic) electron nucleus scattering crosssection data while properly accounting for Coulomb corrections using the so called phase shift model. While this procedure is known since the 60s, the documentation and availability of programs used back then, are very limited and hardly accessable. Furthermore, for a majority of nuclei the extractions of model-independent parameterizations of their charge densities have not been updated since the 80s and do not provide any uncertainty estimates. Since these also applied to the nuclei we were interested in, we build such a framework in python, which calculates elastic electron nucleus scattering crosssections based on input charge densities using the phase shift model, from which this package evolved. With such a more modern implementation and programming language, the calculations can become faster, broader accessable, easier to understand and easier to incorporate into existing programs. This package (in its pre-release version) has already been used successfully in the papers REFs to extract charge densities for nuclei of interest to \( \mu \to e \) conversion a (hypothetical) process of interest in particle physics. The package has been benchmarked for nuclei reaching from ... to ... .

Dependencies

This program is built with python 3. (add Version).

It requires the following packages to operate properly:
numpy (add Version)
scipy (add Version)
mpmath (add Version)
....... (TBD)

Getting Started

Install the package

To install and use the program please follow the following steps for setup. We assume you already have a working python environment setup on your machine (othewise see REF).
The package is listed on PyPI, which is why it can be easily installed via:

 pip install phasr 
If you want the latest (development) version, you can clone it from the phasr GitHub repository .
Note that you need need to also install all the required packages as listed above .

Use the package

You can import the python package the same way as every other package, i.e.:

 import phasr as phr 
You can then use the contained methods in the standard way. For example
 nucleus_Al27 = phr.nucleus(name='Al27',Z=13,A=27,rc=5)
 mmu = phr.masses.mmu
 boundstates_Al27 = phr.boundstates(nucleus_Al27,kappa=-1,lepton_mass=mmu) 
 print('E_1s=',boundstates_Al27.energy_levels[0],'MeV') 
will give you the groundstate energy of a muonic atom for the nucleus parameterization, described used a uniform charge distribution of charge 13 and a radius of 5 fm. We will assume the import with the alias phr throughout this manual.

Basic Modules

The code is structured into two main modules. The nuclei module is home of the nucleus class and manages all the relevant properties of a nucleus/isotope that is considered and calculates other properties based on the information it has available (e.g. the electric field from the charge density). The dirac_solvers module is home of the boundstates and continuumstates classes and manages solving the Dirac equation for a given nucleus based on its nucleus potential and returns the desired consequential quantities like for example electron wave functions, phase shifts or the crosssection. In all applicational scenarios one will first construct a nucleus object, which is then be used to construct a boundstates or continuumstates object. For calculating elastic scattering crosssections we supply the function crosssection_lepton_nucleus_scattering, which can be used to directly calculate the crosssection using the energy, the angle and a nucleus object.

The nucleus class

To get a nucleus parameterization, use the nucleus class to construct a nucleus object. A nucleus object requires a few minimum parameters (e.g. a name, a charge \(Z\) and a nucleon number \(A\)) and and has a lot of optional parameters (e.g. some parameters for a specific parameterization). Based on the parameters supplied the class decides what kind of nucleus it is supposed to construct. You can check what parameterization was choosen based on the attribute nucleus_type

 phr.nucleus(name,Z,A,m=None,**args) 
General parameters are:
  • name (str): Label for the nucleus, part of the filename if things about this nucleus are saved, should be unique.
  • Z (int): Nuclear charge number / proton number of the considered nucleus, necessary for the normalization. Will raise a warning if not consistent with the normalization of a provided charge density in **args.
  • A (int): Atomic mass number / nucleon number of the considered nucleus, necessary for recoil, enables lookup of nucleon mass, spin, parity and isotope abundancy.
  • m (float,optional) mass of the nucleus, necessary for recoil, is looked up based on \(Z\), \(A\) if not supplied.
  • abundance (float,optional) natural abundance of the nucleus, stored for convencience, is looked up based on \(Z\), \(A\) if not supplied.
  • spin (float,optional) spin of the nucleus, restricts what multipoles can exist, is looked up based on \(Z\), \(A\) if not supplied.
  • parity (float,optional) parity of the nucleus, stored for convenience, is looked up based on \(Z\), \(A\) if not supplied.
  • k_barrett (float,optional) value of k used to calculate barrett moments.
  • alpha_barrett (float,optional) value of alpha used to calculate barrett moments.
Parameterization specific parameters are:
nucleus_type   Parameters / **args
'coulomb' None
'fourier-bessel' a (1d-array): Array of parameters \(a_i\) in fm\(^{-3}\). \(N\) given by length of a.
R (float): Cutoff radius R in fm.
'oszillator-basis' Ci_dict (dict): dictionary of 1d-arrays \(C_i\), with multipole names as keys (see TODO).
'fermi' c (float): \(c\) parameter in fm for Fermi parameterization.
z (float): \(z\) parameter in fm for Fermi parameterization.
w (float,optional): \(w\) parameter for Fermi parameterization (default: w=0).
'gauss' b (float): \(b\) parameter in fm for Gauss parameterization.
'uniform' rc (float): \(r_c\) parameter in fm for Uniform parameterization.
'numerical' see below

The parameterizations are based on: \[ \rho_\text{FB}(r) = \begin{cases} \sum_{i=1}^N a_i ~ j_0(q_i r) & r \leq R \\ 0 & r > R \end{cases} \qquad q_i= \frac{i \pi}{R} \] \[ \rho_\text{Fermi}(r) = \rho_0 \frac{1+w \tfrac{r^2}{c^2}}{1 + e^{(r-c)/z}}\] \[ \rho_\text{Gauss}(r) = \frac{Z}{\pi^{3/2} b^3} e^{-(r/b)^2} \] \[ \rho_\text{Uni}(r) = \begin{cases} \frac{3 Z}{4 \pi r_c^3} & r \leq r_c \\ 0 & r > r_c \end{cases} \] \[ F_\text{Osz}(q) = e^{-u(q)/2} \sum_i C_i ~ u(q)^i; \qquad u(q)=\frac{q^2 b^2}{2} \] Alternatively one can also supply general descriptions for charge_density, form_factor or the different multipoles. These get the label of nucleus_type 'numerical'. Here the other quantities can be calculated based on the given ones using numerical methods. Possible parameters are:
Parameters/**args
charge_density (callable,optional): Function describing the charge density.
form_factor (callable,optional): Function describing the form factor.
electric_field (callable,optional): Function describing the electric field.
electric_potential (callable,optional): Function describing the electric potential.
form_factor_dict (dict,optional): Dictionary containing form factor functions for specific multipoles (keys are 'F'+multipole).
density_dict (dict,optional): Dictionary containing density functions for specific multipoles (keys are 'rho'+multipole).

Every nucleus has the following attributes that can be accessed: TODO

The dirac_solvers module

The dirac_solvers module has two main duties, covered by two classes. The boundstates class calculates the boundstate solutions of a given nucleus by solving the Dirac equation in the presence of the nucleus potential, which includes finding the boundstate energy first. The continuumstates class calculates the continuumstate solution for a given energy of a given nucleus by solving the Dirac equation in the presence of the nucleus potential. Finally, there are some functions that make use of the boundstate and continuumstate solutions by calculating elastic electron nucleus scattering crosssections and overlap integrals for \(\mu \to e\) conversion.
The boundstate and continuumstate wavefunctions are both parameterized according to \[ \psi = \psi_\kappa^\mu(\vec{r})= \frac{1}{r} \Bigg(\begin{array}{c} g_\kappa(r) \phi_\kappa^\mu(\hat{r}) \\ i f_\kappa(r) \phi_{-\kappa}^\mu(\hat{r}) \end{array} \Bigg), \] separating off the angular degrees of freedom, which are the same for all radial symmetric potential. The \(\kappa\) label refers to the angular momentum quantum numbers and spin according to \[j = |\kappa| - \frac{1}{2}; \quad l = \begin{cases} \kappa, & \quad \kappa > 0 \\ - \kappa - 1 , & \quad \kappa < 0 \end{cases}\] The numerical calculations then provide descriptions for \(g_\kappa(r)\) and \(f_\kappa(r)\).

The boundstates class:
The boundstates class is used to calculate the boundstates for a given nucleus object, beginning from the ground state for the given partial wave. First the energy for the boundstate is found by iteratively scanning for flips in the asymptotic behaviour of the wave functions, which indicate a boundstate in between, until a set precision is reached. For this energy then the solution of the Dirac equation gives the dirac wavefunction for that state.

 phr.boundstates(nucleus,kappa,lepton_mass,**args) 
  • nucleus (nucleus object): Nucleus, whose electric potential is used for the calculation.
  • kappa (int): selects total spin and angular momentum quantum number.
  • lepton_mass (float): mass of the lepton bound to the nucleus in MeV. Common choices are electron mass phr.masses.me or muon mass phr.masses.mmu. needs to be non-zero.
  • bindingenergy_limit_lower (float,optional): lower limit for the binding energy to look for boundstates. If not provided uses the potential at zero as the lower bound.
  • bindingenergy_limit_upper (float,optional): upper limit for the binding energy to look for boundstates (default: 0).
  • **args (dict,optional): Additional keywords can be used to pass parameters to the underlying numerical routine (including parameters for scipy.integrate.solve_IVP). Default values are set by phr.default_boundstate_settings.
Apart from storing the input paramters the boundstates object has the following attributes, which you might want to access:
  • principal_quantum_numbers (list of int's): principal quantum numbers of the calculated boundstates.
  • energy_levels (list of floats): energy values of the calculated boundstates.
  • find_next_energy_level() (callable): called to find the next boundstate energy (for that /(\kappa/)). The groundstate energy is always calculated at initialisation.
  • solve_IVP_at_current_energy() (callable): called to calculate the boundstate wavefunctions at the current boundstate energy. The groundstate wavefunction is always calculated at initialisation.
  • find_next_solution(**args) (callable): same as calling find_next_energy_level and solve_IVP_at_current_energy consecutively. Can be passed **args to change the solver_setting.
  • wavefunction_g_<state>(r) (callable): returns the upper radial component of the wavefunction evaluated at r in fm for the state <state>. Exists only if find_next_solution() was called enough times. Example: Groundstate is wavefunction_g_1s12 (using \(\kappa=-1\))
  • wavefunction_f_<state>(r) (callable): returns the lower radial component of the wavefunction evaluated at r in fm for the state <state>. Exists only if find_next_solution() was called enough times. Example: Groundstate is wavefunction_f_1s12 (using \(\kappa=-1\))
The states are labeled according to principal quantum numbers, angular momentum quantum number and total spin following the scheme of atomic orbitals here writte according to <state>=1s12, 2s12, ..., 2p12, 2d12, ..., 2p32, 2d32, ... .

The continuumstate class:
 phr.continuumstates(nucleus,kappa,energy,lepton_mass=0,**args) 
  • nucleus (nucleus object): Nucleus, whose electric potential is used for the calculation.
  • kappa (int): selects total spin and angular momentum quantum number.
  • energy (float): energy of the (initial) lepton in MeV.
  • lepton_mass (float,optional): mass of the lepton bound to the nucleus in MeV (default: 0). Can be set to 0 if energy \(\gg\) lepton_mass, which gives significant runtime improvements.
  • **args (dict,optional): Additional keywords can be used to pass parameters to the underlying numerical routine (including parameters for scipy.integrate.solve_IVP). Default values are set by phr.default_continuumstate_settings.
Apart from storing the input paramters the continuumstates object has the following attributes, which you might want to access:
  • solve_IVP() (callable): called to calculate the continuumstate wavefunctions for the given energy. Needs to be called once explicitly, does not happen at initialisation.
  • wavefunction_g(r) (callable): returns the upper radial component of the wavefunction evaluated at r in fm for the given energy. exists only if solve_IVP() was called.
  • wavefunction_f(r) (callable): returns the upper radial component of the wavefunction evaluated at r in fm for the given energy. exists only if solve_IVP() was called.
  • extract_phase_shift() (callable): called to calculate the asymptotic phase shift for the given energy. Needs to be called once explicitly, does not happen at initialisation. Faster than calling solve_IVP().
  • phase_difference (float): difference between the phase shift of a Coulomb potential and the considered nucleus potential.
  • phase_shift (float): phase shift of the considered nucleus potential. Input for elastic scattering crosssection.
The crosssection_lepton_nucleus_scattering function:
Used to calculate the crosssection for elastic lepton nucleus scattering for specific intial lepton energys and scattering angles based on a given nucleus.
 phr.crosssection_lepton_nucleus_scattering(energy,theta,nucleus,lepton_mass=0,recoil=True,subtractions=3,**args) 
  • energy (float): energy of the (initial) lepton in MeV.
  • theta (float / 1d-array): elastic scattering angle in rad (\(\theta=0\) forward direction, \(\theta=2\pi\) backward direction).
  • nucleus (nucleus object): nucleus against which the lepton is scattered.
  • lepton_mass (float,optional): mass of the lepton bound to the nucleus in MeV (default: 0). Can be set to 0 if energy \(\gg\) lepton_mass, which gives significant runtime improvements. No convergent implementation for non-zero lepton_mass yet.
  • recoil (bool,optional): if recoil of the nucleus at leading (kinematic) order should be considered (default: True).
  • subtractions (int, optional): number of subtractions to make the partial wave sum convergent (default: 3). subtractions=0 will not converge. subtractions=3 seems to work most consistently.
  • **args (dict,optional): Additional keywords can be used to pass parameters to the underlying numerical routines. The time-optimal parameters for a given dataset can be accessed by calling optimise_crosssection_precision() (see below).

  • return: crosssection for elastic lepton nucleus scattering of the given nucleus at the given energy and angles in MeV\(^{-2}\)
Possible **args/ parameters for the underlying numerical routines are:
Parameters/**args
N_partial_waves (int,optional): angular momenutum of the partial wave at which the series is terminated. (default: 100)
phase_difference_limit (float,optional): limit at which (numerical) value for the phase_difference attribute of a continuumstates all following partial waves are considered indistinguishable to the coulomb solution. (default:0).
energy_norm (float,optional): scaling factor for the units while solving the inital value problem to increase speed (default: \( \hbar c\) fm\(^{-1}\))
method (str,optional): method keyword passed to scipy.integrate.solve_IVP (default: 'DOP853' / as set in default_continuumstate_settings)
atol (float,optional): atol keyword passed to scipy.integrate.solve_IVP (default: 1e-12 as set in default_continuumstate_settings)
rtol (float,optional): rtol keyword passed to scipy.integrate.solve_IVP (default: 1e-9 as set in default_continuumstate_settings)

If you want to call crosssection_lepton_nucleus_scattering many times for the same energies and angles theta for similar nuclei (for example for a fit), it makes sense to assess the necessary precision for the numerical calculation using optimise_crosssection_precision. This functions can run a few minutes as it is scanning a large **args / parameter space, but can afterwards improve the runtime significantly. Unfortunately the necessary precision is very dependent on how large the energy and angles are (larger energies and angles require significantly more precision and runtime), which is why this is the best solution for a reasonable runtime. Of course you can always also try to finetune the parameters yourself, but note that the convergence might break down, and that you should always check if the calculated crosssection is still reasonable.
 phr.optimise_crosssection_precision(energy,theta,nucleus,lepton_mass=0,recoil=True,subtractions=3,crosssection_precision=1e-3,jump_forward_dist=1) 
  • energy (float): energy of the (initial) lepton in MeV.
  • theta (float / 1d-array): elastic scattering angle in rad (\(\theta=0\) forward direction, \(\theta=2\pi\) backward direction).
  • nucleus (nucleus object): nucleus against which the lepton is scattered.
  • lepton_mass (float,optional): mass of the lepton bound to the nucleus in MeV (default: 0). Can be set to 0 if energy \(\gg\) lepton_mass, which gives significant runtime improvements. No convergent implementation for non-zero lepton_mass yet.

  • return: args with the quickest runtime for crosssection_lepton_nucleus_scattering while still consistent with crosssection calculations with a higher precision up to a relative precision of crosssection_precision.

Common use cases

Calculate the scattering crosssection based on an analytic charge density parametrisation

TBD

Calculate the scattering crosssection based on a numerical charge density parametrisation

TBD

Calculate the left-right asymmetry for parity violating electron scattering

TBD

More Info

... .This website was build based on a template from github user charlyllo.