ler.gw_source_population.prior_functions

Module for JIT-compiled functions used in GW source population simulations.

This module provides Numba JIT-compiled functions for efficient sampling and computation of merger rate densities, star formation rates, and mass distributions for various gravitational wave source populations including BBH, BNS, and NSBH.

Key Features:

  • Merger rate density functions for PopI/II, PopIII, and Primordial BBH

  • Star formation rate models (Madau & Dickinson 2014, Madau & Fragos 2017)

  • Mass distribution samplers (power-law Gaussian, broken power-law, bimodal)

  • JIT-compiled for high performance with Numba

Copyright (C) 2026 Hemanta Ph. Distributed under MIT License.

Module Contents

Functions

merger_rate_density_bbh_oguri2018_function(zs[, R0, ...])

Compute the merger rate density for PopI/II BBH.

merger_rate_density_bbh_popIII_ken2022_function(zs[, ...])

Compute the unnormalized merger rate density for PopIII BBH.

merger_rate_density_madau_dickinson2014_function(zs[, ...])

Compute the merger rate density for BBH using Madau & Dickinson (2014) model.

merger_rate_density_madau_dickinson_belczynski_ng_function(zs)

Compute BBH merger rate density following Ng et al. (2021).

merger_rate_density_bbh_primordial_ken2022_function(zs)

Compute the merger rate density for Primordial BBH.

sfr_madau_fragos2017_with_bbh_td(zs[, R0])

Compute the merger rate density for BBH. This is computed from star formation rate, Madau & Fragos (2017), with an additional time delay. This function is relies on pre-generated data points.

sfr_madau_dickinson2014_with_bbh_td(zs[, R0])

Compute the merger rate density for BBH. This is computed from star formation rate, Madau & Dickinson (2014), with an additional time delay. This function is relies on pre-generated data points.

sfr_madau_fragos2017_with_bns_td(zs[, R0])

Compute the merger rate density for BNS. This is computed from star formation rate, Madau & Fragos (2017), with an additional time delay. This function is relies on pre-generated data points.

sfr_madau_dickinson2014_with_bns_td(zs[, R0])

Compute the merger rate density for BNS. This is computed from star formation rate, Madau & Dickinson (2014), with an additional time delay. This function is relies on pre-generated data points.

sfr_madau_fragos2017(zs[, a, b, c, d])

Compute star formation rate using Madau & Fragos (2017) model.

sfr_madau_dickinson2014(zs[, a, b, c, d])

Compute star formation rate using Madau & Dickinson (2014) model.

binary_masses_BBH_popIII_lognormal_rvs(size[, m_min, ...])

Sample from a lognormal distribution in 2D mass space.

binary_masses_BBH_primordial_lognormal_rvs(size[, ...])

Sample from a lognormal distribution in 2D mass space.

binary_masses_BNS_bimodal_rvs(size[, w, muL, sigmaL, ...])

Sample BNS masses from bimodal Gaussian distribution.

binary_masses_NSBH_broken_powerlaw_rvs([size, mminbh, ...])

Generate NSBH mass samples from broken power-law (BH) and power-law (NS).

binary_masses_BBH_powerlaw_gaussian_rvs(size, mminbh, ...)

Generate BBH mass samples from power-law + Gaussian model with mass ratio.

available_prior_list()

Returns a list of available priors.

ler.gw_source_population.prior_functions.merger_rate_density_bbh_oguri2018_function(zs, R0=19 * 1e-09, b2=1.6, b3=2.1, b4=30)[source]

Compute the merger rate density for PopI/II BBH.

Reference: Oguri et al. (2018). The output is in detector frame and is unnormalized.

Parameters:
zsfloat or numpy.ndarray

Source redshifts.

R0float

Local merger rate density at low redshift (Mpc^-3 yr^-1).

default: 19e-9 (GWTC-4)

b2float

Fitting parameter.

default: 1.6

b3float

Fitting parameter.

default: 2.1

b4float

Fitting parameter.

default: 30

Returns:
rate_densityfloat or numpy.ndarray

Merger rate density.

Examples

>>> from ler.gw_source_population import merger_rate_density_bbh_oguri2018
>>> rate_density = merger_rate_density_bbh_oguri2018(zs=np.array([0.1]))
ler.gw_source_population.prior_functions.merger_rate_density_bbh_popIII_ken2022_function(zs, n0=19.2 * 1e-09, aIII=0.66, bIII=0.3, zIII=11.6)[source]

Compute the unnormalized merger rate density for PopIII BBH.

Reference: Ng et al. (2022). The output is in detector frame and is unnormalized.

Parameters:
zsfloat or numpy.ndarray

