rabacus.atomic package

Submodules

rabacus.atomic.badnell_06 module

A module that stores rate fits from two Badnell 06 papers http://adsabs.harvard.edu/abs/2006ApJS..167..334B http://adsabs.harvard.edu/abs/2006A%26A...447..389B

class rabacus.atomic.badnell_06.Badnell06[source]

Fits to radiative recombination rates from http://adsabs.harvard.edu/abs/2006ApJS..167..334B

Fits from dielectronic recombination rates from http://adsabs.harvard.edu/abs/2006A%26A...447..389B

Attributes:

reH2a(Tin)[source]

Fit to H2 (Z=1,N=0) recombination rate (caseA)

reHe2a(Tin)[source]

Fit to He2 (Z=2, N=1) recombination rate (caseA)

reHe2di(Tin)[source]

Fit to He2 dielectronic recombination rate

reHe3a(Tin)[source]

Fit to He3 (Z=2, N=0) recombination rate (caseA)

rabacus.atomic.chemistry module

A module for loading atomic chemistry rates.

class rabacus.atomic.chemistry.ChemistryRates(T, fcA_H2, fcA_He2, fcA_He3, H1i=None, He1i=None, He2i=None, fit_name='hg97', add_He2di=True, use_badnell=False)[source]

A chemistry rates class.

Provides access to collisional ionization and recombination rates as a function of temperature. The user can optionally supply photoionization rates which will become attributes of the returned instance. These chemistry and photoionization rates determine the evolution of the ionization fractions for primordial cosmological chemistry.

Args:

T (array): Temperatures.

fcA_H2 (array): Case A fraction for HII.

fcA_He2 (array): Case A fraction for HeII.

fcA_He3 (array): Case A fraction for HeIII.

Note

The case A fractions, fcA_H2, fcA_He2, and fcA_He3, interpolate between case A and case B recombination rates (0.0 = case B, 1.0 = case A).

Kwargs:

H1i (array): HI photoionization rate.

He1i (array): HeI photoionization rate.

He2i (array): HeII photoionization rate.

fit_name (string): Source of rate fits {hg97}

add_He2di (bool): If true, reHe2di is added to reHe2.

Attributes:

u (Units): Units.

ciH1 (array): HI collisional ionization rate.

ciHe1 (array): HeI collisional ionization rate.

ciHe2 (array): HeII collisional ionization rate.

reH2 (array): HII radiative recombination rate.

reHe2 (array): HeII radiative recombination rate.

reHe3 (array): HeIII radiative recombination rate.

reHe2di (array): HeII dielectronic recombination rate.

extend_dims(N)[source]

Takes the one dimensional rate arrays (whose size has been determined by the size of the input temperature array) and extends them to 2-D arrays using copies of the original structures such that the new arrays have dimensions (N, self.T.size). This is useful when one wants to use array operations to calculate ionization fractions for all possible combinations of temperature and density.

reduce_dims()[source]

The inverse of extend_dims().

set(T, fcA_H2, fcA_He2, fcA_He3, H1i=None, He1i=None, He2i=None)[source]

Updates rates using new temperatures, case A fractions, and photoionization rates.

Args:

T (array): Temperatures.

fcA_H2 (array): Case A fraction for HII.

fcA_He2 (array): Case A fraction for HeII.

fcA_He3 (array): Case A fraction for HeIII.

Kwargs:

H1i (array): HI photoionization rate.

He1i (array): HeI photoionization rate.

He2i (array): HeII photoionization rate.

rabacus.atomic.cooling module

A module for loading atomic cooling rates.

class rabacus.atomic.cooling.CoolingRates(T, fcA_H2, fcA_He2, fcA_He3, H1h=None, He1h=None, He2h=None, fit_name='hg97', add_He2di=True)[source]

A cooling rates class.

