Photometric Redshifts

frb.halos.photoz.get_des_data(coords: SkyCoord, radius: Quantity = <Quantity 15. arcmin>, starbright: float = 17, starflagval: float = 0.9, gaiacat: str = None, write: bool = False, outfile: str = None) Table[source]

Download photometry for galaxies within an FRB field. :param coords: Coordinates of the center of a cone search. :type coords: SkyCoord :param radius: Radius of cone search. :type radius: Quantity, optional :param starbright: Lower limit of r band mag. Objects brighter

than this will be removed.

Parameters:
  • starflagval (float, optional) – Upper limit for a morphology-based classifier flag. Objects more point-like (i.e. higher value) will be filtered out.

  • gaicat (str, optional) – Optional file with gaia catalog of stars within the same search radius. These stars will be removed. Must contain at least two columns: “ra” and “dec”. The values must be in decimal degrees and the column names are case sensitive.

  • write (bool, optional) – Write output table to file?

  • outfile (str, optional) – Path to the output file. If not given and write is True, the table will be written to “photom_cat_J{coords}_{radius}arcmin.fits” in the current working directory.

Returns:

Table of DES galaxies within the search radius.

Return type:

des_data (Table)

frb.halos.photoz.mhalo_lookup_tables(z_grid: list, datafolder: str = 'data', n_cores: int = 8)[source]

For each z in z_grid, produces a fits file containing m_halo values corresponding to a fixed grid of m_star values. The values are produced by sampling the Moster+13 SHMR relation. The fits files can then be used to produce interpolation functions of the moments of the m_halo distribution (e.g. mean, std.dev) as a function of redshift and log_mstar. :param z_grid: List of redshift values to be sampled. :type z_grid: list or np.ndarray :param datafolder: Path to the directory where the results will be stored. :type datafolder: str, optional :param n_cores: Number of CPU threads used for parallel processing. :type n_cores: int, optional

frb.halos.photoz.dm_grid(frb_z: float, n_z: int = 100, n_o: int = 100, n_m: int = 100, max_log_mhalo: float = 12.8, outdir: str = 'data', outfile: str = None) None[source]

Produce DM estimates for a 3D grid of redshift, offsets and log_halo_masses and write them to disk. :param frb_z: frb redshift :type frb_z: float :param n_z: size of the redshift grid. i.e. np.linspace(0, frb_z, n_z) :type n_z: int, optional :param n_o: size of the offset grid. i.e. np.linspace(0, 600, n_o) :type n_o: int, optional :param n_m: size of the log_halo_mass grid. i.e. np.linspace(8, 16, n_m) :type n_m: int, optional :param max_log_mhalo: DM for halo masses larger than this are currently

set to -99.0 to prevent weirdly large DM contributions from galactic halos.

Parameters:
  • outdir (str, optional) – data directory to store results

  • outfile (str, optional) – name of results .npz file (within outdir).

frb.halos.photoz.dm_for_all_galaxies(frb: FRB, input_catfile: str, datafolder: str, n_cores: int = 8, n_gals: int = None)[source]

Produce DM estimates for all the galaxies provided by the user. Creates two files : “DM_halos_zdraws.npz” which contains all the redshift draws used for the DM realizations and “DM_halos_final.npz” which contains the DM realizations themselves. Each row in each of these files corresponds to one galaxy and each z draw corresponds to 1000 DM realizations for a galaxy. :param frb: The FRB object of interest. :type frb: FRB :param input_catfile: Path to the input catalog of photometry. Assumed

to be from DES for now.

Parameters:
  • datafolder (str) – Path to the folder in which results will be saved.

  • n_cores (int, optional) – Number of CPU threads to be used for computation.

  • n_gals (int, optional) – Limit analysis to n_gals galaxies for testing purposes.

Module for photometric redshift-based halo DM calculations and galaxy analysis.

This module combines photometric redshift estimates with halo models to calculate dispersion measure contributions from galaxy halos along FRB sightlines.

Constants

frb.halos.photoz.DEFAULT_DATA_FOLDER = 'data'

str(object=’’) -> str str(bytes_or_buffer[, encoding[, errors]]) -> str

Create a new string object from the given object. If encoding or errors is specified, then the object must expose a data buffer that will be decoded using the given encoding and error handler. Otherwise, returns the result of object.__str__() (if defined) or repr(object). encoding defaults to sys.getdefaultencoding(). errors defaults to ‘strict’.

Default data folder path: “data”

Functions

get_des_data

frb.halos.photoz.get_des_data(coords: SkyCoord, radius: Quantity = <Quantity 15. arcmin>, starbright: float = 17, starflagval: float = 0.9, gaiacat: str = None, write: bool = False, outfile: str = None) Table[source]