Source redshifts.

n0float

Normalization constant.

default: 19.2e-9

aIIIfloat

Fitting parameter.

default: 0.66

bIIIfloat

Fitting parameter.

default: 0.3

zIIIfloat

Characteristic redshift.

default: 11.6

Returns:
rate_densityfloat or numpy.ndarray

Merger rate density.

Examples

>>> from ler.gw_source_population import merger_rate_density_bbh_popIII_ken2022
>>> rate_density = merger_rate_density_bbh_popIII_ken2022(zs=np.array([0.1]))
ler.gw_source_population.prior_functions.merger_rate_density_madau_dickinson2014_function(zs, R0=19 * 1e-09, a=0.015, b=2.7, c=2.9, d=5.6)[source]

Compute the merger rate density for BBH using Madau & Dickinson (2014) model.

density(zs) ∝ (1 + zs) ** b / (1 + ((1 + zs) / c) ** d)

Reference: Eqn. 15 of https://arxiv.org/pdf/1403.0007

Parameters:
zsfloat or numpy.ndarray

Source redshifts.

R0float

Local merger rate density (Mpc^-3 yr^-1).

default: 19e-9

afloat

Normalization parameter.

default: 0.015

bfloat

Low-redshift power-law slope.

default: 2.7

cfloat

Turnover redshift parameter.

default: 2.9

dfloat

High-redshift power-law slope.

default: 5.6

Returns:
rate_densityfloat or numpy.ndarray

Merger rate density (Mpc^-3 yr^-1).

Examples

>>> from ler.gw_source_population import merger_rate_density_madau_dickinson2014
>>> rate_density = merger_rate_density_madau_dickinson2014(zs=np.array([0.1]))
ler.gw_source_population.prior_functions.merger_rate_density_madau_dickinson_belczynski_ng_function(zs, R0=19 * 1e-09, alpha_F=2.57, beta_F=5.83, c_F=3.36)[source]

Compute BBH merger rate density following Ng et al. (2021).

This model uses a Madau-Dickinson-like functional form to fit the merger rate density of field BHs, accounting for time delays and metallicity effects. Coefficients from Madau & Dickinson (2014) are translated as: B-> alpha_F, D-> beta_F, C-> c_F.

density(zs) ∝ (1 + zs) ** alpha_F / (1 + ((1 + zs) / c_F) ** beta_F)

Parameters:
zsfloat or numpy.ndarray

Source redshifts.

R0float

Local merger rate density (Mpc^-3 yr^-1).

default: 19e-9

alpha_Ffloat

Low-redshift power-law slope.

default: 2.57

beta_Ffloat

High-redshift power-law slope.

default: 5.83

c_Ffloat

Turnover redshift parameter.

default: 3.36

Returns:
rate_densityfloat or numpy.ndarray

Merger rate density (Mpc^-3 yr^-1).

Examples

>>> from ler.gw_source_population import merger_rate_density_madau_dickinson_belczynski_ng
>>> rate_density = merger_rate_density_madau_dickinson_belczynski_ng(zs=np.array([0.1]))
ler.gw_source_population.prior_functions.merger_rate_density_bbh_primordial_ken2022_function(zs, cosmology=None, n0=0.044 * 1e-09, t0=13.786885302009708)[source]

Compute the merger rate density for Primordial BBH.

Reference: Ng et al. (2022). The output is in detector frame and is unnormalized.

Parameters:
zsfloat or numpy.ndarray

Source redshifts.

cosmologyastropy.cosmology or None

Cosmology object for age calculations.

default: LambdaCDM(H0=70, Om0=0.3, Ode0=0.7, Tcmb0=0.0, Neff=3.04, m_nu=None, Ob0=0.0)

n0float

Normalization constant.

default: 0.044e-9

t0float

Present age of the Universe (Gyr).

default: 13.786885302009708

Returns:
rate_densityfloat or numpy.ndarray

Merger rate density.

Examples

>>> from ler.gw_source_population import merger_rate_density_bbh_primordial_ken2022
>>> rate_density = merger_rate_density_bbh_primordial_ken2022(zs=np.array([0.1]))
ler.gw_source_population.prior_functions.sfr_madau_fragos2017_with_bbh_td(zs, R0=19 * 1e-09)[source]

Compute the merger rate density for BBH. This is computed from star formation rate, Madau & Fragos (2017), with an additional time delay. This function is relies on pre-generated data points.

Parameters:
zsfloat or numpy.ndarray

Source redshifts.

R0float

Local merger rate density (Mpc^-3 yr^-1).

default: 19e-9

Returns:
SFRfloat or numpy.ndarray