Provides access to collisional ionization, collisional excitation, recombination, bremsstrahlung, and inverse compton cooling rates as a function of temperature. The user can optionally supply photo-heating rates which will become attributes of the returned instance. These cooling and photo-heating rates determine the evolution of temperature in cosmological gas.

Note

For compton and bremsstrahlung cooling, only the temperature dependent part of rates are stored. See the functions return_cooling() and return_heating() for the proper multiplicative factors.

Note

The case A fractions, fcA_H2, fcA_He2, and fcA_He3, interpolate between case A and case B recombination rates (0.0 = case B, 1.0 = case A).

Args:

T (array): Temperatures.

fcA_H2 (array): Case A fraction for HII.

fcA_He2 (array): Case A fraction for HeII.

fcA_He3 (array): Case A fraction for HeIII.

Kwargs:

H1h (array): HI photoheating rate.

He1h (array): HeI photoheating rate.

He2h (array): HeII photoheating rate.

fit_name (string): Source of rate fits {‘hg97’}

add_He2di (bool): If true, recHe2di is added to recHe2.

Attributes:

u (Units): Units.

cicH1 (array): HI collisional ionization cooling rate.

cicHe1 (array): HeI collisional ionization cooling rate.

cicHe2 (array): HeII collisional ionization cooling rate.

cecH1 (array): HI collisional excitation cooling rate.

cecHe1 (array): HeI collisional excitation cooling rate.

cecHe2 (array): HeII collisional excitation cooling rate.

recH2 (array): HII radiative recombination cooling rate.

recHe2 (array): HeII radiative recombination cooling rate.

recHe3 (array): HeIII radiative recombination cooling rate.

recHe2di (array): HeII dielectronic recombination cooling rate.

compton (array): Inverse Compton scattering cooling rate.

bremss (array): Bremsstrahlung radiation cooling rate.

extend_dims(N)[source]

Takes the one dimensional rate arrays (whose size has been determined by the size of the input temperature array) and extends them to 2-D arrays using copies of the original structures such that the new arrays have dimensions (N, self.T.size). This is useful when one wants to use array operations to calculate cooling/heating for all possible combinations of temperature and density.

reduce_dims()[source]

The inverse of extend_dims().

return_cooling(nH, nHe, x, z, Hz=None, sd93_logZ=None)[source]

Returns a cooling rate in erg/(s cm^3).

Args:

nH (array): Number density of hydrogen.

nHe (array): Number density of helium.

x (ionization fractions): An object which contains ionization fractions. For example, instances of the classes in one_zone.

z (float): Redshift.

Kwargs:

Hz (float): Hubble parameter at z if Hubble cooling is desired.

sd93_logZ (float): set to log metallicity (log Z) to use metal cooling in CIE from Sutherland and Dopita 1993

return_heating(nH, nHe, x)[source]

Returns a heating rate in erg/(s cm^3).

Args:

nH (array): Number density of hydrogen.

nHe (array): Number density of helium.

x (ionization fractions): An object which contains ionization fractions. For example, instances of the classes in one_zone

set(T, fcA_H2, fcA_He2, fcA_He3, H1h=None, He1h=None, He2h=None)[source]

Updates rates using new temperatures, case A fractions, and photo-heating rates.

Args:

T (array): Temperatures.

fcA_H2 (array): Case A fraction for HII.

fcA_He2 (array): Case A fraction for HeII.

fcA_He3 (array): Case A fraction for HeIII.

Kwargs:

H1h (array): HI photoheating rate.

He1h (array): HeI photoheating rate.

He2h (array): HeII photoheating rate.

rabacus.atomic.helium module

Constants and functions relevant for helium atoms.

class rabacus.atomic.helium.Helium(px_fit_type='verner')[source]

Helium atom class.

Attributes:

rabacus.atomic.hui_gnedin_97 module

A module that stores rate fits from Hui & Gnedin 97 http://adsabs.harvard.edu/abs/1997MNRAS.292...27H

class rabacus.atomic.hui_gnedin_97.HuiGnedin97[source]

