pynucastro.screening package

Screening routines

Submodules

pynucastro.screening.screen module

Python implementations of screening routines.

class pynucastro.screening.screen.PlasmaState(*args, **kwargs)[source]

Bases: PlasmaState

Stores precomputed values that are reused for all screening correction factor calculations.

Variables:
  • temp – temperature in K

  • dens – density in g/cm^3

  • qlam0z – TODO: from screen5

  • taufac – TODO: from screen5

  • aa – TODO: from screen5

  • abar – average atomic mass

  • zbar – average ion charge

  • z2bar – average (ion charge)^2

  • n_e – electron number density

  • gamma_e_fac – temperature-independent part of Gamma_e

class_type = jitclass.PlasmaState#7fa611f388d0<temp:float64,dens:float64,qlam0z:float64,taufac:float64,aa:float64,abar:float64,zbar:float64,z2bar:float64,n_e:float64,gamma_e_fac:float64>
class pynucastro.screening.screen.ScreenFactors(*args, **kwargs)[source]

Bases: ScreenFactors

Stores values that will be used to calculate the screening correction factor for a specific pair of nuclei.

Variables:
  • z1 – atomic number of first nucleus

  • z2 – atomic number of second nucleus

  • a1 – atomic mass of first nucleus

  • a2 – atomic mass of second nucleus

  • zs13 – (z1+z2)**(1/3)

  • zhat – combination of z1 and z2 raised to the 5/3 power

  • zhat2 – combination of z1 and z2 raised to the 5/12 power

  • lzav – log of effective charge

  • aznut – combination of a1, z1, a2, z2 raised to 1/3 power

  • ztilde – effective ion radius factor for a MCP

class_type = jitclass.ScreenFactors#7fa5fe494fd0<z1:int64,z2:int64,a1:int64,a2:int64,zs13:float64,zhat:float64,zhat2:float64,lzav:float64,aznut:float64,ztilde:float64>
pynucastro.screening.screen.chugunov_2007(state, scn_fac)[source]

Calculates screening factors based on Chugunov et al. [2007].

Follows the approach in Yakovlev et al. [2006] to extend to a multi-component plasma.

Parameters:
  • state (PlasmaState) – the precomputed plasma state factors

  • scn_fac (ScreenFactors) – the precomputed ion pair factors

Returns:

screening correction factor

pynucastro.screening.screen.chugunov_2009(state, scn_fac)[source]

Calculates screening factors based on Chugunov and DeWitt [2009].

Parameters:
  • state (PlasmaState) – the precomputed plasma state factors

  • scn_fac (ScreenFactors) – the precomputed ion pair factors

Returns:

screening correction factor

pynucastro.screening.screen.make_plasma_state(temp, dens, molar_fractions)[source]

Construct a PlasmaState object from simulation data.

Parameters:
  • temp – temperature in K

  • dens – density in g/cm^3

  • molar_fractions – dictionary of molar abundances for each ion, as returned by Composition.get_molar()

pynucastro.screening.screen.make_screen_factors(n1, n2)[source]

Construct a ScreenFactors object from a pair of nuclei.

Parameters:
pynucastro.screening.screen.potekhin_1998(state, scn_fac)[source]

Calculates screening factors based on Chabrier and Potekhin [1998].

Parameters:
  • state (PlasmaState) – the precomputed plasma state factors

  • scn_fac (ScreenFactors) – the precomputed ion pair factors

Returns:

screening correction factor

pynucastro.screening.screen.screen5(state: PlasmaState, scn_fac)[source]

Calculates screening factors following the appendix of Wallace et al. [1982].

Based on Graboske et al. [1973] for weak screening. Based on Alastuey and Jancovici [1978] with plasma parameters from Itoh et al. [1979], for strong screening.

pynucastro.screening.screening_util module

Some helper functions for determining which rates need screening

class pynucastro.screening.screening_util.ScreeningPair(name, nuc1, nuc2, rate=None)[source]

Bases: object

a pair of nuclei that will have rate screening applied. We store a list of all rates that match this pair of nuclei

add_rate(rate)[source]
pynucastro.screening.screening_util.get_screening_map(rates, *, symmetric_screening=False)[source]

a screening map is just a list of ScreeningPar objects containing the information about nuclei pairs for screening If symmetric_screening=True, then for reverse rates, we screen using the forward rate nuclei (assuming that we got here via detailed balance).