-
Notifications
You must be signed in to change notification settings - Fork 14
Feature/cored nfw #266
New issue
Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.
By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.
Already on GitHub? Sign in to your account
Feature/cored nfw #266
Changes from all commits
39aab61
3397bd4
c0e0e98
8748cc5
19fb509
efbb0d9
b9d974e
a1a3e15
a779a76
0579846
bb41d60
36750ab
189d401
700f4fc
728bcb8
e2791a5
1505abf
980e655
16d2dc2
1dcf3d0
cfea2c0
7e8bdda
5daf88d
791d7f9
e719e44
ec3f05c
0a82598
fec2847
45162d1
1de95e6
c50d7f0
bfde04c
ee1bc3c
7dea916
3f1f8cc
ed50465
34efa98
aa9ac34
d867d6a
dec498c
File filter
Filter by extension
Conversations
Jump to
Diff view
Diff view
There are no files selected for viewing
| Original file line number | Diff line number | Diff line change |
|---|---|---|
| @@ -0,0 +1,51 @@ | ||
| cNFWSph: | ||
| centre_0: | ||
| type: Gaussian | ||
| mean: 0.0 | ||
| sigma: 0.1 | ||
| width_modifier: | ||
| type: Absolute | ||
| value: 0.05 | ||
| limits: | ||
| lower: -inf | ||
| upper: inf | ||
| centre_1: | ||
| type: Gaussian | ||
| mean: 0.0 | ||
| sigma: 0.1 | ||
| width_modifier: | ||
| type: Absolute | ||
| value: 0.05 | ||
| limits: | ||
| lower: -inf | ||
| upper: inf | ||
| kappa_s: | ||
| type: Uniform | ||
| lower_limit: 0.0 | ||
| upper_limit: 1.0 | ||
| width_modifier: | ||
| type: Relative | ||
| value: 0.2 | ||
| limits: | ||
| lower: 0.0 | ||
| upper: inf | ||
| scale_radius: | ||
| type: Uniform | ||
| lower_limit: 0.0 | ||
| upper_limit: 30.0 | ||
| width_modifier: | ||
| type: Relative | ||
| value: 0.2 | ||
| limits: | ||
| lower: 0.0 | ||
| upper: inf | ||
| core_radius: | ||
| type: Uniform | ||
| lower_limit: 0.0 | ||
| upper_limit: 15.0 | ||
| width_modifier: | ||
| type: Relative | ||
| value: 0.2 | ||
| limits: | ||
| lower: 0.0 | ||
| upper: inf | ||
| Original file line number | Diff line number | Diff line change |
|---|---|---|
| @@ -0,0 +1,61 @@ | ||
| cNFWMCRLudlowSph: | ||
| centre_0: | ||
| type: Gaussian | ||
| mean: 0.0 | ||
| sigma: 0.1 | ||
| width_modifier: | ||
| type: Absolute | ||
| value: 0.05 | ||
| limits: | ||
| lower: -inf | ||
| upper: inf | ||
| centre_1: | ||
| type: Gaussian | ||
| mean: 0.0 | ||
| sigma: 0.1 | ||
| width_modifier: | ||
| type: Absolute | ||
| value: 0.05 | ||
| limits: | ||
| lower: -inf | ||
| upper: inf | ||
| mass_at_200: | ||
| type: LogUniform | ||
| lower_limit: 100000000.0 | ||
| upper_limit: 1000000000000000.0 | ||
| width_modifier: | ||
| type: Relative | ||
| value: 0.5 | ||
| limits: | ||
| lower: 0.0 | ||
| upper: inf | ||
| f_c: | ||
| type: Uniform | ||
| lower_limit: 0.0001 | ||
| upper_limit: 0.5 | ||
| width_modifier: | ||
| type: Relative | ||
| value: 0.2 | ||
| limits: | ||
| lower: 0.0001 | ||
| upper: inf | ||
| redshift_object: | ||
| type: Uniform | ||
| lower_limit: 0.0 | ||
| upper_limit: 1.0 | ||
| width_modifier: | ||
| type: Relative | ||
| value: 0.5 | ||
| limits: | ||
| lower: 0.0 | ||
| upper: inf | ||
| redshift_source: | ||
| type: Uniform | ||
| lower_limit: 0.0 | ||
| upper_limit: 1.0 | ||
| width_modifier: | ||
| type: Relative | ||
| value: 0.5 | ||
| limits: | ||
| lower: 0.0 | ||
| upper: inf |
| Original file line number | Diff line number | Diff line change |
|---|---|---|
| @@ -0,0 +1,175 @@ | ||
| import numpy as np | ||
|
|
||
| from typing import Tuple | ||
|
|
||
| from autogalaxy.profiles.mass.dark.abstract import AbstractgNFW | ||
|
|
||
| import autoarray as aa | ||
|
|
||
| class cNFWSph(AbstractgNFW): | ||
| def __init__( | ||
| self, | ||
| centre: Tuple[float, float] = (0.0, 0.0), | ||
| kappa_s: float = 0.05, | ||
| scale_radius: float = 1.0, | ||
| core_radius: float = 0.01, | ||
| ): | ||
| """ | ||
| Represents a spherical cored NFW density distribution | ||
|
|
||
| Parameters | ||
| ---------- | ||
| centre | ||
| The (y,x) arc-second coordinates of the profile centre. | ||
| kappa_s | ||
| The overall normalization of the dark matter halo \| | ||
| (kappa_s = (rho_0 * D_d * scale_radius)/lensing_critical_density) | ||
| scale_radius | ||
| The cored NFW scale radius `theta_s`, as an angle on the sky in arcseconds. | ||
| core_radius | ||
| The cored NFW core radius `theta_c`, as an angle on the sky in arcseconds. | ||
| """ | ||
|
|
||
| super().__init__( | ||
| centre=centre, | ||
| ell_comps=(0.0, 0.0)) | ||
|
|
||
| self.kappa_s = kappa_s | ||
| self.scale_radius = scale_radius | ||
| self.core_radius = core_radius | ||
|
|
||
|
|
||
| @aa.grid_dec.to_vector_yx | ||
| @aa.grid_dec.transform | ||
| def deflections_yx_2d_from(self, grid: aa.type.Grid2DLike, xp=np, **kwargs): | ||
| """ | ||
| Calculate the deflection angles on a grid of (y,x) arc-second coordinates. | ||
|
|
||
| The input grid of (y,x) coordinates are transformed to a coordinate system centred on the profile centre with | ||
| and rotated based on the position angle defined from its `ell_comps` (this is described fully below). | ||
|
|
||
| The numerical backend can be selected via the ``xp`` argument, allowing this | ||
| method to be used with both NumPy and JAX (e.g. inside ``jax.jit``-compiled | ||
| code). This is described fully later in this example. | ||
|
|
||
| Parameters | ||
| ---------- | ||
| grid | ||
| The grid of (y,x) arc-second coordinates the deflection angles are computed on. | ||
| xp | ||
| The numerical backend to use, either `numpy` or `jax.numpy`. | ||
| """ | ||
| theta = self.radial_grid_from(grid=grid, xp=xp, **kwargs).array | ||
| theta = xp.maximum(theta, 1e-8) | ||
|
|
||
| factor = ( | ||
| 4.0 | ||
| * self.kappa_s | ||
| * self.scale_radius**2 | ||
| ) | ||
|
|
||
| deflection_r = ( | ||
| factor | ||
| * (self.F_func(theta, self.scale_radius, xp=xp) - self.F_func(theta, self.core_radius, xp=xp) | ||
| - (self.scale_radius - self.core_radius) * self.dev_F_func(theta, self.scale_radius, xp=xp) | ||
| ) | ||
| / (theta * (self.scale_radius - self.core_radius)**2) | ||
| ) | ||
|
Comment on lines
+71
to
+77
|
||
|
|
||
|
|
||
| return self._cartesian_grid_via_radial_from( | ||
| grid=grid, | ||
| radius=deflection_r, | ||
| xp=xp, | ||
| **kwargs, | ||
| ) | ||
|
|
||
| def F_func(self, theta, radius, xp=np): | ||
|
|
||
| F = theta * 0.0 | ||
|
|
||
| # theta < radius | ||
| mask1 = (theta > 0) & (theta < radius) | ||
|
|
||
| # theta > radius | ||
| mask2 = theta > radius | ||
|
|
||
| F = xp.where( | ||
| mask1, | ||
| ( | ||
| radius / 2 * xp.log(2 * radius / theta) | ||
| - xp.sqrt(radius ** 2 - theta ** 2) | ||
| * xp.arctanh(xp.sqrt((radius - theta) / (radius + theta))) | ||
| ), | ||
| F, | ||
| ) | ||
|
|
||
| F = xp.where( | ||
| mask2, | ||
| ( | ||
| radius / 2 * xp.log(2 * radius / theta) | ||
| + xp.sqrt(theta ** 2 - radius ** 2) | ||
| * xp.arctan(xp.sqrt((theta - radius) / (theta + radius))) | ||
| ), | ||
| F, | ||
| ) | ||
|
|
||
| return 2 * radius * F | ||
|
|
||
| def dev_F_func(self, theta, radius, xp=np): | ||
|
|
||
| dev_F = theta * 0.0 | ||
|
|
||
| mask1 = (theta > 0) & (theta < radius) | ||
| mask2 = theta == radius | ||
| mask3 = theta > radius | ||
|
|
||
| dev_F = xp.where( | ||
| mask1, | ||
| ( | ||
| radius * xp.log(2 * radius / theta) | ||
| - (2 * radius ** 2 - theta ** 2) / xp.sqrt(radius ** 2 - theta ** 2) | ||
| * xp.arctanh(xp.sqrt((radius - theta) / (radius + theta))) | ||
| ), | ||
| dev_F, | ||
| ) | ||
|
|
||
| dev_F = xp.where( | ||
| mask2, | ||
| radius * (xp.log(2) - 1 / 2), | ||
| dev_F, | ||
| ) | ||
|
|
||
| dev_F = xp.where( | ||
| mask3, | ||
| ( | ||
| radius * xp.log(2 * radius / theta) | ||
| + (theta ** 2 - 2 * radius ** 2) / xp.sqrt(theta ** 2 - radius ** 2) | ||
| * xp.arctan(xp.sqrt((theta - radius) / (theta + radius))) | ||
| ), | ||
| dev_F, | ||
| ) | ||
|
|
||
| return 2 * dev_F | ||
|
|
||
| @aa.grid_dec.to_array | ||
| def convergence_2d_from(self, grid: aa.type.Grid2DLike, xp=np, **kwargs): | ||
| """ | ||
| Convergence (dimensionless surface mass density) for the cored NFW profile. | ||
| This is not yet implemented for `cNFWSph`. | ||
| """ | ||
| raise NotImplementedError( | ||
| "convergence_2d_from is not implemented for cNFWSph; a physical cNFW " | ||
| "convergence expression must be added before use." | ||
| ) | ||
|
|
||
| @aa.grid_dec.to_array | ||
| def potential_2d_from(self, grid: aa.type.Grid2DLike, xp=np, **kwargs): | ||
| """ | ||
| Lensing potential for the cored NFW profile. | ||
| This is not yet implemented for `cNFWSph`. | ||
| """ | ||
| raise NotImplementedError( | ||
| "potential_2d_from is not implemented for cNFWSph; a physical cNFW " | ||
| "potential expression must be added before use." | ||
| ) | ||
| Original file line number | Diff line number | Diff line change |
|---|---|---|
| @@ -0,0 +1,34 @@ | ||
| from typing import Tuple | ||
|
|
||
| from autogalaxy.profiles.mass.dark.cnfw import cNFWSph | ||
|
|
||
| from autogalaxy.profiles.mass.dark import mcr_util | ||
|
|
||
| class cNFWMCRLudlowSph(cNFWSph): | ||
| def __init__( | ||
| self, | ||
| centre: Tuple[float, float] = (0.0, 0.0), | ||
| mass_at_200: float = 1e9, | ||
| f_c = 0.01, | ||
| redshift_object: float = 0.5, | ||
| redshift_source: float = 1.0, | ||
| ): | ||
| self.mass_at_200 = mass_at_200 | ||
| self.f_c = f_c | ||
| self.redshift_object = redshift_object | ||
| self.redshift_source = redshift_source | ||
|
|
||
| ( | ||
| kappa_s, | ||
| scale_radius, | ||
| core_radius, | ||
| radius_at_200, | ||
| ) = mcr_util.kappa_s_scale_radius_and_core_radius_for_ludlow( | ||
| mass_at_200=mass_at_200, | ||
| scatter_sigma=0.0, | ||
| f_c=f_c, | ||
| redshift_object=redshift_object, | ||
| redshift_source=redshift_source, | ||
| ) | ||
|
|
||
| super().__init__(centre=centre, kappa_s=kappa_s, scale_radius=scale_radius, core_radius=core_radius) |
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
The prior for
core_radiusallowslower_limit: 0.0, but the cNFW deflection implementation callsF_func(theta, radius=self.core_radius)whereF_funccontains terms likeradius * log(2 * radius / theta); settingcore_radius=0will drivelog(0)and produce invalid values. To avoid this singular case, consider enforcingcore_radius > 0(or handling thecore_radius -> 0limit analytically) instead of allowing exactly zero.