Star formation rate (Mpc^-3 yr^-1).

ler.gw_source_population.prior_functions.sfr_madau_dickinson2014_with_bbh_td(zs, R0=19 * 1e-09)[source]

Compute the merger rate density for BBH. This is computed from star formation rate, Madau & Dickinson (2014), with an additional time delay. This function is relies on pre-generated data points.

Parameters:
zsfloat or numpy.ndarray

Source redshifts.

R0float

Local merger rate density (Mpc^-3 yr^-1).

default: 19e-9

Returns:
SFRfloat or numpy.ndarray

Star formation rate (Mpc^-3 yr^-1).

ler.gw_source_population.prior_functions.sfr_madau_fragos2017_with_bns_td(zs, R0=89 * 1e-09)[source]

Compute the merger rate density for BNS. This is computed from star formation rate, Madau & Fragos (2017), with an additional time delay. This function is relies on pre-generated data points.

Parameters:
zsfloat or numpy.ndarray

Source redshifts.

R0float

Local merger rate density (Mpc^-3 yr^-1).

default: 89e-9

Returns:
SFRfloat or numpy.ndarray

Star formation rate (Mpc^-3 yr^-1).

ler.gw_source_population.prior_functions.sfr_madau_dickinson2014_with_bns_td(zs, R0=89 * 1e-09)[source]

Compute the merger rate density for BNS. This is computed from star formation rate, Madau & Dickinson (2014), with an additional time delay. This function is relies on pre-generated data points.

Parameters:
zsfloat or numpy.ndarray

Source redshifts.

R0float

Local merger rate density (Mpc^-3 yr^-1).

default: 89e-9

Returns:
SFRfloat or numpy.ndarray

Star formation rate (Mpc^-3 yr^-1).

ler.gw_source_population.prior_functions.sfr_madau_fragos2017(zs, a=0.01, b=2.6, c=3.2, d=6.2)[source]

Compute star formation rate using Madau & Fragos (2017) model.

Reference: https://arxiv.org/pdf/1606.07887.pdf

Parameters:
zsfloat or numpy.ndarray

Source redshifts.

afloat

Normalization parameter.

default: 0.01

bfloat

Low-redshift power-law slope.

default: 2.6

cfloat

Turnover redshift parameter.

default: 3.2

dfloat

High-redshift power-law slope.

default: 6.2

Returns:
SFRfloat or numpy.ndarray

Star formation rate (Msun yr^-1 Mpc^-3).

ler.gw_source_population.prior_functions.sfr_madau_dickinson2014(zs, a=0.015, b=2.7, c=2.9, d=5.6)[source]

Compute star formation rate using Madau & Dickinson (2014) model.

Reference: Eqn. 15 of https://arxiv.org/pdf/1403.0007

Parameters:
zsfloat or numpy.ndarray

Source redshifts.

afloat

Normalization parameter.

default: 0.015

bfloat

Low-redshift power-law slope.

default: 2.7

cfloat

Turnover redshift parameter.

default: 2.9

dfloat

High-redshift power-law slope.

default: 5.6

Returns:
SFRfloat or numpy.ndarray

Star formation rate (Msun yr^-1 Mpc^-3).

Examples

>>> from ler.gw_source_population import sfr_madau_dickinson2014
>>> sfr = sfr_madau_dickinson2014(zs=np.array([0.1]))
ler.gw_source_population.prior_functions.binary_masses_BBH_popIII_lognormal_rvs(size, m_min=1.0, m_max=100.0, Mc=20.0, sigma=0.3, chunk_size=10000)[source]

Sample from a lognormal distribution in 2D mass space.

Reference: Ng et al. (2022). This is a helper function for PopIII BBH and primordial BBH merger rate density distribution functions.

Parameters:
sizeint

Number of samples to draw.

m_minfloat

Minimum mass (Msun).

default: 1.0

m_maxfloat

Maximum mass (Msun).

default: 100.0

Mcfloat

Characteristic mass scale (Msun).

default: 20.0

sigmafloat

Width of the distribution.

default: 0.3

chunk_sizeint

Number of samples per rejection sampling chunk.

default: 10000

Returns:
m1_samplenumpy.ndarray

Primary mass samples (Msun).

m2_samplenumpy.ndarray

Secondary mass samples (Msun).

Examples

>>> from ler.gw_source_population import binary_masses_BBH_popIII_lognormal
>>> m1, m2 = binary_masses_BBH_popIII_lognormal(size=1000)
ler.gw_source_population.prior_functions.binary_masses_BBH_primordial_lognormal_rvs(size, m_min=1.0, m_max=100.0, Mc=20.0, sigma=0.3, chunk_size=10000)[source]

