galsbi.ucat_sps.galaxy_population_models package

Submodules

galsbi.ucat_sps.galaxy_population_models.galaxy_agn module

galsbi.ucat_sps.galaxy_population_models.galaxy_agn.sample_ratio_agn_to_galaxy_bolometric_luminosity(par, redshifts, log10_stellar_masses, galaxy_type='blue', seed=None)[source]

We first compute the probability that a galaxy hosts an AGN as function of redshift and stellar mass by taking the ratio between the redshift- evolving SMFs of blue and red galaxies and the AGN host galaxy mass function dervied in Bongiorno et al 2016. Then, we statistically assign an AGN to a host galaxy with a Bernoulli trial proportional to the computed probability. As in López-López et al. (2024), we assume fixed AGN fractions outside the redshift ranges specified by Bongiorno et al. (2016). In particular, the probability of hosting an AGN is 1% at all stellar masses below z = 0.55, while at z > 2 we maintain a constant fraction taken from the last redshift bin in Bongiorno et al. (2016). We then assign AGN flux contributions to all galaxies, regardless of whether they are classified as AGN hosts. However, fAGN is sampled from different distributions depending on the AGN classification. Galaxies hosting AGNs draw from the high-log (fAGN) distribution, while non-AGN galaxies draw from the low-log (fAGN) one. We adopt a threshold of log(fAGN) = −1 to separate these two regimes. The ratio of the AGN to galaxy bolometric luminosity fagn is sampled from the sum of two Normal distributions that matches the fagn distribution computed with Prospect using the best-fitting parameters of DEVILS D10 field used in Thorne+2022. fagn = L_bol_AGN / L_bol_galaxy.

Parameters:
  • par – (obj) par objects containing the Ucat parameters.

  • redshifts – (array_like[n_gal,]) redshifts drawn from the stellar mass function.

  • log10_stellar_masses – (array_like[n_gal,]) log10 of the stellar masses sampled from the stellar mass function.

  • galaxy_type – (str) “blue” or “red” galaxy population.

  • seed – (int) random seed.

Return fagn:

(array_like[n_gal,]) AGN to galaxy bolometric luminosity fraction.

galsbi.ucat_sps.galaxy_population_models.galaxy_dust_attenuation module

galsbi.ucat_sps.galaxy_population_models.galaxy_dust_attenuation.sample_dust_attenuation_parameters(par, log10_stellar_mass, log_specific_star_formation_rate, redshift, galaxy_type='blue', seed=None)[source]

This function samples dust attenuation parameters depending on the specific user- defined dust attenuation law.

Parameters:
  • par – (obj) par objects containing the Ucat parameters.

  • log10_stellar_mass

  • log_specific_star_formation_rate – (array_like[n_gal,]) log10 of the specific star formation rate.

  • redshift – (array_like[n_gal,]) galaxy redshifts.

  • galaxy_type – (str) “blue” or “red” galaxy population.

  • seed – (int) random seed.

Return dust_attenuation_params:

(array_like[n_gal,n_params]) dust attenuation parameters for specific attenuation law.

galsbi.ucat_sps.galaxy_population_models.galaxy_dust_attenuation.sample_dust_attenuation_parameters_charlotfall2000(par, log10_stellar_masses, log_specific_star_formation_rates, redshifts, galaxy_type='blue', seed=None)[source]

The Charlot&Fall2000 dust attenuation prescription is used in this case. It is implemented in ProSpect by means of the tau_birth, tau_screen and Eb parameters. We use DEVILS data to model the dependence of tau_screen and tau_birth on redshift, stellar mass and specific star-formation rate. We find that tau_birth does not depend on any of those quantities. tau_screen depends instead on redshift, stellar mass and specific star-formation rate. tau_birth is sampled from a Gaussian distribution with different parameter values for blue and red galaxies. The Eb of the 2175 Å bump in the Drude profile is assumed to have the typical value of 1.

Parameters:
  • par – (obj) par objects containing the Ucat parameters.

  • log10_stellar_masses – (array_like[n_gal,]) log10 of the stellar masses drawn from the stellar mass function.

  • log_specific_star_formation_rates – (array_like[n_gal,]) log10 of the specific star formation rate.

  • redshifts – (array_like[n_gal,]) galaxy redshifts.

  • galaxy_type – (str) “blue” or “red” galaxy population.

  • seed – (int) random seed.