Download photometry for galaxies within an FRB field. :param coords: Coordinates of the center of a cone search. :type coords: SkyCoord :param radius: Radius of cone search. :type radius: Quantity, optional :param starbright: Lower limit of r band mag. Objects brighter

than this will be removed.

Parameters:
  • starflagval (float, optional) – Upper limit for a morphology-based classifier flag. Objects more point-like (i.e. higher value) will be filtered out.

  • gaicat (str, optional) – Optional file with gaia catalog of stars within the same search radius. These stars will be removed. Must contain at least two columns: “ra” and “dec”. The values must be in decimal degrees and the column names are case sensitive.

  • write (bool, optional) – Write output table to file?

  • outfile (str, optional) – Path to the output file. If not given and write is True, the table will be written to “photom_cat_J{coords}_{radius}arcmin.fits” in the current working directory.

Returns:

Table of DES galaxies within the search radius.

Return type:

des_data (Table)

Download photometry for galaxies within an FRB field.

Parameters:

  • coords : SkyCoord - Center coordinates for cone search

  • radius : Quantity, optional - Search radius (default: 15 arcmin)

  • starbright : float, optional - Lower r-band magnitude limit (default: 17)

  • starflagval : float, optional - Star classification upper limit (default: 0.9)

  • gaiacat : str, optional - Gaia catalog file for star removal

  • write : bool, optional - Write output table to file (default: False)

  • outfile : str, optional - Output filename

Returns:

  • des_data : Table - DES galaxies within search radius

dm_grid

frb.halos.photoz.dm_grid(frb_z: float, n_z: int = 100, n_o: int = 100, n_m: int = 100, max_log_mhalo: float = 12.8, outdir: str = 'data', outfile: str = None) None[source]

Produce DM estimates for a 3D grid of redshift, offsets and log_halo_masses and write them to disk. :param frb_z: frb redshift :type frb_z: float :param n_z: size of the redshift grid. i.e. np.linspace(0, frb_z, n_z) :type n_z: int, optional :param n_o: size of the offset grid. i.e. np.linspace(0, 600, n_o) :type n_o: int, optional :param n_m: size of the log_halo_mass grid. i.e. np.linspace(8, 16, n_m) :type n_m: int, optional :param max_log_mhalo: DM for halo masses larger than this are currently

set to -99.0 to prevent weirdly large DM contributions from galactic halos.

Parameters:
  • outdir (str, optional) – data directory to store results

  • outfile (str, optional) – name of results .npz file (within outdir).

Produce DM estimates for a 3D grid of redshift, offsets and halo masses.

Parameters:

  • frb_z : float - FRB redshift

  • n_z : int, optional - Size of redshift grid (default: 100)

  • n_o : int, optional - Size of offset grid (default: 100)

  • n_m : int, optional - Size of halo mass grid (default: 100)

  • max_log_mhalo : float, optional - Maximum log halo mass (default: 12.8)

  • outdir : str, optional - Output directory (default: DEFAULT_DATA_FOLDER)

  • outfile : str, optional - Output .npz filename

Creates a 3D interpolation grid for DM calculations with dimensions:

  • Redshift: np.linspace(0, frb_z, n_z)

  • Offsets: np.linspace(0, 600, n_o) kpc

  • Halo masses: np.linspace(8, 16, n_m) log solar masses

mhalo_lookup_tables

frb.halos.photoz.mhalo_lookup_tables(z_grid: list, datafolder: str = 'data', n_cores: int = 8)[source]

For each z in z_grid, produces a fits file containing m_halo values corresponding to a fixed grid of m_star values. The values are produced by sampling the Moster+13 SHMR relation. The fits files can then be used to produce interpolation functions of the moments of the m_halo distribution (e.g. mean, std.dev) as a function of redshift and log_mstar. :param z_grid: List of redshift values to be sampled. :type z_grid: list or np.ndarray :param datafolder: Path to the directory where the results will be stored. :type datafolder: str, optional :param n_cores: Number of CPU threads used for parallel processing. :type n_cores: int, optional

For each redshift in z_grid, produces files containing halo mass values corresponding to stellar masses.

Parameters:

  • z_grid : list or ndarray - Redshift values to sample

  • datafolder : str, optional - Storage directory (default: DEFAULT_DATA_FOLDER)

  • n_cores : int, optional - CPU threads for parallel processing (default: 8)

Values are produced by sampling the Moster+13 stellar-halo mass relation (SHMR).

_mhalo_lookup_table

frb.halos.photoz._mhalo_lookup_table(z: float, npz_out: str = 'm_halo_realizations', n_cores: int = 8)[source]

For a given z, produce realizations of m_halo for relevant m_star values using only the uncertainty in the SHMR relation. Internal function. Use directly if you know what you’re doing. :param z: redshift :type z: float :param npz_out: output .npz file path. :type npz_out: str, optional :param n_cores: Number of CPU threads used for parallel processing. :type n_cores: int, optional