Fits to atomic data from http://adsabs.harvard.edu/abs/1997MNRAS.292...27H

Ferland3_1e9 are fits to the data from Ferland 1992 and are accurate to 2% from 3K - 1e9K http://adsabs.harvard.edu/abs/1992ApJ...387...95F

BurgessSeaton5e3_5e5 are fits to the data from Burgess and Seaton 1960 and are accurate to 10% from 5e3K - 5e5K http://adsabs.harvard.edu/abs/1960MNRAS.121..471B

Lotz1e4_1e9 are fits to the data from Lotz 1967 and are accurate to 3% from 1e4K to 1e9K http://adsabs.harvard.edu/abs/1967ApJS...14..207L

AP3e4_1e6 are fits to the data from Aldrovandi and Pequignot 1973 and are accurate to 5% from 3e4K to 1e6K http://adsabs.harvard.edu/abs/1973A%26A....25..137A

Black5e3_5e5 are from Black 1981 with the correction from Cen 1992 and are accurate to 10% from 5e3 to 5e5 http://adsabs.harvard.edu/abs/1981MNRAS.197..553B

Attributes:

u (Units)

pc (PhysicalConstants)

T_H1 (float): HI ionization threshold in temperature units.

T_He1 (float): HeI ionization threshold in temperature units.

T_He2 (float): HeII ionization threshold in temperature units.

bremss(Tin)[source]

Fit to temperature dependent part of bremsstrahlung cooling.

cecH1(Tin)[source]

Fit to H1 collisional excitation cooling rate (Lotz1e4_1e9)

cecHe2(Tin)[source]

Fit to He2 collisional excitation cooling rate (Lotz1e4_1e9)

ciH1(Tin)[source]

Fit to H1 collisional ionization rate (Lotz1e4_1e9).

ciHe1(Tin)[source]

Fit to He1 collisional ionization rate (Lotz1e4_1e9)

ciHe2(Tin)[source]

Fit to He2 collisional ionization rate (Lotz1e4_1e9)

cicH1(Tin)[source]

Fit to H1 collisional ionization cooling rate (Lotz1e4_1e9)

cicHe1(Tin)[source]

Fit to He1 collisional ionization cooling rate (Lotz1e4_1e9)

cicHe2(Tin)[source]

Fit to He2 collisional ionization cooling rate (Lotz1e4_1e9)

compton(Tin)[source]

Fit to temperature dependent part of compton cooling rate.

reH2a(Tin)[source]

Fit to H2 recombination rate (caseA) (Ferland3_1e9)

reH2b(Tin)[source]

Fit to H2 recombination rate (caseB) (Ferland3_1e9)

reHe2a(Tin)[source]

Fit to He2 recombination rate (caseA) (BurgessSeaton5e3_5e5)

reHe2b(Tin)[source]

Fit to He2 recombination rate (caseB) (BurgessSeaton5e3_5e5)

reHe2di(Tin)[source]

Fit to He2 dielectronic recombination rate (AP3e4_1e6)

reHe3a(Tin)[source]

Fit to He3 recombination rate (caseA) (Ferland3_1e9)

reHe3b(Tin)[source]

Fit to He3 recombination rate (caseB) (Ferland3_1e9)

recH2a(Tin)[source]

Fit to H2 recombination cooling rate (caseA) (Ferland3_1e9)

recH2b(Tin)[source]

Fit to H2 recombination cooling rate (caseB) (Ferland3_1e9)

recHe2a(Tin)[source]

Fit to He2 recombination cooling rate (caseA) (BurgessSeaton5e3_5e5)

recHe2b(Tin)[source]

Fit to He2 recombination cooling rate (caseB) (BurgessSeaton5e3_5e5)

recHe2di(Tin)[source]

Fit to He2 dielectronic recombination cooling rate (AP3e4_1e6)

recHe3a(Tin)[source]

Fit to He3 recombination cooling rate (caseA) (Ferland3_1e9)

