Skip to content
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

Implements Friends of Friends Sampling function #460

Merged
merged 9 commits into from
Dec 1, 2023
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension


Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
109 changes: 109 additions & 0 deletions btk/sampling_functions.py
Original file line number Diff line number Diff line change
Expand Up @@ -4,7 +4,9 @@

import astropy
import numpy as np
from astropy.coordinates import SkyCoord
from astropy.table import Table
from fast3tree import find_friends_of_friends
ismael-mendoza marked this conversation as resolved.
Show resolved Hide resolved

from btk.utils import DEFAULT_SEED

Expand Down Expand Up @@ -429,3 +431,110 @@
blend_table["dec"] *= 3600

return blend_table


class FriendsOfFriendsSampling(SamplingFunction):
"""Randomly selects galaxies using Friends Of Friends clustering algorithm.

Friends of friends clustering algorithm assigns a group of object same index if one can
each member of the group is within link_distance of any other member in the group.
This sampling function explicitly uses the spatial information in the input catalog to
generate scenes of galaxies. However, blends might not always be returned as a result
if the provided `link_distance` is too small.
"""

def __init__(
self,
max_number: int = 10,
min_number: int = 2,
link_distance: int = 2.5,
stamp_size: float = 24.0,
seed: int = DEFAULT_SEED,
min_mag: float = -np.inf,
max_mag: float = 25.3,
mag_name: str = "i_ab",
):
"""Initializes the FriendsOfFriendsSampling sampling function.

Args:
max_number: Defined in parent class
min_number: Defined in parent class
link_distance: Minimum linkage distance to form a group (arcsec).
stamp_size: Size of the desired stamp (arcsec).
seed: Seed to initialize randomness for reproducibility.
min_mag: Minimum magnitude allowed in samples
max_mag: Maximum magnitude allowed in samples.
mag_name: Name of the magnitude column in the catalog.
"""
super().__init__(max_number=max_number, min_number=min_number, seed=seed)
self.stamp_size = stamp_size
self.link_distance = link_distance
self.max_number = max_number
self.min_number = min_number
self.max_mag = max_mag
self.min_mag = min_mag
self.mag_name = mag_name

Check warning on line 476 in btk/sampling_functions.py

View check run for this annotation

Codecov / codecov/patch

btk/sampling_functions.py#L469-L476

Added lines #L469 - L476 were not covered by tests

@staticmethod
def _arcsec2dist(sep, r=1.0):
ismael-mendoza marked this conversation as resolved.
Show resolved Hide resolved
"""Converts arcsec separation to a cartesian distance."""
return np.sin(np.deg2rad(sep / 3600.0 / 2.0)) * 2.0 * r

Check warning on line 481 in btk/sampling_functions.py

View check run for this annotation

Codecov / codecov/patch

btk/sampling_functions.py#L481

Added line #L481 was not covered by tests

def _precompute_friends_of_friends_index(self, table: Table):
"""Computes galaxy groups using friends of friends algorithm on cartesian coordinates."""
points = SkyCoord(table["ra"], table["dec"], unit=("deg", "deg")).cartesian.xyz.value.T
group_ids = find_friends_of_friends(

Check warning on line 486 in btk/sampling_functions.py

View check run for this annotation

Codecov / codecov/patch

btk/sampling_functions.py#L485-L486

Added lines #L485 - L486 were not covered by tests
points, self._arcsec2dist(self.link_distance), reassign_group_indices=False
)
table["friends_of_friends_id"] = group_ids

Check warning on line 489 in btk/sampling_functions.py

View check run for this annotation

Codecov / codecov/patch

btk/sampling_functions.py#L489

Added line #L489 was not covered by tests

def __call__(self, table: Table):
"""Samples galaxies from input catalog to make scene.

We assume the input catalog has `ra` and `dec` in degrees, like CATSIM does.
"""
# group galaxies by indices
if "friends_of_friends_id" not in table.colnames:
self._precompute_friends_of_friends_index(table)

