class TwentyOne(LymanAlpha): def __init__(self, cosmo=None, z_reio=7.82, z_min=10.0, z_max=150.0, f_star_L=0.5, f_star_X=0.004, T_vir_cut=1e4 * Kelv, use_hirata_fits=True, sed_X="PL", sed_X_kwargs={}, hmf_kwargs={ "mdef": "vir", "model": "despali16" }): """ Class to calculate Lyman-alpha heating :param cosmo: The cosmology, specified as an astropy.cosmology.FlatLambdaCDM instance :param z_reio: Redshift at reionization, passed to CLASS :param z_min: Minimum redshift (for interpolation tables) :param z_max: Maximum redshift (for interpolation tables) :param f_star_L: Efficiency of star formation, passed to Ly-A/X-ray classes :param f_star_X: Efficiency of star formation, X-ray, passed to Ly-A/X-ray classes :param T_vir_cut: Minimum virial temperature of star-forming halos, passed to Ly-A/X-ray classes :param use_hirata_fits: Whether to use fitting functions from Hirata for S_alpha and T_c :param sed_X: X-ray luminosity, "PL" or "Hirata" :param sed_X_kwargs: Parameters for X-ray luminosity """ LymanAlpha.__init__(self, cosmo=cosmo, z_min=z_min, z_max=z_max, f_star_L=f_star_L, f_star_X=f_star_X, T_vir_cut=T_vir_cut, sed_X=sed_X, sed_X_kwargs=sed_X_kwargs, hmf_kwargs=hmf_kwargs) self.z_reio = z_reio # Redshift at reionization self.use_hirata_fits = use_hirata_fits # Whether to use fitting functions from Hirata for S_alpha and T_c self.data_path = str(Path(__file__).parent / "../data/") self.load_constants() # Set class-specific constants self.load_interpolations() # Load interpolation tables self.initialize_class_inst() # Initialize CLASS instance def load_constants(self): """ Load class-specific constants """ self.nu_Lya = 2466 * 1e12 / Sec # Ly-A absorption frequency, originally in THz self.gamma = 50 * 1e6 / Sec # HWHM of the 21-cm resonance, after Eq. (11) of astro-ph/0507102, originally in MHz self.T_21 = 68.2e-3 # 21-cm temperature, in K self.nu_21 = 1420405751.7667 / Sec # Frequency of 21-cm transition self.A_21 = 2.86e-15 / Sec # Spontaneous emission (Einstein-A) coefficient of 21-cm transition, after Eq. (3) of 0802.2102 self.lambda_Lya = 121.567e-9 * Meter # Ly-A absorption wavelength self.sigma_T = 0.66524587158 * barn # Thomson scattering cross-section self.EHth = 13.6 * eV # Hydrogen ground state energy def initialize_class_inst(self): """ Get electron ionization fraction from CLASS """ class_parameters = { "H0": self.cosmo.H0.value, "Omega_b": self.cosmo.Ob0, "N_ur": self.cosmo.Neff, "Omega_cdm": self.cosmo.Odm0, "YHe": self.Y_p, "z_reio": self.z_reio } self.CLASS_inst = Class() self.CLASS_inst.set(class_parameters) self.CLASS_inst.compute() def x_e(self, z): """ Electron ionization fraction, from CLASS instance """ return self.CLASS_inst.ionization_fraction(z) def T_b(self, z): """ Baryon temperature at given redshift, from CLASS instance """ return self.CLASS_inst.baryon_temperature(z) def load_interpolations(self): """ Load interpolation tables """ ## Load table from https://github.com/ntveem/lyaheating heffs = np.load(self.data_path + "/heffs.npy") # heffs[:, :, :, 1, 0][np.where(heffs[:, :, :, 1, 0] < 0)] = 1e-15 # Argument arrays for T_k, T_s,... l10_t_ary = np.linspace(np.log10(0.1), np.log10(100.0), num=175) # ... and tau_GP (Gunn-Peterson optical depth) l10_tau_gp_ary = np.linspace(4.0, 7.0) # Net energy loss efficiency, defined in Eq. (37) of 1804.02406 self.E_c_interp = RegularGridInterpolator( points=[l10_t_ary, l10_t_ary, l10_tau_gp_ary], values=((heffs[:, :, :, 0, 0])), bounds_error=False, fill_value=None) self.E_i_interp = RegularGridInterpolator( points=[l10_t_ary, l10_t_ary, l10_tau_gp_ary], values=((heffs[:, :, :, 1, 0])), bounds_error=False, fill_value=None) # Energy loss to spins, defined in Eq. (32) of astro-ph/0507102 self.S_alpha_c_interp = RegularGridInterpolator( points=[l10_t_ary, l10_t_ary, l10_tau_gp_ary], values=(np.log10(heffs[:, :, :, 0, 2])), bounds_error=False, fill_value=None) self.S_alpha_i_interp = RegularGridInterpolator( points=[l10_t_ary, l10_t_ary, l10_tau_gp_ary], values=(np.log10(heffs[:, :, :, 1, 2])), bounds_error=False, fill_value=None) # Effective colour temperature, defined in Eq. (32) of astro-ph/0507102 self.T_c_c_interp = RegularGridInterpolator( points=[l10_t_ary, l10_t_ary, l10_tau_gp_ary], values=(np.log10(heffs[:, :, :, 0, 3])), bounds_error=False, fill_value=None) self.T_c_i_interp = RegularGridInterpolator( points=[l10_t_ary, l10_t_ary, l10_tau_gp_ary], values=(np.log10(heffs[:, :, :, 1, 3])), bounds_error=False, fill_value=None) ## Rate coefficients # From astro-ph/0608067, Table 1 kappa_10_eH_ary = np.loadtxt(self.data_path + "/kappa_10_eH_tab.dat") # From Zygelman (2005), http://adsabs.harvard.edu/abs/2005ApJ...622.1356Z, Table 2 kappa_10_HH_ary = np.loadtxt(self.data_path + "/kappa_10_HH_tab.dat") self.l10_kappa_10_eH_interp = interp1d( np.log10(kappa_10_eH_ary)[:, 0], np.log10(kappa_10_eH_ary * Centimeter**3 / Sec)[:, 1], bounds_error=False, fill_value="extrapolate") self.l10_kappa_10_HH_interp = interp1d( np.log10(kappa_10_HH_ary)[:, 0], np.log10(kappa_10_HH_ary * Centimeter**3 / Sec)[:, 1], bounds_error=False, fill_value="extrapolate") def S_alpha_Hirata(self, Tk, Ts, tauGP): """ Hirata fitting functional form for energy loss to spins """ xi = (1e-7 * tauGP)**(1.0 / 3.0) * Tk**(-2.0 / 3.0) a = 1.0 - 0.0631789 / Tk + 0.115995 / Tk**2 - 0.401403 / Ts / Tk + 0.336463 / Ts / Tk**2 b = 1.0 + 2.98394 * xi + 1.53583 * xi**2 + 3.85289 * xi**3 return a / b def T_c_Hirata(self, Tk, Ts): """ Hirata fitting functional form for effective colour temperature """ T_c_inv = Tk**(-1.0) + 0.405535 * Tk**(-1.0) * (Ts**(-1.0) - Tk** (-1.0)) return 1 / T_c_inv def T_c_c(self, T_k, T_s, x_e, z): """ Effective color temperature from interpolation, continuum """ if T_k <= 100.0: return 10**self.T_c_c_interp( [np.log10(T_k), np.log10(T_s), np.log10(self.tau_GP(x_e, z))]) else: return T_k / 100 * 10**self.T_c_c_interp([ np.log10(100.0), np.log10(T_s), np.log10(self.tau_GP(x_e, z)) ]) # self.T_b(z) def T_c_i(self, T_k, T_s, x_e, z): """ Effective color temperature from interpolation, injected """ if T_k <= 100.0: return 10**self.T_c_i_interp( [np.log10(T_k), np.log10(T_s), np.log10(self.tau_GP(x_e, z))]) else: return T_k / 100 * 10**self.T_c_i_interp([ np.log10(100.0), np.log10(T_s), np.log10(self.tau_GP(x_e, z)) ]) # self.T_b(z) def x_HI(self, x_e): """ IGM neutral fraction, from electron ionization fraction """ # return np.max(1 - x_e, 0) return np.where(x_e <= 1 + self.Y_p / 4 * (1 - self.Y_p), 1 - x_e * (1 - self.Y_p / (4 - 3 * self.Y_p)), 0) def T_CMB(self, z): """ CMB temperature, in K """ return self.T_CMB_0 * (1 + z) def tau_GP(self, x_e, z): """ Gunn-Peterson optical depth, Eq. (35) of astro-ph/0507102 """ H = self.cosmo.H(z).value * Kmps / Mpc return (3 * self.n_H(z) * self.x_HI(x_e) * self.gamma) / (2 * H * self.nu_Lya**3) def tau_21(self, T_s, x_e, z): """ 21-cm optical depth, Eq. (2) of 1804.02406 """ H = self.cosmo.H(z).value * Kmps / Mpc return 3 / (32 * np.pi) * (self.n_H(z) * self.x_HI(x_e) * self.A_21 ) / (self.nu_21**3 * H) * self.T_21 / T_s def x_CMB(self, T_s, x_e, z): """ Spin-flip rate due to CMB heating, Eq. (14) of 1804.02406 """ return 1 / self.tau_21(T_s, x_e, z) * (1 - np.exp(-self.tau_21(T_s, x_e, z))) def x_c(self, T_k, T_gamma, x_e, z): """ Collisional coupling spin-flip rate coefficient, Eq (3) of 0802.2102 """ return 4 * self.T_21 / (3 * self.A_21 * T_gamma) * self.n_H(z) * ( 10**self.l10_kappa_10_eH_interp(np.log10(T_k)) * x_e + 10**self.l10_kappa_10_HH_interp(np.log10(T_k))) def x_alpha_c(self, T_k, T_s, T_gamma, x_e, J_c_o_J_0, z): """ Wouthuysen-Field spin-flip rate coefficient, continuum, Eq. (11) of astro-ph/0507102 """ if self.use_hirata_fits: S_alpha_c = self.S_alpha_Hirata(T_k, T_s, self.tau_GP(x_e, z)) else: S_alpha_c = 10**self.S_alpha_c_interp( np.log10([T_k, T_s, self.tau_GP(x_e, z)])) return 8 * np.pi * self.lambda_Lya**2 * self.gamma * self.T_21 / ( 9 * self.A_21 * T_gamma) * S_alpha_c * J_c_o_J_0 * self.J_0(z) def x_alpha_i(self, T_k, T_s, T_gamma, x_e, J_i_o_J_0, z): """ Wouthuysen-Field spin-flip rate coefficient, injected, Eq. (11) of astro-ph/0507102 """ if self.use_hirata_fits: S_alpha_i = self.S_alpha_Hirata(T_k, T_s, self.tau_GP(x_e, z)) else: S_alpha_i = 10**self.S_alpha_i_interp( np.log10([T_k, T_s, self.tau_GP(x_e, z)])) return 8 * np.pi * self.lambda_Lya**2 * self.gamma * self.T_21 / ( 9 * self.A_21 * T_gamma) * S_alpha_i * J_i_o_J_0 * self.J_0(z) def T_s_inv(self, T_s, T_k, T_gamma, x_e, J, z): """ Spin temperature (inverse), e.g. Eq (13) of 1804.02406 """ x_CMB = self.x_CMB(T_s, x_e, z) x_alpha_c = self.x_alpha_c(T_k, T_s, T_gamma, x_e, J[0], z) x_alpha_i = self.x_alpha_i(T_k, T_s, T_gamma, x_e, J[1], z) x_c = self.x_c(T_k, T_gamma, x_e, z) if self.use_hirata_fits: T_c_c = T_c_i = self.T_c_Hirata(T_k, T_s) else: T_c_c = self.T_c_c(T_k, T_s, x_e, z) T_c_i = self.T_c_i(T_k, T_s, x_e, z) return (x_CMB * T_gamma**-1 + x_alpha_c * T_c_c**-1 + x_alpha_i * T_c_i**-1 + x_c * T_k**-1) / (x_CMB + x_alpha_c + x_alpha_i + x_c) def T_s_solve(self, T_k, T_gamma, x_e, J, z): """ Solve for spin temperature """ T_s = (root( lambda T_s: self.T_s_inv(T_s[0], T_k, T_gamma, x_e, J, z) - T_s** -1, np.min([T_k, T_gamma])).x)[0] return T_s def E_CMB(self, T_s, T_gamma, T_k, x_e, z): """ Heating efficiency due to CMB, from Eq. (17) of 1804.02406 """ H = self.cosmo.H(z).value * Kmps / Mpc return self.x_HI(x_e) * self.A_21 / (2 * H) * self.x_CMB( T_s, x_e, z) * (T_gamma / T_s - 1) * self.T_21 / T_k def E_Compton(self, T_k, x_e, z): """ Compton heating efficiency, from Eq. (22) of 1312.4948 (TODO: but is it) """ H = self.cosmo.H(z).value * Kmps / Mpc a_r = np.pi**2 / 15.0 return 8 * self.sigma_T * a_r * (self.T_CMB(z) * Kelv)**4 * x_e / ( 3 * m_e * H) * (self.T_CMB(z) / T_k - 1) def dT_k_dz(self, T_s, T_k, T_gamma, x_e, J, z): """ Kinetic temperature evolution, from Eq (18) of 1804.02406 """ E_c = self.E_c_interp(np.log10([T_k, T_s, self.tau_GP(x_e, z)])) E_i = self.E_i_interp(np.log10([T_k, T_s, self.tau_GP(x_e, z)])) dT_k_dz = 1 / (1 + z) * ( 2 * T_k - 1 / (1 + self.f_He + x_e) * (E_c * J[0] + E_i * J[1] + self.E_CMB(T_s, T_gamma, T_k, x_e, z) + self.E_Compton(T_k, x_e, z)) * T_k) return dT_k_dz - 1 / (1 + z) * self.heat_coeff(z) def alpha_A(self, T): """ Case-A recombination coefficient, from Pequignot et al (1991), Eq. (1) and Table 1 """ zi = 1.0 a = 5.596 b = -0.6038 c = 0.3436 d = 0.4479 t = 1e-4 * T / zi**2 return zi * (a * t**b) / (1 + c * t**d) * 1e-13 * Centimeter**3 / Sec def alpha_B(self, T): """ Case-B recombination coefficient, from `DarkHistory` """ return alpha_recomb( (k_B * T * Kelv) / eV, species="HI") * Centimeter**3 / Sec def rec_rate(self, z, x_e, T_k, case="B"): """ Recombination rate (Eq. 29 of 1312.4948) """ Gamma_rec = -peebles_C(1 - self.x_HI(x_e), 1 + z) * self.alpha_B(T_k) * x_e**2 * self.n_H(z) return self.dz_dt(z)**-1 * Gamma_rec
class TransitionProbabilities: def __init__(self, cosmo=None, z_reio=7.82, z_recomb=1089.80): """ Container class to compute dark photon transition probabilities. :param cosmo: Cosmology as an astropy object. If `None`, defaults to Planck18 cosmology. :param z_reio: Redshift of reionization. By default consistent with Planck18 cosmology. :param z_recomb: Redshift of recombination. By default consistent with Planck18 cosmology. """ # Set constants self.xi_3 = 1.20205 # Apéry's constant, from Wikipedia self.eta = 6.129e-10 # Baryon-to-photon ratio, from 1912.01132 self.Y_p = 0.247 # Helium-to-hydrogen mass fraction, from 1912.01132 self.T_0 = 2.725 # CMB temperature today, from 0911.1955 self.T_0_n = (c.k_B * self.T_0 * u.Kelvin).to(u.eV).value * \ eV # CMB temperature today, in natural units # Set cosmology if cosmo is None: self.cosmo = self.get_Planck18_cosmology() else: self.cosmo = cosmo # Set redshift of reionization self.z_reio = z_reio self.z_recomb = z_recomb # Initialize CLASS instance to compute ionization fraction from self.initialize_class_inst() def get_Planck18_cosmology(self): """ Stolen from https://github.com/abatten/fruitbat/blob/c074abff432c3b267d00fbb49781a0e0c6eeab75/fruitbat/cosmologies.py Planck 2018 paper VI Table 2 Final column (68% confidence interval) This is the Planck 2018 cosmology that will be added to Astropy when the paper is accepted. :return: astropy.cosmology.FlatLambdaCDM instance describing Planck18 cosmology """ planck18_cosmology = { 'Oc0': 0.2607, 'Ob0': 0.04897, 'Om0': 0.3111, 'H0': 67.66, 'n': 0.9665, 'sigma8': 0.8102, 'tau': 0.0561, 'z_reion': 7.82, 't0': 13.787, 'Tcmb0': 2.7255, 'Neff': 3.046, 'm_nu': [0., 0., 0.06], 'z_recomb': 1089.80, 'reference': "Planck 2018 results. VI. Cosmological Parameters, " "A&A, submitted, Table 2 (TT, TE, EE + lowE + lensing + BAO)" } Planck18 = FlatLambdaCDM(H0=planck18_cosmology['H0'], Om0=planck18_cosmology['Om0'], Tcmb0=planck18_cosmology['Tcmb0'], Neff=planck18_cosmology['Neff'], Ob0=planck18_cosmology['Ob0'], name="Planck18", m_nu=u.Quantity(planck18_cosmology['m_nu'], u.eV)) return Planck18 def initialize_class_inst(self): """ Get electron ionization fraction from CLASS """ class_parameters = { 'H0': self.cosmo.H0.value, 'Omega_b': self.cosmo.Ob0, 'N_ur': self.cosmo.Neff, 'Omega_cdm': self.cosmo.Odm0, 'YHe': self.Y_p, 'z_reio': self.z_reio } self.CLASS_inst = Class() self.CLASS_inst.set(class_parameters) self.CLASS_inst.compute() # z_ary = np.logspace(-3, 5, 100000) z_ary = np.linspace(0, 33000., 300000) x_e_ary = [self.CLASS_inst.ionization_fraction(z) for z in z_ary] self.x_e = interp1d(z_ary, x_e_ary) def m_A_sq(self, z, omega, x_e=None): """ Effective photon plasma mass squared, in natural units, from 1507.02614 :param z: Redshift :param omega: Photon frequency :param x_e: Free electron fraction if not default (optional) :return: Effective photon plasma mass squared, in natural units """ if x_e is None: x_e = self.x_e(z) m_A_sq = 1.4e-21 * (x_e - 6e-3 * (omega / eV)**2 * (1 - x_e)) * (self.n_H(z) / Centimeter**-3) return m_A_sq * eV**2 # Convert to natural units def n_p(self, z): """ Proton density at redshift `z` in cm^3, from 1507.02614 """ return (1 - self.Y_p / 2.) * self.eta * 2 * self.xi_3 / np.pi**2 * ( self.T_0_n * (1 + z))**3 def n_H(self, z): """ Number density of hydrogen nuclei at redshift `z` in cm^3 """ return (1 - self.Y_p) * self.eta * 2 * self.xi_3 / np.pi**2 * (self.T_0_n * (1 + z))**3 def n_b(self, z): """ Baryon density at redshift `z` in cm^3. """ return self.eta * 2 * self.xi_3 / np.pi**2 * (self.T_0_n * (1 + z))**3 def dz_dt(self, z): """ dz/dt """ return -self.cosmo.H(z).value * Kmps / Mpc * (1 + z) def omega(self, omega_0, z, evolve_z=True): """ Frequency corresponding to present-day omega_0 evolved to `z` is `evolve_z` is `True`, otherwise just return `omega_0`. """ if evolve_z: return omega_0 * (1 + z) else: return omega_0 def get_z_crossings(self, m_A, omega_0, evolve_z=True): """ Find redshifts at which resonance occurs :param m_A: Dark photon mass :param omega_0: Present-day frequency :param evolve_z: Whether to evolve frequency in redshift. :return: Array of redshifts at which resonance occurs """ z_ary = np.logspace(-3, 4.5, 20000) m_A_ary = np.nan_to_num(np.sqrt( self.m_A_sq(z_ary, self.omega(omega_0, z_ary, evolve_z))), nan=1e-18 * eV) where_ary = np.where( np.logical_or((m_A_ary[:-1] < m_A) * (m_A_ary[1:] > m_A), (m_A_ary[:-1] > m_A) * (m_A_ary[1:] < m_A))) def m_A_sq(z): return np.nan_to_num(np.sqrt( self.m_A_sq(z, self.omega(omega_0, z, evolve_z))), nan=1e-18 * eV) - m_A z_cross_ary = [] for i in range(len(where_ary[0])): z_cross_ary.append( brentq(m_A_sq, z_ary[where_ary[0][i]], z_ary[where_ary[0][i] + 1])) return np.array(z_cross_ary) def P_trans(self, m_A, z_res_ary, omega_0, eps, evolve_z=True, approx_linearize=True): """ Photon transition probability :param m_A: Dark photon mass :param z_res_ary: Array of resonance redshifts :param omega_0: Photon frequency (present day if `approx_linearize`, otherwise absolute) :param eps: Kinetic mixing coupling :param evolve_z: Whether to evolve `omega_0` in redshift :param approx_linearize: Linearize probability in `epsilon` :return: Transition probability array at redshifts `z_res_ary` """ d_log_m_A_sq_dz = np.array([ derivative(lambda z: np.log( self.m_A_sq(z=z, omega=self.omega(omega_0, z, evolve_z))), x0=z, dx=1e-7) for z in z_res_ary ]) omega_res_ary = self.omega(omega_0, z_res_ary, evolve_z) if approx_linearize: P_homo = np.pi * m_A ** 2 * eps ** 2 / omega_res_ary * \ np.abs((d_log_m_A_sq_dz * self.dz_dt(z_res_ary))) ** -1 else: r = np.abs((d_log_m_A_sq_dz * self.dz_dt(z_res_ary)))**-1 k = m_A**2 / (2 * omega_res_ary) P_homo = 1 - np.exp(-2 * np.pi * r * k * np.sin(eps)**2) return np.nan_to_num(P_homo) def P_tot(self, omega_0, eps, m_A, approx_linearize=True, evolve_z=True, sum_probs=False, **kwargs): """ Total conversion probability in the homogeneous limit :param omega_0: Present-day photon frequency :param eps: Dark photon coupling :param m_A: Dark photon mass :param approx_linearize: Whether to use linearized probability approximation :param evolve_z: Whether to evolve frequency in redshift. :param sum_probs: Whether to sum over probabilities associated with different z :return: Redshift resonance array, transition probability array """ # Find redshift at which resonance occurs z_res_ary = [ self.get_z_crossings(m_A, omega, evolve_z) for omega in omega_0 ] # Get transition probabilities at resonance if sum_probs: P_ary = np.array([ np.nansum( self.P_trans(m_A, z, omega, eps, approx_linearize=approx_linearize, evolve_z=evolve_z)) for z, omega in zip(z_res_ary, omega_0) ]) else: P_ary = np.array([(self.P_trans(m_A, z, omega, eps, approx_linearize=approx_linearize, evolve_z=evolve_z)) for z, omega in zip(z_res_ary, omega_0)]) return z_res_ary, 1., P_ary def B_CMB(self, omega, T): """ CMB spectral intensity at frequency `omega` (in natural units) for temperature `T` (in Kelvin) """ T_N = (c.k_B * T * u.Kelvin).to(u.eV).value * eV return omega**3 / (2 * np.pi**2) * (np.exp(omega / T_N) - 1)**-1