Return dust_attenuation_params:

(array_like[n_gal,3]) dust attenuation parameters tau_birth, tau_screen, Eb.

galsbi.ucat_sps.galaxy_population_models.galaxy_dust_emission module

galsbi.ucat_sps.galaxy_population_models.galaxy_dust_emission.sample_Dale2014_dust_emission_parameters(par, n_galaxies, galaxy_type='blue', seed=None)[source]

Sample Dale 2014 dust emission parameters from modelling the distribution in Thorne+2022. These parameters are used by Prospect to simulate dust emission which has a contribution of the order of 10^-4 in griZ and 10^-3 in JHKs. DaleNormSFR only for Dale2014 templates including star formation and not agn component which is instead given by Temple templates.

Parameters:
  • par – (obj) par objects containing the Ucat parameters.

  • n_galaxies – (int) number of galaxies to sample.

  • galaxy_type – (str) “blue” or “red” galaxy population.

  • seed – (int) random seed.

Return dust_emission_params:

(array_like[n_gal,2]) dust emission parameters for Dale 2014 implementation in ProSpect.

galsbi.ucat_sps.galaxy_population_models.galaxy_gas_ionization module

galsbi.ucat_sps.galaxy_population_models.galaxy_gas_ionization.sample_gas_ionization_Kashino19(Zgas_final, star_formation_rates, log10_stellar_masses, par, seed=None)[source]

This function computes the gas ionization parameter using the relation in equation 12 of Kashino+19. We keep the electron density fixed to the value with which the CLOUDY templates of FSPS have been generated, i.e. ne~100 cm^-3. We sample logU from a Normal distribution with mean given by the Kashino+19 relation and scatter of 0.1 dex, as highlighted in Kashino+19. The Normal distribution is truncated within the range [-4, -1].

Parameters:
  • Zgas_final – (array_like[n_gal,]) final gas-phase metallicity in absolute units.

  • star_formation_rates – (array_like[n_gal,]) SFR averaged over the last 100 Myr of the galaxy SFH. Units are Msun/yr.

  • log10_stellar_masses – (array_like[n_gal,]) log10 of the stellar masses sampled from the stellar mass function.

  • par – (obj) par objects containing the Ucat parameters.

  • seed – (int) random seed for the truncnorm sampling.

Return logU:

(array_like[n_gal,]) gas ionization parameters.

galsbi.ucat_sps.galaxy_population_models.light_profile module

galsbi.ucat_sps.galaxy_population_models.galaxy_light_profile.sample_sersic_for_galaxy_type(n_gal, galaxy_type, par, log10_stellar_mass)[source]

This function samples Sérsic indices conditioned on the galaxy stellar masses for different galaxy types.

Parameters:
  • n_gal – (int) number of galaxies.

  • galaxy_type – (str) “blue” or “red” galaxy population.

  • par – (obj) par objects containing the Ucat parameters.

  • log10_stellar_mass – (array_like[n_gal,]) log10 of the stellar masses sampled from the stellar mass function.

Return sersic_n:

(array_like[n_gal,]) Sérsic indices.

galsbi.ucat_sps.galaxy_population_models.galaxy_light_profile.sample_sersic_stellar_mass(n_gal, n_0, n_1, n_2, n_scatter, log10_stellar_mass, min_n, max_n)[source]

This function samples Sérsic indices conditioned on the galaxy stellar masses. The stellar mass dependence of the Sérsic index is modelled using the GAMA data in Kelvin+12.

Parameters:
  • n_gal – (int) number of galaxies.

  • n_0 – (float) Sérsic index - stellar mass relation intercept.

  • n_1 – (float) Sérsic index - stellar mass relation slope.

  • n_2 – (float) Sérsic index - stellar mass relation exponent.

  • n_scatter – (float) Sérsic index - stellar mass relation scatter.

  • log10_stellar_mass – (array_like[n_gal,]) log10 of the stellar masses sampled from the stellar mass function.

  • min_n – (float) minimum Sérsic index for the truncated distribution.

  • max_n – (float) maximum Sérsic index for the truncated distribution.

Return sersic_n:

(array_like[n_gal,]) Sérsic indices.

galsbi.ucat_sps.galaxy_population_models.galaxy_metallicity_history module