Sample from a lognormal distribution in 2D mass space.

Based on Eqn. 1 and 4 of Ng et al. 2022 for primordial black holes.

Parameters:
sizeint

Number of samples to draw.

m_minfloat

Minimum mass (Msun).

default: 1.0

m_maxfloat

Maximum mass (Msun).

default: 100.0

Mcfloat

Characteristic mass scale (Msun).

default: 20.0

sigmafloat

Width of the distribution.

default: 0.3

chunk_sizeint

Number of samples per rejection sampling chunk.

default: 10000

Returns:
m1_samplenumpy.ndarray

Primary mass samples (Msun).

m2_samplenumpy.ndarray

Secondary mass samples (Msun).

Examples

>>> from ler.gw_source_population import binary_masses_BBH_primordial_lognormal
>>> m1, m2 = binary_masses_BBH_primordial_lognormal(size=1000)
ler.gw_source_population.prior_functions.binary_masses_BNS_bimodal_rvs(size, w=0.643, muL=1.352, sigmaL=0.08, muR=1.88, sigmaR=0.3, mmin=1.0, mmax=2.3, resolution=500)[source]

Sample BNS masses from bimodal Gaussian distribution.

Based on Will M. Farr et al. 2020 Eqn. 6 for neutron star mass distribution combining two Gaussian peaks.

Parameters:
sizeint

Number of samples to draw.

wfloat

Weight of the left (low-mass) peak.

default: 0.643

muLfloat

Mean of the left peak (Msun).

default: 1.352

sigmaLfloat

Standard deviation of the left peak (Msun).

default: 0.08

muRfloat

Mean of the right peak (Msun).

default: 1.88

sigmaRfloat

Standard deviation of the right peak (Msun).

default: 0.3

mminfloat

Minimum mass (Msun).

default: 1.0

mmaxfloat

Maximum mass (Msun).

default: 2.3

resolutionint

Number of points to use for the CDF.

default: 500

Returns:
m1numpy.ndarray

Primary mass samples (Msun).

m2numpy.ndarray

Secondary mass samples (Msun).

ler.gw_source_population.prior_functions.binary_masses_NSBH_broken_powerlaw_rvs(size=1000, mminbh=26.0, mmaxbh=125.0, alpha_1=6.75, alpha_2=0.0, b=0.5, delta_m=5.0, mminns=1.0, mmaxns=3.0, alphans=0.0, normalization_size=1000)[source]

Generate NSBH mass samples from broken power-law (BH) and power-law (NS).

Parameters:
sizeint

Number of samples to draw.

default: 1000

mminbhfloat

Minimum BH mass (Msun).

default: 26.0

mmaxbhfloat

Maximum BH mass (Msun).

default: 125.0

alpha_1float

BH power-law index below break.

default: 6.75

alpha_2float

BH power-law index above break.

default: 0.0

bfloat

Break location parameter (0-1).

default: 0.5

delta_mfloat

Smoothing width (Msun).

default: 5.0

mminnsfloat

Minimum NS mass (Msun).

default: 1.0

mmaxnsfloat

Maximum NS mass (Msun).

default: 3.0

alphansfloat

NS power-law index.

default: 0.0

normalization_sizeint

Grid size for CDF computation.

default: 1000

Returns:
m1_samplesnumpy.ndarray

BH mass samples (Msun).

m2_samplesnumpy.ndarray

NS mass samples (Msun).

ler.gw_source_population.prior_functions.binary_masses_BBH_powerlaw_gaussian_rvs(size, mminbh, mmaxbh, alpha, mu_g, sigma_g, lambda_peak, delta_m, beta, normalization_size=1000)[source]

Generate BBH mass samples from power-law + Gaussian model with mass ratio.

Parameters:
sizeint

Number of samples to draw.

mminbhfloat

Minimum BH mass (Msun).

mmaxbhfloat

Maximum BH mass (Msun).

alphafloat

Power-law spectral index for m1.

mu_gfloat

Mean of the Gaussian peak (Msun).

sigma_gfloat

Standard deviation of the Gaussian peak (Msun).

lambda_peakfloat

Fraction in Gaussian component (0-1).

delta_mfloat

Low-mass smoothing width (Msun).

betafloat

Power-law index for mass ratio distribution.

normalization_sizeint

Grid size for CDF computation.

default: 1000

Returns:
m1numpy.ndarray

Primary mass samples (Msun).

m2numpy.ndarray

Secondary mass samples (Msun).

ler.gw_source_population.prior_functions.available_prior_list()[source]

Returns a list of available priors.