Added help strings and ksz support
This commit is contained in:
parent
ff01447ab6
commit
247d8512e4
7 changed files with 333 additions and 3 deletions
3
python_sample/ksz/__init__.py
Normal file
3
python_sample/ksz/__init__.py
Normal file
|
@ -0,0 +1,3 @@
|
|||
from .constants import *
|
||||
from .gal_prof import KSZ_Profile, KSZ_Isothermal, KSZ_NFW
|
||||
|
31
python_sample/ksz/constants.py
Normal file
31
python_sample/ksz/constants.py
Normal file
|
@ -0,0 +1,31 @@
|
|||
import numpy as np
|
||||
|
||||
Mpc=3.08e22
|
||||
rhoc = 1.8864883524081933e-26 # m^(-3)
|
||||
sigmaT = 6.6524e-29
|
||||
mp = 1.6726e-27
|
||||
lightspeed = 299792458.
|
||||
v_unit = 1e3 # Unit of 1 km/s
|
||||
T_cmb=2.725
|
||||
h = 0.71
|
||||
Y = 0.245 #The Helium abundance
|
||||
Omega_matter = 0.26
|
||||
Omega_baryon=0.0445
|
||||
|
||||
G=6.67e-11
|
||||
MassSun=2e30
|
||||
frac_electron = 1.0 # Hmmmm
|
||||
frac_gas_galaxy = 0.14
|
||||
mu = 1/(1-0.5*Y)
|
||||
|
||||
tmpp_cat={'Msun':3.29,'alpha':-0.7286211634758224,'Mstar':-23.172904033796893,'PhiStar':0.0113246633636846,'lbar':393109973.22508669}
|
||||
|
||||
baryon_fraction = Omega_baryon / Omega_matter
|
||||
|
||||
ksz_normalization = T_cmb*sigmaT*v_unit/(lightspeed*mu*mp) * baryon_fraction
|
||||
rho_mean_matter = Omega_matter * (3*(100e3/Mpc)**2/(8*np.pi*G))
|
||||
Lbar = tmpp_cat['lbar'] / Mpc**3
|
||||
M_over_L_galaxy = rho_mean_matter / Lbar
|
||||
|
||||
|
||||
del np
|
149
python_sample/ksz/gal_prof.py
Normal file
149
python_sample/ksz/gal_prof.py
Normal file
|
@ -0,0 +1,149 @@
|
|||
from .constants import *
|
||||
|
||||
# -----------------------------------------------------------------------------
|
||||
# Generic profile generator
|
||||
# -----------------------------------------------------------------------------
|
||||
|
||||
class KSZ_Profile(object):
|
||||
R_star= 0.050 # 15 kpc/h
|
||||
L_gal0 = 10**(0.4*(tmpp_cat['Msun']-tmpp_cat['Mstar']))
|
||||
|
||||
def __init__(self):
|
||||
self.rGalaxy = 1.0
|
||||
|
||||
def evaluate_profile(r):
|
||||
raise NotImplementedError("Abstract function")
|
||||
|
||||
def projected_profile(cos_theta,angularDistance):
|
||||
|
||||
idx = np.where(cos_theta > 0)[0]
|
||||
tan_theta_2 = 1/(cos_theta[idx]**2) - 1
|
||||
tan_theta_2_max = (self.rGalaxy/angularDistance)**2
|
||||
tan_theta_2_min = (R_star/angularDistance)**2
|
||||
|
||||
idx0 = np.where((tan_theta_2 < tan_theta_2_max))
|
||||
idx = idx[idx0]
|
||||
tan_theta_2 = tan_theta_2[idx0]
|
||||
tan_theta = np.sqrt(tan_theta_2)
|
||||
|
||||
r = (tan_theta*DA)
|
||||
|
||||
m,idx_mask = self.evaluate_profile(r)
|
||||
|
||||
idx_mask = np.append(idx_mask,np.where(tan_theta_2<tan_theta_2_min)[0])
|
||||
idx_mask = np.append(idx_mask,[tan_theta_2.argmin()])
|
||||
|
||||
return idx,idx_mask,m
|
||||
|
||||
|
||||
# -----------------------------------------------------------------------------
|
||||
# Isothermal profile generator
|
||||
# -----------------------------------------------------------------------------
|
||||
|
||||
class KSZ_Isothermal(KSZ_Profile):
|
||||
sigma_FP=160e3 #m/s
|
||||
R_innergal = 0.226
|
||||
|
||||
def __init__(self, Lgal, x, y=0.0):
|
||||
"Support for Isothermal profile"
|
||||
|
||||
super(KSZ_Isothermal,self).__init__()
|
||||
|
||||
self.R_gal = 0.226 * x
|
||||
self.R_innergal *= y
|
||||
|
||||
self.rho0 = self.sigma_FP**2/(2*np.pi*G) # * (Lgal/L_gal0)**(2./3)
|
||||
self.rGalaxy = self.R_gal*(Lgal/self.L_gal0)**(1./3)
|
||||
self.rInnerGalaxy = self.R_innergal*(Lgal/self.L_gal0)**(1./3)
|
||||
# self._prepare()
|
||||
|
||||
def _prepare(self, x_min, x_max):
|
||||
from scipy.integrate import quad
|
||||
from numpy import sqrt, log10
|
||||
|
||||
lmin = log10(x_min)
|
||||
lmax = log10(x_max)
|
||||
|
||||
x = 10**(np.arange(100)*(lmax-lmin)/100.+lmin)
|
||||
profile = np.empty(x.size)
|
||||
|
||||
nu_tilde = lambda u: (1/(u**2*(1+u)))
|
||||
|
||||
for i in range(x.size):
|
||||
profile[i] = 2*quad(lambda y: (nu_tilde(sqrt(x[i]**2+y**2))), 0, args.x)[0]
|
||||
|
||||
self._table = x,profile
|
||||
|
||||
|
||||
def evaluate_profile(r):
|
||||
rho0, rGalaxy, rInner = self.rho0, self.rGalaxy, self.rInnerGalaxy
|
||||
|
||||
Q=rho0*2/r*np.arctan(np.sqrt((rGalaxy/r)**2 - 1))/Mpc
|
||||
# Q[r<rInner] = 0
|
||||
return Q,np.where(r<rInner)[0]
|
||||
|
||||
|
||||
# -----------------------------------------------------------------------------
|
||||
# NFW profile generator
|
||||
# -----------------------------------------------------------------------------
|
||||
|
||||
class KSZ_NFW(KSZ_Profile):
|
||||
""" Support for NFW profile
|
||||
"""
|
||||
|
||||
def __init__(self,x,y=0.0):
|
||||
from numpy import log, pi
|
||||
|
||||
if 'pre_nfw' not in self:
|
||||
self._prepare()
|
||||
|
||||
kiso = KSZ_Isothermal(x,y)
|
||||
r_is = kiso.rGalaxy
|
||||
rho_is = kiso.rho0
|
||||
r_inner = kiso.rInnerGalaxy
|
||||
|
||||
self.Mgal = rho_is*4*pi*(r_is/args.x)*Mpc #Lgal*M_over_L_galaxy
|
||||
self.Rvir = r_is/x
|
||||
|
||||
cs = self._get_concentration(Mgal)
|
||||
|
||||
self.rs = Rvir/cs
|
||||
b = (log(1.+cs)-cs/(1.+cs))
|
||||
|
||||
self.rho_s = Mgal/(4*pi*b*(rs*Mpc)**3)
|
||||
|
||||
|
||||
def _prepare(self, _x_min=1e-4, _x_max=1e4):
|
||||
from scipy.integrate import quad
|
||||
from numpy import sqrt, log10
|
||||
from scipy.interpolate import interp1d
|
||||
|
||||
lmin = log10(x_min)
|
||||
lmax = log10(x_max)
|
||||
|
||||
x = 10**(np.arange(100)*(lmax-lmin)/100.+lmin)
|
||||
profile = np.empty(x.size)
|
||||
|
||||
nu_tilde = lambda u: (1/(u*(1+u)**2))
|
||||
|
||||
for i in range(x.size):
|
||||
if x[i] < args.x:
|
||||
profile[i] = 2*quad(lambda y: (nu_tilde(sqrt(x[i]**2+y**2))), 0, np.sqrt((args.x)**2-x[i]**2))[0]
|
||||
else:
|
||||
profile[i] = 0
|
||||
|
||||
# Insert the interpolator into the class definition
|
||||
KSZ_NFW.pre_nfw = self.pre_nfw = interp1d(x,prof)
|
||||
|
||||
def _get_concentration(self, Mvir):
|
||||
from numpy import exp, log
|
||||
|
||||
return exp(0.971 - 0.094*log(Mvir/(1e12*MassSun)))
|
||||
|
||||
def evaluate_profile(r):
|
||||
cs = self._get_concentration(self.Mvir)
|
||||
rs = self.Rvir/cs
|
||||
|
||||
return self.rho_s*rs*Mpc*self.pre_nfw(r/rs),np.array([],dtype=int)
|
||||
|
||||
|
Loading…
Add table
Add a link
Reference in a new issue