Internal function to create halo mass lookup tables for a single redshift.

Parameters:

  • z : float - Redshift

  • npz_out : str, optional - Output .npz file path

  • n_cores : int, optional - CPU threads (default: 8)

Note

This is an internal function. Use mhalo_lookup_tables() directly if you know what you’re doing.

_instantiate_intepolators

frb.halos.photoz._instantiate_intepolators(datafolder: str = 'data', dmfilename: str = None, frb_name: str = 'FRB180924') list[source]

Produce interpolator functions for key quantities required for the analysis. :param datfolder: Folder where the interpolation data files exist :type datfolder: str, optional :param dmfilename: file name (within datafolder) for the DM interpolation data. :type dmfilename: str, optional :param frb_name: Assumes “FRB180924” by default. :type frb_name: str, optional

Returns:

DM(z, offset_kpc, log_mhalo) mean_interp (interp2d): <log_mhalo(log_mstar, z)> (based on SHMR) stddev_interp (interp2d): std.dev. log_mhalo(log_mstar, z) (based on SHMR) ang_dia_interp (interp1d): angular_diameter_distance(z) (default Repo cosmology)

Return type:

dm_interpolator (RegularGridInterpolator)

Produce interpolator functions for key quantities required for the analysis.

Parameters:

  • datafolder : str, optional - Folder with interpolation data files (default: DEFAULT_DATA_FOLDER)

  • dmfilename : str, optional - DM interpolation data filename

  • frb_name : str, optional - FRB name (default: “FRB180924”)

Returns:

  • dm_interpolator : RegularGridInterpolator - DM(z, offset_kpc, log_mhalo)

  • mean_interp : interp2d - <log_mhalo(log_mstar, z)> based on SHMR

  • stddev_interp : interp2d - std.dev. log_mhalo(log_mstar, z) based on SHMR

  • ang_dia_interp : interp1d - angular_diameter_distance(z) for default cosmology

_mhalo_realizations

frb.halos.photoz._mhalo_realizations(log_mstar: float, log_mstar_err: float, z: float, mean_interp: interp2d, stddev_interp: interp2d, n_mstar: int = 100, n_norm: int = 10, max_log_mhalo: float = 12.8) ndarray[source]

Using the lookup tables generated (see function mhalo_lookup_tables), produce realiztions of mhalo. This takes into account both the stellar mass uncertainty and the uncertainty in the SMHR relation from Moster+13. :param log_mstar: log stellar mass in M_sun. :type log_mstar: float :param log_mstar_err: log error in log_mstar :type log_mstar_err: float :param z: redshift :type z: float :param mean_interp: <log_mhalo(log_mstar, z)> (based on SHMR) :type mean_interp: interp2d :param stddev_interp: std.dev. log_mhalo(log_mstar, z) (based on SHMR) :type stddev_interp: interp2d :param n_mstrar: Number of m_star samples to be produced. :type n_mstrar: int, optional :param n_norm: Number of m_halo samples for each m_star sample. :type n_norm: int, optional :param max_log_mhalo: Maximum allowed log halo mass. log halo masses

are capped artificially to this value if any exceed.

Returns:

log_mhalo realizations.

Return type:

mhalo_reals (np.ndarray)

Generate halo mass realizations using lookup tables, accounting for both stellar mass uncertainty and SHMR scatter.

Parameters:

  • log_mstar : float - log stellar mass in M_sun

  • log_mstar_err : float - log error in log_mstar

  • z : float - Redshift

  • mean_interp : interp2d - <log_mhalo(log_mstar, z)> interpolator

  • stddev_interp : interp2d - std.dev. log_mhalo(log_mstar, z) interpolator

  • n_mstar : int, optional - Number of stellar mass samples (default: 100)

  • n_norm : int, optional - Number of normal distribution samples (default: 10)

  • max_log_mhalo : float, optional - Maximum log halo mass (default: 12.8)

Returns:

  • ndarray - Halo mass realizations

_dm_pdf

frb.halos.photoz._dm_pdf(cigale_tab: Table, eazy_outdir: str, mean_interp: interp2d, stddev_interp: interp2d, ang_dia_interp: interp1d, dm_interpolator: RegularGridInterpolator, n_cores: int = 8)[source]

For a given galaxy, compute its PDF of DM from the CIGALE and EAZY inputs. :param cigale_tab: On of the groups

from the full cigale result. This group contains data on only one galaxy at various assumed redshifts.

Parameters:
  • eazy_outdir (str) – Path to the directory with EAZY output

  • mean_interp (interp2d) – <log_mhalo(log_mstar, z)> (based on SHMR)

  • stddev_interp (interp2d) – std.dev. log_mhalo(log_mstar, z) (based on SHMR)

  • ang_dia_interp (interp1d) – angular_diameter_distance(z) (default Repo cosmology)

  • dm_interpolator (RegularGridInterpolator) – DM(z, offset_kpc, log_mhalo)

  • n_cores (int, optional) – Number of CPU threads to use.