galsbi.ucat_sps.galaxy_population_models.galaxy_metallicity_history.generate_gas_metallicity_history_snorm_trunc(Zgas_final, star_formation_histories, par, Zgas_init=0.0001)[source]

This function creates, for each galaxy, a gas metallicity history under the assumption that the metallicity enrichment follows the mass build up. This prescription assumes the same low initial gas-phase metallicity for every galaxy, while different final Zgas, drawn from the sample_gas_metallicity_function. In order to compute the stellar mass formed in each time step, we use the star formation histories as function of lookback time. par.time_step needs to be equal to the time step in par.time_grid.

Parameters:
  • Zgas_final – (array_like[n_gal,]) final gas-phase metallicity in absolute units.

  • star_formation_histories – (array_like[n_gal,n_time_bins]) galaxy SFHs along a grid in lookback time. Units are Msun/yr.

  • par – (obj) par objects containing the Ucat parameters.

  • Zgas_init – (float) The initial value of the gas-phase metallicity history. Best to set this value to the lowest metallicity available in the stellar/gas templates.

Return Z_gas:

(array_like[n_gal,n_time_bins]) galaxy gas-phase metallicity histories along a grid in time since the beginning of the Universe.

galsbi.ucat_sps.galaxy_population_models.galaxy_metallicity_history.integrate_sfh(sfh, time_grid, time_step)[source]

This function computes the cumulative formed stellar mass based on the input star formation history.

Parameters:
  • sfh – (array_like[n_time_bins, ]) galaxy SFHs along a grid in increasing time. Units are Msun/yr.

  • time_grid – (array_like[n_time_bins, ]) time grid in units of yr.

  • time_step – (float) time grid step in units of yr.

Return cumulative_mass:

(array_like[n_time_bins, ]) formed stellar mass in units of Msun.

galsbi.ucat_sps.galaxy_population_models.galaxy_metallicity_history.sample_gas_metallicity(log10_stellar_masses, redshifts, star_formation_rates, par, cosmology, galaxy_type='blue', seed=None)[source]

This function samples the final gas metallicities of the galaxies based on equation 10 in Bellstedt+21. The coefficients of the equation evolves quadratically with lookback time of the observer (not that of the galaxy). The values in Bellstedt+21 have been obtained only for a sample of star-forming galaxies. The scatter around the plane also evolves quadratically with lookback time. The lookback times should be in units of Gyr.

The final gas metallicity in absolute units is obtained following equation 1 in Bellstedt+21, $Z_{gas} = Z_{odot} imes 10^{[12+log{(O/H)}] - [12+log{(O/H)}]_{odot}}$ It represents the current gas-phase metallicity of the galaxy.

Parameters:
  • log10_stellar_masses – (array_like[n_gal,]) log10 of the stellar masses drawn from the stellar mass function.

  • redshifts – (array_like[n_gal,]) redshifts drawn from the stellar mass function.

  • star_formation_rates – (array_like[n_gal,]) SFR averaged over the last 100 Myr of the galaxy SFH. Units are Msun/yr.

  • par – (obj) par objects containing the Ucat parameters.

  • cosmology – (obj) PyCosmo class object.

  • galaxy_type – (str) “blue” or “red” galaxy population.

  • seed – (int) random seed.

Return Zgas_final:

(array_like[n_gal,]) final gas-phase metallicity in absolute units.

galsbi.ucat_sps.galaxy_population_models.galaxy_metallicity_history.vectorized_truncnorm_sample(mean, std, lower_bound, upper_bound, size=None, seed=None)[source]

galsbi.ucat_sps.galaxy_population_models.galaxy_size module

galsbi.ucat_sps.galaxy_population_models.galaxy_size.sample_r50_from_stellar_mass(redshifts, log10_stellar_mass, cosmology, par, galaxy_type='blue', seed=None)[source]

This function samples galaxy sizes conditioned on their stellar masses and types.

Parameters:
  • redshifts – (array_like[n_gal,]) redshifts drawn from the stellar mass function.

  • log10_stellar_mass – (array_like[n_gal,]) log10 of the stellar masses drawn from the stellar mass function.

  • cosmology – (obj) PyCosmo class object.

  • par – (obj) par objects containing the Ucat parameters.

  • galaxy_type – (str) “blue” or “red” galaxy population.

  • seed – (int) random seed.

Return r50_ang:

(array_like[n_gal,]) angular sizes in pixel.

Return r50_phys:

