holodeck.hardening

Binary evolution hardening submodules.

To-Do (Hardening)

  • Dynamical_Friction_NFW
    • Allow stellar-density profiles to also be specified (instead of using a hard-coded Dehnen profile)

    • Generalize calculation of stellar characteristic radius. Make self-consistent with stellar-profile, and user-specifiable.

  • Evolution
    • _sample_universe() : sample in comoving-volume instead of redshift

  • Sesana_Scattering
    • Allow stellar-density profile (or otherwise the binding-radius) to be user-specified and flexible. Currently hard-coded to Dehnen profile estimate.

  • _SHM06
    • Interpolants of hardening parameters return 1D arrays.

  • Fixed_Time_2PL
    • Handle rchar better with respect to interpolation. Currently not an interpolation variable, which restricts it’s usage.

    • This class should be separated into a generic _Fixed_Time class that can use any functional form, and then a 2-power-law functional form that requires a specified normalization. When they’re combined, it will produce the same effect. Another good functional form to implement would be GW + log-uniform hardening time, the same as the current phenomenological model but with both power-laws set to 0.

References

class holodeck.hardening.CBD_Torques(f_edd=0.1, subpc=True)[source]

Binary Orbital Evolution based on Hydrodynamic Simulations by Siwek+23.

This module uses data from Siwek+23, which supplies rates of change of binary semi-major axis a_b and binary eccentricity e_b. The calculation of a_b and e_b versus time requires accretion rates (for scale).

dadt_dedt(evo, step)[source]

Circumbinary Disk Torque hardening rate.

Parameters:
  • evo (Evolution) – Evolution instance providing binary parameters at the given intergration step.

  • step (int) – Integration step at which to calculate hardening rates.

Returns:

  • dadt (array_like) – Binary hardening rates in units of [cm/s], defined to be negative.

  • dedt (array_like) – Binary rate-of-change of eccentricity in units of [1/sec].

class holodeck.hardening.Dynamical_Friction_NFW(mmbulge=None, msigma=None, smhm=None, coulomb=10.0, attenuate=True, rbound_from_density=True)[source]

Dynamical Friction (DF) hardening module assuming an NFW dark-matter density profile.

This class calculates densities and orbital velocities based on a NFW profile with parameters based on those of each MBH binary. The holodeck.observations.NFW class is used for profile calculations, and the halo parameters are calculated from Stellar-mass–Halo-mass relations (see ‘arguments’ below). The ‘effective-mass’ of the inspiralling secondary is modeled as a power-law decreasing from the sum of secondary MBH and its stellar-bulge (calculated using the mmbulge - Mbh-Mbulge relation), down to just the bare secondary MBH after 10 dynamical times. This is to model tidal-stripping of the secondary host galaxy.

Attenuation of the DF hardening rate is typically also included, to account for the inefficiency of DF once the binary enters the hardened regime. This is calculated using the prescription from [BBR1980]. The different characteristic radii, needed for the attenuation calculation, currently use a fixed Dehnen stellar-density profile as in [Chen2017], and a fixed scaling relationship to find the characteristic stellar-radius.

Notes

  • This module does not evolve eccentricity.

  • The hardening rate (da/dt) is not allowed to be larger than the orbital/virial velocity of the halo (as a function of radius).

dadt_dedt(evo, step, attenuate=None)[source]

Calculate DF hardening rate given Evolution instance, and an integration step.

Parameters:
  • evo (Evolution instance) – The evolutionary tracks of the binary population, providing binary parameters.

  • step (int,) – Integration step at which to calculate hardening rates.

Returns:

  • dadt ((N,) np.ndarray of scalar,) – Binary hardening rates in units of [cm/s].

  • dedt ((N,) np.ndarray or None) – Rate-of-change of eccentricity, which is not included in this calculation, it is zero. None is returned if the input eccen is None.

class holodeck.hardening.Fixed_Time_2PL(time, mtot, mrat, redz, sepa_init, rchar=3.0856775814913674e+20, gamma_inner=-1.0, gamma_outer=1.5, progress=False, interpolate_norm=False)[source]

