__author__ = 'sibirrer'
# this file contains a class to compute the truncated Navaro-Frank-White function (Baltz et al 2009)in mass/kappa space
# the potential therefore is its integral
import numpy as np
import warnings
from lenstronomy.LensModel.Profiles.base_profile import LensProfileBase
__all__ = ['TNFW']
[docs]class TNFW(LensProfileBase):
"""
this class contains functions concerning the truncated NFW profile with a truncation function (r_trunc^2)*(r^2+r_trunc^2)
relation are: R_200 = c * Rs
"""
param_names = ['Rs', 'alpha_Rs', 'r_trunc', 'center_x', 'center_y']
lower_limit_default = {'Rs': 0, 'alpha_Rs': 0, 'r_trunc': 0, 'center_x': -100, 'center_y': -100}
upper_limit_default = {'Rs': 100, 'alpha_Rs': 10, 'r_trunc': 100, 'center_x': 100, 'center_y': 100}
def __init__(self):
"""
:param interpol: bool, if True, interpolates the functions F(), g() and h()
"""
self._s = 0.001
super(LensProfileBase, self).__init__()
[docs] def function(self, x, y, Rs, alpha_Rs, r_trunc, center_x=0, center_y=0):
"""
:param x: angular position
:param y: angular position
:param Rs: angular turn over point
:param alpha_Rs: deflection at Rs
:param center_x: center of halo
:param center_y: center of halo
:return:
"""
rho0_input = self._alpha2rho0(alpha_Rs=alpha_Rs, Rs=Rs)
#if Rs < 0.0001:
# Rs = 0.0001
x_ = x - center_x
y_ = y - center_y
R = np.sqrt(x_ ** 2 + y_ ** 2)
R = np.maximum(R, self._s * Rs)
f_ = self.nfwPot(R, Rs, rho0_input, r_trunc)
return f_
[docs] def L(self, x, tau):
"""
Logarithm that appears frequently
:param x: r/Rs
:param tau: t/Rs
:return:
"""
x = np.maximum(x, self._s)
return np.log(x * (tau + np.sqrt(tau ** 2 + x ** 2)) ** -1)
[docs] def F(self, x):
"""
Classic NFW function in terms of arctanh and arctan
:param x: r/Rs
:return:
"""
x = np.maximum(x, self._s)
if isinstance(x, np.ndarray):
nfwvals = np.zeros_like(x)
inds1 = np.where(x < 1)
inds2 = np.where(x > 1)
inds3 = np.where(x == 1)
nfwvals[inds1] = (1 - x[inds1] ** 2) ** -.5 * np.arctanh((1 - x[inds1] ** 2) ** .5)
nfwvals[inds2] = (x[inds2] ** 2 - 1) ** -.5 * np.arctan((x[inds2] ** 2 - 1) ** .5)
nfwvals[inds3] = 1
return nfwvals
elif isinstance(x, float) or isinstance(x, int):
if x == 1:
return 1
elif x == 0:
return 0
elif x < 1:
return (1 - x ** 2) ** -.5 * np.arctanh((1 - x ** 2) ** .5)
else:
return (x ** 2 - 1) ** -.5 * np.arctan((x ** 2 - 1) ** .5)
[docs] def derivatives(self, x, y, Rs=None, alpha_Rs=None, r_trunc=None, center_x=0, center_y=0):
"""
returns df/dx and df/dy of the function (integral of TNFW), which are the deflection angles
:param x: angular position (normally in units of arc seconds)
:param y: angular position (normally in units of arc seconds)
:param Rs: turn over point in the slope of the NFW profile in angular unit
:param alpha_Rs: deflection (angular units) at projected Rs
:param r_trunc: truncation radius (angular units)
:param center_x: center of halo (in angular units)
:param center_y: center of halo (in angular units)
:return: deflection angle in x, deflection angle in y
"""
rho0_input = self._alpha2rho0(alpha_Rs=alpha_Rs, Rs=Rs)
#if Rs < 0.0000001:
# Rs = 0.0000001
x_ = x - center_x
y_ = y - center_y
R = np.sqrt(x_ ** 2 + y_ ** 2)
R = np.maximum(R, self._s * Rs)
f_x, f_y = self.nfwAlpha(R, Rs, rho0_input, r_trunc, x_, y_)
return f_x, f_y
[docs] def hessian(self, x, y, Rs, alpha_Rs, r_trunc, center_x=0, center_y=0):
"""
returns d^2f/dx^2, d^2f/dxdy, d^2f/dydx, d^2f/dy^2 of the TNFW potential f
:param x: angular position (normally in units of arc seconds)
:param y: angular position (normally in units of arc seconds)
:param Rs: turn over point in the slope of the NFW profile in angular unit
:param alpha_Rs: deflection (angular units) at projected Rs
:param r_trunc: truncation radius (angular units)
:param center_x: center of halo (in angular units)
:param center_y: center of halo (in angular units)
:return: Hessian matrix of function d^2f/dx^2, d^f/dy^2, d^2/dxdy
"""
rho0_input = self._alpha2rho0(alpha_Rs=alpha_Rs, Rs=Rs)
#if Rs < 0.0001:
# Rs = 0.0001
x_ = x - center_x
y_ = y - center_y
R = np.sqrt(x_ ** 2 + y_ ** 2)
R = np.maximum(R, self._s * Rs)
kappa = self.density_2d(x_, y_, Rs, rho0_input, r_trunc)
gamma1, gamma2 = self.nfwGamma(R, Rs, rho0_input, r_trunc, x_, y_)
f_xx = kappa + gamma1
f_yy = kappa - gamma1
f_xy = gamma2
return f_xx, f_xy, f_xy, f_yy
[docs] def density(self, R, Rs, rho0, r_trunc):
"""
three dimenstional truncated NFW profile
:param R: radius of interest
:type R: float/numpy array
:param Rs: scale radius
:type Rs: float
:param rho0: density normalization (characteristic density)
:type rho0: float
:return: rho(R) density
"""
return (r_trunc ** 2 * (r_trunc ** 2 + R ** 2) ** -1) * rho0 / (R / Rs * (1 + R / Rs) ** 2)
[docs] def density_2d(self, x, y, Rs, rho0, r_trunc, center_x=0, center_y=0):
"""
projected two dimensional NFW profile (kappa*Sigma_crit)
:param R: radius of interest
:type R: float/numpy array
:param Rs: scale radius
:type Rs: float
:param rho0: density normalization (characteristic density)
:type rho0: float
:param r200: radius of (sub)halo
:type r200: float>0
:return: Epsilon(R) projected density at radius R
"""
x_ = x - center_x
y_ = y - center_y
R = np.sqrt(x_ ** 2 + y_ ** 2)
x = R * Rs ** -1
tau = float(r_trunc) * Rs ** -1
Fx = self._F(x, tau)
return 2 * rho0 * Rs * Fx
[docs] def mass_3d(self, R, Rs, rho0, r_trunc):
"""
mass enclosed a 3d sphere or radius r
:param r:
:param Ra:
:param Rs:
:return:
"""
x = R * Rs ** -1
x = np.maximum(x, self._s)
func = (r_trunc ** 2 * (-2 * x * (1 + r_trunc ** 2) + 4 * (1 + x) * r_trunc * np.arctan(x / r_trunc) -
2 * (1 + x) * (-1 + r_trunc ** 2) * np.log(Rs) + 2 * (1 + x) * (-1 + r_trunc ** 2) * np.log(Rs * (1 + x)) +
2 * (1 + x) * (-1 + r_trunc ** 2) * np.log(Rs * r_trunc) -
(1 + x) * (-1 + r_trunc ** 2) * np.log(Rs ** 2 * (x ** 2 + r_trunc ** 2)))) / (2. * (1 + x) * (1 + r_trunc ** 2) ** 2)
m_3d = 4*np.pi*Rs ** 3 * rho0 * func
return m_3d
[docs] def nfwPot(self, R, Rs, rho0, r_trunc):
"""
lensing potential of NFW profile
:param R: radius of interest
:type R: float/numpy array
:param Rs: scale radius
:type Rs: float
:param rho0: density normalization (characteristic density)
:type rho0: float
:return: Epsilon(R) projected density at radius R
"""
x = R / Rs
x = np.maximum(x, self._s)
tau = float(r_trunc) / Rs
hx = self._h(x, tau)
return 2 * rho0 * Rs ** 3 * hx
[docs] def nfwAlpha(self, R, Rs, rho0, r_trunc, ax_x, ax_y):
"""
deflection angel of NFW profile along the projection to coordinate axis
:param R: radius of interest
:type R: float/numpy array
:param Rs: scale radius
:type Rs: float
:param rho0: density normalization (characteristic density)
:type rho0: float
:param r200: radius of (sub)halo
:type r200: float>0
:param axis: projection to either x- or y-axis
:type axis: same as R
:return: Epsilon(R) projected density at radius R
"""
R = np.maximum(R, self._s * Rs)
#if isinstance(R, int) or isinstance(R, float):
# R = max(R, 0.00001)
#else:
# R[R <= 0.00001] = 0.00001
x = R / Rs
x = np.maximum(x, self._s)
tau = float(r_trunc) / Rs
gx = self._g(x, tau)
a = 4 * rho0 * Rs * gx / x ** 2
return a * ax_x, a * ax_y
[docs] def nfwGamma(self, R, Rs, rho0, r_trunc, ax_x, ax_y):
"""
shear gamma of NFW profile (times Sigma_crit) along the projection to coordinate 'axis'
:param R: radius of interest
:type R: float/numpy array
:param Rs: scale radius
:type Rs: float
:param rho0: density normalization (characteristic density)
:type rho0: float
:param r200: radius of (sub)halo
:type r200: float>0
:param axis: projection to either x- or y-axis
:type axis: same as R
:return: Epsilon(R) projected density at radius R
"""
c = 0.000001
#if isinstance(R, int) or isinstance(R, float):
# R = max(R, c)
#else:
# R[R <= c] = c
R = np.maximum(R, self._s * Rs)
x = R / Rs
#x = np.maximum(x, self._s)
#R = np.maximum(R, self._s * Rs)
tau = float(r_trunc) * Rs ** -1
gx = self._g(x, tau)
Fx = self._F(x, tau)
a = 2*rho0*Rs*(2*gx/x**2 - Fx) # /x #2*rho0*Rs*(2*gx/x**2 - Fx)*axis/x
return a * (ax_y ** 2 - ax_x ** 2) / R ** 2, -a * 2 * (ax_x * ax_y) / R ** 2
[docs] def mass_2d(self, R, Rs, rho0, r_trunc):
"""
analytic solution of the projection integral
(convergence)
:param x: R/Rs
:type x: float >0
"""
x = R / Rs
x = np.maximum(x, self._s)
tau = r_trunc / Rs
gx = self._g(x,tau)
m_2d = 4 * rho0 * Rs * R ** 2 * gx / x ** 2 * np.pi
return m_2d
def _F(self, X, tau):
"""
analytic solution of the projection integral
(convergence)
:param x: R/Rs
:type x: float >0
"""
t2 = tau ** 2
#Fx = self.F(X)
X = np.maximum(X, self._s )
_F = self.F(X)
a = t2*(t2+1)**-2
if isinstance(X, np.ndarray):
#b = (t2 + 1) * (X ** 2 - 1) ** -1 * (1 - _F)
b = np.ones_like(X)
b[X == 1] = (t2+1) * 1./3
b[X != 1] = (t2 + 1) * (X[X != 1] ** 2 - 1) ** -1 * (1 - _F[X != 1])
elif isinstance(X, float) or isinstance(X, int):
if X == 1:
b = (t2+1)* 1./3
else:
b = (t2+1)*(X**2-1)**-1*(1-_F)
else:
raise ValueError("The variable type is not compatible with the function, please use float, int or ndarray's.")
c = 2*_F
d = -np.pi*(t2+X**2)**-0.5
e = (t2-1)*(tau*(t2+X**2)**0.5)**-1*self.L(X,tau)
result = a * (b + c + d + e)
return result
def _g(self, x, tau):
"""
analytic solution of integral for NFW profile to compute deflection angel and gamma
:param x: R/Rs
:type x: float >0
"""
x = np.maximum(x, self._s)
return tau ** 2 * (tau ** 2 + 1) ** -2 * (
(tau ** 2 + 1 + 2 * (x ** 2 - 1)) * self.F(x) + tau * np.pi + (tau ** 2 - 1) * np.log(tau) +
np.sqrt(tau ** 2 + x ** 2) * (-np.pi + self.L(x, tau) * (tau ** 2 - 1) * tau ** -1))
@staticmethod
def _cos_function(x):
if isinstance(x, np.ndarray) or isinstance(x, list):
out = np.empty_like(x)
inds1 = np.where(x < 1)
inds2 = np.where(x >= 1)
out[inds1] = -np.arccosh(1 / x[inds1]) ** 2
out[inds2] = np.arccos(1 / x[inds2]) ** 2
elif isinstance(x, float) or isinstance(x, int):
if x < 1:
out = -np.arccosh(1 / x) ** 2
else:
out = np.arccos(1 / x) ** 2
else:
raise Exception('x data type '+type(x)+' not recognized.')
return out
def _h(self, x, tau):
"""
expression for the integral to compute potential
:param x: R/Rs
:param tau: r_trunc/Rs
:type x: float >0
"""
x = np.maximum(x, self._s)
u = x ** 2
t2 = tau ** 2
Lx = self.L(x, tau)
Fx = self.F(x)
return (t2 + 1) ** -2 * (
2 * t2 * np.pi * (tau - (t2 + u) ** .5 + tau * np.log(tau + (t2 + u) ** .5))
+
2 * (t2 - 1) * tau * (t2 + u) ** .5 * Lx
+
t2 * (t2 - 1) * Lx ** 2
+
4 * t2 * (u - 1) * Fx
+
t2 * (t2 - 1) * self._cos_function(x)
+
t2 * ((t2 - 1) * np.log(tau) - t2 - 1) * np.log(u)
-
t2 * (
(t2 - 1) * np.log(tau) * np.log(4 * tau) + 2 * np.log(0.5 * tau) - 2 * tau * (
tau - np.pi) * np.log(
tau * 2)))
def _alpha2rho0(self, alpha_Rs, Rs):
"""
convert angle at Rs into rho0; neglects the truncation
"""
rho0 = alpha_Rs / (4. * Rs ** 2 * (1. + np.log(1. / 2.)))
return rho0
def _rho02alpha(self, rho0, Rs):
"""
neglects the truncation
convert rho0 to angle at Rs
:param rho0:
:param Rs:
:return:
"""
alpha_Rs = rho0 * (4 * Rs ** 2 * (1 + np.log(1. / 2.)))
return alpha_Rs