(array_like[n_gal,]) physical sizes of galaxies in kpc.

galsbi.ucat_sps.galaxy_population_models.galaxy_size.sample_r50_phys_from_stellar_mass(redshifts, log10_stellar_mass, par, galaxy_type, seed=None)[source]

This function samples the physical galaxy sizes conditioned on the stellar mass. The mean size has a power-law dependence from stellar mass. A single power-law correctly represents blue galaxy sizes, while a flattend power-law at low sizes is required for red galaxies. The parameters of the power-laws evolve linearly with redshift. The actual galaxy sizes are drawn from a log-normal distribution centred on the mean size from the power-law and with scatter value taken from Shen+2003. The functional forms are from Nedkova+2021.

Parameters:
  • redshifts – (array_like[n_gal,]) redshifts drawn from the stellar mass function.

  • log10_stellar_mass – (array_like[n_gal,]) log10 of the stellar masses drawn from the stellar mass function.

  • par – (obj) par objects containing the Ucat parameters.

  • galaxy_type – (str) “blue” or “red” galaxy population.

  • seed – (int) random seed.

Return r50_phys:

(array_like[n_gal,]) physical sizes of galaxies in kpc.

galsbi.ucat_sps.galaxy_population_models.galaxy_star_formation_history module

galsbi.ucat_sps.galaxy_population_models.galaxy_star_formation_history.compute_surviving_stellar_mass(redshifts, log10_stellar_masses, mSFR, mpeak, mperiod, mskew, Zgas_final, par)[source]

This function uses an emulator model to compute the surviving stellar mass of individual galaxies conditioned on their star formation histories and metallicities. Galaxies are split in bins of 0.5 in redshift and then passed to the emulator for log10 surviving stellar mass predictions.

Parameters:
  • redshifts – (array_like[n_gal, ]) redshifts drawn from the stellar mass function.

  • log10_stellar_masses – (array_like[n_gal, ]) log10 of the stellar masses drawn from the stellar mass function.

  • mSFR – (array_like[n_gal, ]) The values of the SFR at the peak of the SFH in

units of Msun/yr. :param mpeak: (array_like[n_gal, ]) The times at which the SFH peaks, i.e. the

time at which the SFR is equal to mSFR, in units of Gyr.

Parameters:
  • mperiod – (array_like[n_gal, ]) The widths of the normal distribution in units of Gyr.

  • mskew – (array_like[n_gal, ]) The skewnesses of the normal distribution in units of Gyr.

  • Zgas_final – (array_like[n_gal, ]) final gas-phase metallicity in absolute units.

  • par – (obj) par objects containing the Ucat parameters.

Return log10_surviv_stellar_masses:

(array_like[n_gal,]) log10 of surviving stellar masses

galsbi.ucat_sps.galaxy_population_models.galaxy_star_formation_history.get_age_of_the_universe_from_redshift(redshifts, cosmology)[source]

This function returns the age of the Universe at redshift z in years.

Parameters:
  • redshifts – (array_like[n_gal,]) galaxy redshifts.

  • cosmology – (obj) PyCosmo cosmology object.

Return t:

(array_like[n_time_bins,]) age of the Universe at the galaxy redshifts in years.

galsbi.ucat_sps.galaxy_population_models.galaxy_star_formation_history.get_lookbacktime_from_redshift(redshifts, cosmology)[source]

This function returns the lookback time of the observer at redshift z in years.

Parameters:
  • redshifts – (array_like[n_time_bins,]) galaxy redshifts.

  • cosmology – (obj) PyCosmo cosmology object.

Return lookbacktime:

(array_like[n_time_bins,]) observer lookback time at the galaxy redshifts in years.

galsbi.ucat_sps.galaxy_population_models.galaxy_star_formation_history.get_mean_sfh_shape_params(log10_stellar_masses, redshifts, par, galaxy_type='blue')[source]

This function returns the mean values of the SFH massfunc_snorm_trunc shape parameters conditioned on galaxy redshifts and stellar masses. We assume that mskew and logmperiod evolve linearly with stellar mass, while mpeak with both redshift and stellar mass. The stellar mass and redshift evolutions are assumed to be linear in log10(stellar_mass) and log10(1+z).