Returns:

Array containing DM realizations for the galaxy. z_draws (np.ndarray): Array containing redshift draws from which dm_values were produced.

Return type:

dm_values (np.ndarray)

Calculate DM realizations for a galaxy using photometric redshift and stellar mass estimates.

Parameters:

  • cigale_tab : Table - CIGALE results for the galaxy

  • eazy_outdir : str - EAZY output directory

  • mean_interp : interp2d - Mean SHMR interpolator

  • stddev_interp : interp2d - SHMR scatter interpolator

  • ang_dia_interp : interp1d - Angular diameter distance interpolator

  • dm_interpolator : RegularGridInterpolator - DM interpolator

  • n_cores : int, optional - CPU threads

Returns:

  • dm_values : ndarray - DM realizations for the galaxy

  • z_draws : ndarray - Redshift draws used for DM calculations

full_analysis

Dependencies

Required packages for full functionality:

  • pathos - For multiprocessing: pip install pathos

  • progressbar2 - For progress tracking: pip install progressbar2

  • threedhst - For EAZY output: pip install threedhst

Examples

Basic DES data retrieval:

from frb.halos.photoz import get_des_data
from astropy.coordinates import SkyCoord
from astropy import units as u

# Define search parameters
center = SkyCoord('12h34m56s', '+12d34m56s')
radius = 10 * u.arcmin

# Get DES photometry
des_cat = get_des_data(center, radius=radius, write=True)
print(f"Found {len(des_cat)} galaxies")

Create DM interpolation grid:

from frb.halos.photoz import dm_grid

# Create 3D DM grid for FRB at z=0.5
frb_redshift = 0.5
dm_grid(frb_redshift, n_z=50, n_o=50, n_m=50,
       outdir='./halo_data/', outfile='dm_grid.npz')

Full analysis workflow:

from frb.frb import FRB
from frb.halos.photoz import full_analysis

# Load FRB
frb_obj = FRB.by_name('FRB20180924B')

# Run complete halo DM analysis
full_analysis(frb_obj,
              input_catfile='field_photometry.fits',
              datafolder='./analysis_results/',
              n_cores=8,
              n_gals=100)  # Limit for testing

Setup interpolators:

from frb.halos.photoz import _instantiate_intepolators

# Create interpolation functions
dm_interp, mean_interp, std_interp, ang_interp = _instantiate_intepolators(
    datafolder='./halo_data/',
    dmfilename='dm_grid.npz'
)

# Use interpolators for DM calculations
import numpy as np
z_test = 0.3
offset_test = 100.0  # kpc
mhalo_test = 12.0   # log solar masses

dm_value = dm_interp((z_test, offset_test, mhalo_test))
print(f"DM contribution: {dm_value:.2f} pc/cm³")

Generate halo mass lookup tables:

from frb.halos.photoz import mhalo_lookup_tables
import numpy as np

# Create lookup tables for multiple redshifts
z_array = np.linspace(0.1, 1.0, 10)
mhalo_lookup_tables(z_array,
                   datafolder='./lookup_tables/',
                   n_cores=4)

Advanced usage with custom parameters:

from frb.halos.photoz import get_des_data, dm_grid
from astropy.coordinates import SkyCoord
from astropy import units as u

# Custom DES query with star removal
coords = SkyCoord(ra=123.45, dec=-23.45, unit='deg')

des_data = get_des_data(
    coords,
    radius=20*u.arcmin,
    starbright=16.0,        # Remove bright stars
    starflagval=0.8,        # More aggressive star removal
    gaiacat='gaia_stars.csv', # Additional star catalog
    write=True,
    outfile='frb_field_photometry.fits'
)

# High-resolution DM grid
dm_grid(frb_z=1.2,
       n_z=200, n_o=150, n_m=120,  # Higher resolution
       max_log_mhalo=13.5,          # Include more massive halos
       outdir='./high_res_grid/',
       outfile='dm_grid_hires.npz')

Integration with FRB analysis:

from frb.frb import FRB
from frb.halos.photoz import get_des_data, full_analysis
from astropy import units as u

# Load FRB and get field data
frb = FRB.by_name('FRB20180924B')

# Get photometry around FRB location
field_data = get_des_data(frb.coord, radius=15*u.arcmin)

# Save for analysis
field_data.write('frb_field.fits', overwrite=True)

# Run complete halo DM analysis
full_analysis(frb, 'frb_field.fits', './results/', n_cores=8)

print("Halo DM analysis complete!")
print("Results saved in ./results/DM_halos_final.npz")