Provide a binary hardening rate such that the total lifetime matches a given value.

This class uses a phenomenological functional form (defined in Fixed_Time_2PL.function()) to model the hardening rate (\(da/dt\)) of binaries. The functional form is,

\[\dot{a} = - A * (1.0 + x)^{-g_2 - 1} / x^{g_1 - 1},\]

where \(x \equiv a / r_\mathrm{char}\) is the binary separation scaled to a characteristic transition radius (\(r_\mathrm{char}\)) between two power-law indices \(g_1\) and \(g_2\). There is also an overall normalization \(A\) that is calculated to yield the desired binary lifetimes.

The normalization for each binary, to produce the desired lifetime, is calculated as follows:

  1. A set of random test binary parameters are chosen.

  2. The normalization constants are determined, using least-squares optimization, to yield the desired lifetime.

  3. Interpolants are constructed to interpolate between the test binary parameters.

  4. The interpolants are called on the provided binary parameters, to calculate the interpolated normalization constants to reach the desired lifetimes.

Construction/Initialization: note that in addition to the standard Fixed_Time_2PL.__init__() constructor, there are two additional constructors are provided:

#! Using a callable for rchar probably doesnt work - _calculate_norm_interpolant looks like #! it only accepts a scalar value.

dadt_dedt(evo, step)[source]

Calculate hardening rate at the given integration step, for the given population.

Parameters:
  • evo (Evolution instance) – The evolutionary tracks of the binary population, providing binary parameters.

  • step (int,) – Integration step at which to calculate hardening rates.

Returns:

  • dadt ((N,) np.ndarray) – Binary hardening rates in units of [cm/s].

  • dedt ((N,) np.ndarray or None) – Rate-of-change of eccentricity, which is not included in this calculation, it is zero. None is returned if the input eccen is None.

classmethod from_pop(pop, time, **kwargs)[source]

Initialize a Fixed_Time_2PL instance using a provided _Discrete_Population instance.

Parameters:
  • pop (_Discrete_Population) – Input population, from which to use masses, redshifts and separations.

  • time (float, callable or array_like) –

    Total merger time of binaries, units of [sec], specifiable in the following ways:

    • float : uniform merger time for all binaries

    • callable : function time(mtot, mrat, redz) which returns the total merger time

    • array_like : (N,) matching the shape of mtot (etc) giving the merger time for each binary

  • **kwargs (dict) – Additional keyword-argument pairs passed to the Fixed_Time_2PL initialization method.

Returns:

Instance configured for the given binary population.

Return type:

Fixed_Time_2PL

classmethod from_sam(sam, time, sepa_init=3.0856775814913676e+22, **kwargs)[source]

Initialize a Fixed_Time_2PL instance using a provided Semi_Analytic_Model instance.

Parameters:
  • sam (holodeck.sam.Semi_Analytic_Model) – Input population, from which to use masses, redshifts and separations.

  • time (float, callable or array_like) –

    Total merger time of binaries, units of [sec], specifiable in the following ways:

    • float : uniform merger time for all binaries

    • callable : function time(mtot, mrat, redz) which returns the total merger time

    • array_like : (N,) matching the shape of mtot (etc) giving the merger time for each binary

  • sepa_init (float or array_like) –

    Initial binary separation. Units of [cm].

    • float : initial separation applied to all binaries,

    • array_like : initial separations for all binaries, shaped (N,) matching the number binaries.

  • **kwargs (dict) – Additional keyword-argument pairs passed to the Fixed_Time_2PL initialization method.

Returns:

Instance configured for the given binary population.

Return type:

Fixed_Time_2PL

classmethod function(norm, xx, gamma_inner, gamma_outer)[source]

Hardening rate given the parameters for this hardening model.

The functional form is,

\[\dot{a} = - A * (1.0 + x)^{-g_out + g_in} / x^{g_in - 1},\]

Where \(A\) is an overall normalization, and \(x = a / r_\textrm{char}\) is the binary separation scaled to a characteristic transition radius (\(r_\textrm{char}\)) between two power-law indices \(g_inner\) and \(g_outer\).