Parameters:
  • log10_stellar_masses – (array_like[n_gal,]) log10 of the stellar masses drawn from the stellar mass function.

  • redshifts – (array_like[n_gal,]) redshifts drawn from the stellar mass function.

  • par – (obj) par objects containing the Ucat parameters.

  • galaxy_type – either blue or red.

Return mean_sfh_shape_params:

(array_like[n_gal, 3]) mean values of the SFH shape parameters obtained as function of redshift and stellar mass.

galsbi.ucat_sps.galaxy_population_models.galaxy_star_formation_history.massfunc_snorm_trunc(age, mSFR=10.0, mpeak=10.0, mperiod=1.0, mskew=0.5, magemax=13.8, mtrunc=2.0)[source]

This functions is a Python translation of the R code in ProSpect that generates SFHs from the massfunc_snorm_trunc parametrisation described in Robotham+20.

Parameters:
  • age – (array_like[n_time_bins,]) Lookback time in units of years. 0 means now from the point of view of the galaxy.

  • mSFR – (float) The value of the SFR at the peak of the SFH in Msun/yr.

  • mpeak – (float) The time at which the SFH peaks, i.e. the time at which the SFR is equal to mSFR, in units of Gyr.

  • mperiod – (float) The width of the normal distribution in units of Gyr.

  • mskew – (float) The skewness of the normal distribution in units of Gyr.

  • magemax – (float) The age beyond which the SFR is set to 0, it is defined as 13.4 Gyr - lookbacktime(z), where magemax_z0=13.4 Gyr represents the epoch at which the highest redshift galaxies are known to exist.

  • mtrunc – (float) This parameter controls truncation between mpeak and magemax. 0 means no extra truncation, 1 is about half truncated and 2 is softly but fairly fully truncated.

Return sfh:

(array_like[n_time_bins,]) The SFH of the galaxy in units of Msun/yr. The array as the same length as age.

galsbi.ucat_sps.galaxy_population_models.galaxy_star_formation_history.sample_sfh_shape_params(log10_stellar_masses, redshifts, par, galaxy_type='blue', seed=None)[source]

This function samples the SFH shape parameters from a Multivariate Truncated Normal distribution.

Parameters:
  • log10_stellar_masses – (array_like[n_gal, ]) log10 of the stellar masses drawn from the stellar mass function.

  • redshifts – (array_like[n_gal, ]) redshifts drawn from the stellar mass function.

  • par – (obj) par objects containing the Ucat parameters.

  • galaxy_type – either blue or red.

  • seed – (int) random seed for the TruncatedMVN sampling.

Return mpeak_norm:

(array_like[n_gal, ]) fraction of the galaxy age at which SFH peaks.

Return mperiod:

(array_like[n_gal, ]) The width of the normal distribution in units of Gyr.

Return mskew:

(array_like[n_gal, ]) The skewness of the normal distribution in units of Gyr.

galsbi.ucat_sps.galaxy_population_models.galaxy_star_formation_history.sample_sfh_snorm_trunc(log10_stellar_masses, redshifts, par, cosmology, galaxy_type='blue', seed=None)[source]

This function computes the galaxy star formation histories using the parametric star formation history massfunc_snorm_trunc implemented in Robotham+20 and used in Bellstedt+20,21. This parametrization is based on the four parameters mSFR, mpeak, mperiod, mskew. We do not model their distribution from GAMA and DEVILS data directly, otherwise we would end up imposing their redshift and stellar mass functions. Instead, we model the stellar mass and redshift dependence of the shape parameters mpeak, mperiod and mskew, since observationally we expect that more massive galaxies assembled their mass earlier and more rapidly than less massive galaxies. mpeak has been redefined to be the fractional age at which the SFH peaks, i.e., when mSFR occurs, where 0 means at the beginning of the galaxy life, and 1 means at the galaxy redshift. mpeak depends on both redshift and stellar mass, while mperiod and mskew depend only on stellar mass as we did not find a well-defined trend in data. We model the dependency on these physical quantities as linear relations. Once sampled from a truncated multivariate Normal distribution, we compute the shape of the SFH by imposing mSFR=1.0. We then compute the correct value of mSFR such that the integration of the star formation history yields the formed mass sampled from the stellar mass function.