recHe3b(Tin)[source]

Fit to He3 recombination cooling rate (caseB) (Ferland3_1e9)

rabacus.atomic.hydrogen module

Constants and functions relevant for hydrogen atoms.

class rabacus.atomic.hydrogen.Hydrogen(px_fit_type='verner')[source]

Hydrogen atom class.

Attributes:

analytic_slab_soltn_tau(nH, k, xH1)[source]

Analytic inverse solution for monochromatic radiation incident onto a constant density and temperature slab. Given a neutral fraction this function returns the optical depth into the slab.

Args:

nH (array): H number density

k (rates object): must contain the attributes ciH1, reH2, and H1i, see ChemistryRates

xH1 (array): neutral hydrogen fraction

Returns:

tauH (array): depth into the slab for given neutral fraction, tau = NH * sigma
analytic_soltn_H1i(nH, y, k, xH1)[source]

Analytic equilibrium solution for the photoionization rate given the neutral fraction xH1 = nH1/nH

Args:

nH (array): H number density

y (array): electrons from elements heavier than hydrogen, ne = nH * (xH2 + y)

k (rates object): must contain the attributes ciH1 and reH2, see ChemistryRates

xH1 (array): neutral hydrogen fraction

Returns:

H1i (array): hydrogen photoionization rate
analytic_soltn_xH1(nH, y, k)[source]

Analytic equilibrium solution for the neutral fraction xH1 = nH1/nH

Args:

nH (array): H number density

y (array): electrons from elements heavier than hydrogen, ne = nH * (xH2 + y)

k (rates object): must contain the attributes ciH1, reH2, and H1i, see ChemistryRates

Returns:

xH1 (array): neutral hydrogen fraction

rabacus.atomic.photo_xsection module

Generic access to all photoionization cross section routines. Ions are labeled with Z (atomic number) and N (electron number). For convenience, we provide named routines for H1, He1, and He2. All access to photoionization cross sections should be through this module.

class rabacus.atomic.photo_xsection.PhotoXsections(fit_type='verner')[source]

General photoionization cross sections class.

Kwargs:

fit_type (string): fits come from Verner {verner}

Attributes:

Eth_H1 (float): H1 ionization threshold energy

Eth_He1 (float): He1 ionization threshold energy

Eth_He2 (float): He2 ionization threshold energy

Eth(Z, N)[source]

Returns threshold ionization energy for ions defined by Z and N.

Args:

Z (int): atomic number (number of protons)

N (int): electron number (number of electrons)

Returns:

Eth (float): ionization energy
sigma(Z, N, E)[source]

Returns a photo-ionization cross section for an ion defined by Z and N at energies E.

Args:

Z (int): atomic number (number of protons)

N (int): electron number (number of electrons)

E (array): calculate cross section at these energies

Returns:

sigma (array): photoionization cross sections
sigma_H1(E)[source]

Convenience function that calls sigma() with Z and N set for neutral hydrogen.

Args:

E (array): calculate cross section at these energies

Returns:

sigma (array): photoionization cross sections
sigma_He1(E)[source]

Convenience function that calls sigma() with Z and N set for neutral helium.

Args:

E (array): calculate cross section at these energies

Returns:

sigma (array): photoionization cross sections
sigma_He2(E)[source]

Convenience function that calls sigma() with Z and N set for singly ionized helium.

Args:

E (array): calculate cross section at these energies

Returns:

sigma (array): photoionization cross sections
sigma_th(Z, N)[source]

Returns photo-ionization cross sections for ions defined by Z and N at their ionization thresholds.

Args:

Z (int): atomic number (number of protons)

N (int): electron number (number of electrons)

Returns:

sigma (float): photoionization cross sections

rabacus.atomic.sd93 module

A module that stores rate fits from Sutherland and Dopita 93 http://adsabs.harvard.edu//abs/1993ApJS...88..253S

class rabacus.atomic.sd93.SD93[source]