Parameters:
  • norm (array_like) – Hardening rate normalization, units of [cm/s].

  • xx (array_like) – Dimensionless binary separation, the semi-major axis in units of the characteristic (i.e. transition) radius of the model rchar.

  • gamma_inner (scalar) – Power-law of hardening timescale in the inner regime (small separations: r < rchar).

  • gamma_outer (scalar) – Power-law of hardening timescale in the outer regime (large separations: r > rchar).

class holodeck.hardening.Fixed_Time_2PL_SAM(sam, time, sepa_init=3.0856775814913673e+21, rchar=3.0856775814913675e+19, gamma_inner=-1.0, gamma_outer=1.5, num_steps=300)[source]

SAM-Optimized version of Fixed_Time_2PL: binary evolution for a fixed total lifetime.

class holodeck.hardening.Hard_GW[source]

Gravitational-wave driven binary hardening.

static dadt(mtot, mrat, sepa, eccen=None)[source]

Calculate GW Hardening rate of semi-major-axis vs. time.

See [Peters1964], Eq. 5.6

Parameters:
  • mtot (array_like) – Total mass of each binary system. Units of [gram].

  • mrat (array_like) – Mass ratio of each binary, defined as \(q \equiv m_1/m_2 \leq 1.0\).

  • sepa (array_like) – Binary semi-major axis (separation), in units of [cm].

  • eccen (array_like or None) – Binary eccentricity, None is the same as zero eccentricity (circular orbit).

Returns:

dadt – Hardening rate in semi-major-axis, result is negative, units [cm/s].

Return type:

np.ndarray

static dadt_dedt(evo, step)[source]

Calculate GW binary evolution (hardening rate) in semi-major-axis and eccentricity.

Parameters:
  • evo (Evolution) – Evolution instance providing the binary parameters for calculating hardening rates.

  • step (int) – Evolution integration step index from which to load binary parameters. e.g. separations are loaded as evo.sepa[:, step].

Returns:

  • dadt (np.ndarray) – Hardening rate in semi-major-axis, returns negative value, units [cm/s].

  • dedt (np.ndarray) – Hardening rate in eccentricity, returns negative value, units [1/s].

static deda(sepa, eccen)[source]

Rate of eccentricity change versus separation change.

See [Peters1964], Eq. 5.8

Parameters:
  • sepa (array_like,) – Binary semi-major axis (i.e. separation) [cm].

  • eccen (array_like,) – Binary orbital eccentricity.

Returns:

rv – Binary deda rate [1/cm] due to GW emission. Values are always positive.

Return type:

array_like,

static dedt(mtot, mrat, sepa, eccen=None)[source]

Calculate GW Hardening rate of eccentricity vs. time.

See [Peters1964], Eq. 5.7

If eccen is None, zeros are returned.

Parameters:
  • mtot (array_like) – Total mass of each binary system. Units of [gram].

  • mrat (array_like) – Mass ratio of each binary, defined as \(q \equiv m_1/m_2 \leq 1.0\).

  • sepa (array_like) – Binary semi-major axis (separation), in units of [cm].

  • eccen (array_like or None) – Binary eccentricity, None is the same as zero eccentricity (circular orbit).

Returns:

dedt – Hardening rate in eccentricity, result is <= 0.0, units [1/s]. Zero values if eccen is None.

Return type:

np.ndarray

class holodeck.hardening.Sesana_Scattering(gamma_dehnen=1.0, mmbulge=None, msigma=None)[source]

Binary-Hardening Rates calculated based on the Sesana stellar-scattering model.

This module uses the stellar-scattering rate constants from the fits in [Sesana2006] using the _SHM06 class. Scattering is assumed to only be effective once the binary is bound. An exponential cutoff is imposed at larger radii.

dadt_dedt(evo, step)[source]

Stellar scattering hardening rate.

Parameters:
  • evo (Evolution) – Evolution instance providing binary parameters at the given intergration step.

  • step (int) – Integration step at which to calculate hardening rates.

Returns:

  • dadt (array_like) – Binary hardening rates in units of [cm/s], defined to be negative.

  • dedt (array_like) – Binary rate-of-change of eccentricity in units of [1/sec].