The time grid used here refers to the lookback time with respect to the galaxy, so t=0 effectively means now for the galaxy, not the observer. The grid goes, however, up to 13.8 Gyr to have a common time grid for every galaxy. The SFH is set to be 0 beyond magemax = 13.4 - lookbacktime(z), i.e., age of the Universe at the galaxy redshift. We limit the sampled parameters to the observational limits set in Bellstedt+21. We also compute the star formation rate averaged over the last 100 Myr of the galaxy. magemax_z0=13.4 Gyr is the age at which the earliest galaxies are known to exist. mtrunc=2 ensures smooth truncation at SFR=0 at the beginning of the galaxy life.

Parameters:
  • log10_stellar_masses – (array_like[n_gal, ]) log10 of the stellar masses drawn from the stellar mass function.

  • redshifts – (array_like[n_gal, ]) redshifts drawn from the stellar mass function.

  • par – (obj) par objects containing the Ucat parameters.

  • cosmology – (obj) PyCosmo class object.

  • galaxy_type – Either blue or red.

  • seed – (int) random seed for the TruncatedMVN sampling.

Return sfhs:

(array_like[n_gal, n_time_bins]) galaxy star formation histories in units of Msun/yr on the linear lookback time grid.

Return sfrs:

(array_like[n_gal, ]) star formation rates averaged over the last 100 Myr of the galaxy star formation history in units of Msun/yr.

Return mSFR:

(float) The value of the SFR at the peak of the SFH in units of Msun/yr after multiplying by formed stellar mass.

Return mpeak_reconstrs:

(float) The time at which the SFH peaks, i.e. the time at which the SFR is equal to mSFR, in units of Gyr.

Return mperiod:

(float) The width of the normal distribution in units of Gyr.

Return mskew:

(float) The skewness of the normal distribution in units of Gyr.

galsbi.ucat_sps.galaxy_population_models.galaxy_star_formation_history.single_galaxy_sfh_snorm_trunc(time_grid, mSFR_norm, mpeak, mperiod, mskew, magemax, mtrunc)[source]

Evaluate the truncated skewed normal SFH for a single set of parameters.

This is a light-weight convenience wrapper around massfunc_snorm_trunc() that makes the intent explicit in the public API and keeps scalar inputs behaving like scalars (which is relied upon in the tests).

galsbi.ucat_sps.galaxy_population_models.galaxy_stellar_mass_function module

class galsbi.ucat_sps.galaxy_population_models.galaxy_stellar_mass_function.NumGalCalculator(z_max, sm_min, parametrization, sm_alpha, sm_m_star_par, sm_phi_star_par, cosmo, pixarea, ngal_multiplier=1)[source]

Bases: object

Computes galaxy number counts by integrating the galaxy stellar mass function (Schechter function).

The integral over stellar masses can be done analytically, while the integral over redshifts is computed numerically. See also docs/jupyter_notebooks/sample_redshift_magnitude.ipynb.

class galsbi.ucat_sps.galaxy_population_models.galaxy_stellar_mass_function.RedshiftSMSampler(z_res, z_max, sm_res, sm_min, parametrization, sm_alpha, sm_m_star_par, sm_phi_star_par, cosmo)[source]

Bases: object

Samples redshifts and stellar masses from the galaxy stellar mass function.

The sampling is done by first drawing redshifts from the redshift-pdf obtained by integrating out stellar masses. Then, we sample stellar masses from the conditional pdfs obtained by conditioning the stellar mass function on the sampled redshifts (the conditional pdf is different for each redshift). See also docs/jupyter_notebooks/sample_redshift_magnitude.ipynb and docs/jupyter_notebooks/test_self_consistency.ipynb.

class galsbi.ucat_sps.galaxy_population_models.galaxy_stellar_mass_function.StellarMassFunction(name, sm_fct_parametrization, sm_m_star_0, sm_m_star_1, sm_m_star_2, sm_phi_star_amp, sm_phi_star_exp, sm_alpha, cosmo, pixarea, galaxy_type, seed_ngal, z_res=0.001, sm_res=0.01, z_max=3, sm_min=8, z_sm_intp=None, ngal_multiplier=1)[source]

Bases: object

Stellar Mass function.

sample_z_sm_and_apply_cut(seed_ngal, n_gal_max=inf, size_chunk=10000)[source]

This function gets the abs mag and z using chunking, which uses less memory than the original method.

It does not give exactly the same result as before due to different order of random draws in z_mabs_sampler, but it’s the same sample.