Check warning on line 498 in btk/sampling_functions.py

View check run for this annotation

Codecov / codecov/patch

btk/sampling_functions.py#L497-L498

Added lines #L497 - L498 were not covered by tests

# filter by magnitude
if self.mag_name not in table.colnames:
raise ValueError(f"Catalog must have '{self.mag_name}' column.")
cond1 = table[self.mag_name] <= self.max_mag
cond2 = table[self.mag_name] > self.min_mag
cond = cond1 & cond2
blend_table = table[cond]

Check warning on line 506 in btk/sampling_functions.py

View check run for this annotation

Codecov / codecov/patch

btk/sampling_functions.py#L501-L506

Added lines #L501 - L506 were not covered by tests

# filter by number of galaxies
num_galaxies = blend_table.to_pandas()["friends_of_friends_id"].value_counts()
cond1 = num_galaxies >= self.min_number
cond2 = num_galaxies <= self.max_number
if (cond1 & cond2).sum() == 0:
raise ValueError(

Check warning on line 513 in btk/sampling_functions.py

View check run for this annotation

Codecov / codecov/patch

btk/sampling_functions.py#L509-L513

Added lines #L509 - L513 were not covered by tests
f"No groups with number of galaxies in [{self.min_number}, {self.max_number}], "
"found using a `link_distance` of {self.link_distance}"
)
group_id = np.random.choice(num_galaxies[cond1 & cond2].index)
indices = blend_table["friends_of_friends_id"] == group_id

Check warning on line 518 in btk/sampling_functions.py

View check run for this annotation

Codecov / codecov/patch

btk/sampling_functions.py#L517-L518

Added lines #L517 - L518 were not covered by tests

# check if galaxies fit in the window
ra = blend_table["ra"]
dec = blend_table["dec"]
ra_out_of_bound = (ra[indices].max() - ra[indices].min()) * 3600 > self.stamp_size
dec_out_of_bounds = (dec[indices].max() - dec[indices].min()) * 3600 > self.stamp_size
if ra_out_of_bound or dec_out_of_bounds:
raise ValueError(

Check warning on line 526 in btk/sampling_functions.py

View check run for this annotation

Codecov / codecov/patch

btk/sampling_functions.py#L521-L526

Added lines #L521 - L526 were not covered by tests
"sampled galaxies don't fit in the stamp size, increase the stamp size, "
"decrease `max_number` or `link_distance`."
)

# recenter catalog 'ra' and 'dec' to the center of the stamp
blend_table = blend_table[indices]
blend_table["ra"] = ra[indices] - ra[indices].mean()
blend_table["dec"] = dec[indices] - dec[indices].mean()

Check warning on line 534 in btk/sampling_functions.py

View check run for this annotation

Codecov / codecov/patch

btk/sampling_functions.py#L532-L534

Added lines #L532 - L534 were not covered by tests

# finally, convert to arcsec
blend_table["ra"] *= 3600
blend_table["dec"] *= 3600

Check warning on line 538 in btk/sampling_functions.py

View check run for this annotation

Codecov / codecov/patch

btk/sampling_functions.py#L537-L538

Added lines #L537 - L538 were not covered by tests

return blend_table

Check warning on line 540 in btk/sampling_functions.py

View check run for this annotation

Codecov / codecov/patch

btk/sampling_functions.py#L540

Added line #L540 was not covered by tests
23 changes: 18 additions & 5 deletions poetry.lock

Some generated files are not rendered by default. Learn more about how customized files appear on GitHub.

1 change: 1 addition & 0 deletions pyproject.toml
Original file line number Diff line number Diff line change
Expand Up @@ -44,6 +44,7 @@ galsim = ">=2.4.9"
python = "^3.8.1,<3.12"
pre-commit = "^3.3.3"
h5py = "^3.9.0"
fast3tree = "^0.4.1"

[tool.poetry.dev-dependencies]
black = ">=23.3.0"
Expand Down
Loading