Fits to atomic data from http://adsabs.harvard.edu//abs/1993ApJS...88..253S Provides acess to Lambda_N from Table 6

Read and clean data.

Zcool(Tin, logZ)[source]

Interpolate metal cooling function at a fixed metallicity and many temperatures.

Args:

Tin (array like): temperatures

logZ (float): log metallicity

cool(Tin, logZ)[source]

Interpolate cooling function at a fixed metallicity and many temperatures.

Args:

Tin (array like): temperatures

logZ (float): log metallicity

Module contents

Classes to provide atomic data. In particular,

Note

Both chemistry and cooling objects take three case A fraction arguments, fcA_H2, fcA_He2, and fcA_He3. These variables interpolate between case A and case B recombination rates (0.0 = case B, 1.0 = case A).

Chemistry Examples:

The following example will return chemistry rates at 128 temperatures logarithmically spaced between 10^4 and 10^5 K using case A recombination rates.

import numpy as np
import rabacus as ra

NT = 128
T = 10**np.linspace( 4.0, 5.0, NT ) * ra.U.K
fcA_H2 = 1.0; fcA_He2 = 1.0; fcA_He3 = 1.0
kchem = ra.atomic.ChemistryRates( T, fcA_H2, fcA_He2, fcA_He3 )

The object kchem can now be used in calls to collisional equilibrium solvers such as Solve_CE. If you are interested in solutions which include non zero photoionization rates, such as Solve_PCE, then those rates must be passed to the chemistry object when it is created,

H1i = np.ones( NT ) * 1.0e-12 / ra.U.s
He1i = np.ones( NT ) * 1.0e-13 / ra.U.s
He2i = np.ones( NT ) * 1.0e-14 / ra.U.s
kchem = ra.atomic.ChemistryRates( T, fcA_H2, fcA_He2, fcA_He3,
                                  H1i=H1i, He1i=He1i, He2i=He2i )

Cooling Examples:

The cooling rates object behaves analogously to the chemistry rates object. The following example will return cooling rates at 128 temperatures logarithmically spaced between 10^4 and 10^5 K using recombination rates half way between case A and case B.

import numpy as np
import rabacus as ra

NT = 128
T = 10**np.linspace( 4.0, 5.0, NT ) * ra.U.K
fcA_H2 = 0.5; fcA_He2 = 0.5; fcA_He3 = 0.5
kcool = ra.atomic.CoolingRates( T, fcA_H2, fcA_He2, fcA_He3 )

If photoheating rates are needed in a particular solution, they should be passed to the cooling object when it is created,

H1h = np.ones( NT ) * 1.0e-12 * ra.U.eV / ra.U.s
He1h = np.ones( NT ) * 1.0e-12 * ra.U.eV / ra.U.s
He2h = np.ones( NT ) * 1.0e-14 * ra.U.eV / ra.U.s
kcool = ra.atomic.CoolingRates( T, fcA_H2, fcA_He2, fcA_He3,
                                H1h=H1h, He1h=He1h, He2h=He2h )

Note that the heating and cooling rates per volume for a given density of ions and redshift is calculated using the functions return_cooling() and return_heating().

Photoionization Cross Section Examples:

The photoionization cross section class contains functions to calculate the cross section as well as the threshold ionization energies of common ions. The following code will store the photoionization cross sections of HI, HeI and HeII at 200 frequencies between 1 and 100 Rydbergs:

import numpy as np
import rabacus as ra

Nnu = 200
px = ra.atomic.PhotoXsections()
q_min = 1.0
q_max = 1.0e2
log_E_min = np.log10( q_min * px.Eth_H1.magnitude )
log_E_max = np.log10( q_max * px.Eth_H1.magnitude )
E = 10**np.linspace( log_E_min, log_E_max, Nnu ) * px.Eth_H1.units
sigma_H1 = px.sigma_H1( E )
sigma_He1 = px.sigma_He1( E )
sigma_He2 = px.sigma_He2( E )