galsbi.ucat_sps.galaxy_population_models.galaxy_stellar_mass_function.double_stellar_mass_function_parametrisation(log10_stellar_masses, log_Mstar, phi_star_low_mass, phi_star_high_mass, alpha_low_mass, alpha_high_mass)[source]

This function parametrises the stellar mass function as a double Schechter function. It computes the number density of galaxies as a function of their stellar mass.

Parameters:
  • log10_stellar_masses – (array_like[n_mass_bins,]) log10 of the stellar masses drawn from the SMF.

  • log_Mstar – (float) knee of the Schechter function. It corresponds to the stellar mass at which the Schechter transitions from a simple power law with slope alpha at lower masses into an exponential function at higher masses. It has the same value for both the low and high mass SMFs.

  • phi_star_low_mass – (float) number density of galaxies at stellar mass M* (log_Mstar) for the low mass SMF.

  • phi_star_high_mass – (float) number density of galaxies at stellar mass M* (log_Mstar) for the high mass SMF.

  • alpha_low_mass – (float) faint-end slope of the low mass SMF.

  • alpha_high_mass – (float) faint-end slope of the high mass SMF.

Return stellar_mass_function:

(array_like[n_mass_bins,]) the number density of galaxies as a function of their stellar mass.

galsbi.ucat_sps.galaxy_population_models.galaxy_stellar_mass_function.find_closest_ind(grid, vals)[source]
galsbi.ucat_sps.galaxy_population_models.galaxy_stellar_mass_function.initialize_stellar_mass_functions(par, pixarea, cosmo, z_sm_intp=None, mass_key=None)[source]
galsbi.ucat_sps.galaxy_population_models.galaxy_stellar_mass_function.maximum_redshift(z_sm_intp, sm_min, z_max, parametrization, alpha, sm_m_star_par, seed_ngal)[source]

Computes the maximum redshift up to which we sample objects from the luminosity function.

The cutoff is based on the criterion that the CDF for absolute magnitudes is larger than 1e-5, i.e. that there is a reasonable probability of actually obtaining objects at this redshift and absolute magnitude which still pass the cut on par.gals_mag_max.

galsbi.ucat_sps.galaxy_population_models.galaxy_stellar_mass_function.single_stellar_mass_function_parametrisation(log10_stellar_masses, log_Mstar, phi_star, alpha)[source]

This function parametrises the stellar mass function as a single Schechter function. It computes the number density of galaxies as a function of their stellar mass.

Parameters:
  • log10_stellar_masses – (array_like[n_mass_bins,]) log10 of the stellar masses drawn from the SMF.

  • log_Mstar – (float) knee of the Schechter function. It corresponds to the stellar mass at which the Schechter transitions from a simple power law with slope alpha at lower masses into an exponential function at higher masses.

  • phi_star – (float) number density of galaxies at stellar mass M* (log_Mstar).

  • alpha – (float) faint-end slope.

Return stellar_mass_function:

(array_like[n_mass_bins,]) the number density of galaxies as a function of their stellar mass.

galsbi.ucat_sps.galaxy_population_models.galaxy_stellar_mass_function.stellar_mass_function_m_star(z, parametrization, intercept, slope, quadratic_term=0)[source]
galsbi.ucat_sps.galaxy_population_models.galaxy_stellar_mass_function.stellar_mass_function_phi_star(z, parametrization, amplitude, exp)[source]
galsbi.ucat_sps.galaxy_population_models.galaxy_stellar_mass_function.upper_inc_gamma(a, x)[source]

galsbi.ucat_sps.galaxy_population_models.galaxy_velocity_dispersion module

galsbi.ucat_sps.galaxy_population_models.galaxy_velocity_dispersion.sample_velocity_dispersion_zahid(par, log10_stellar_masses, seed=None)[source]

This function samples the velocity dispersion of each galaxy conditioned on its stellar mass. We use equation 5 in Zahid+2016 to determine the mean velocity dispersion as function of stellar mass. Then we sample the actual velocity dispersion from a Gaussian having width from Figure 9 in Zahid+2016. We also set a lower limit equal to 10 km/s for the velocity dispersion.

Parameters:
  • par – (obj) par objects containing the Ucat parameters.

  • log10_stellar_masses – (array_like[n_gal,]) log10 of the stellar masses drawn from the stellar mass function.

  • seed – (int) random seed.

Return vel_disp:

(array_like[n_gal,]) velocity dispersions in units of km/s.

Module contents