evol_ecc

legwork.evol.evol_ecc(ecc_i, t_evol=None, n_step=100, timesteps=None, beta=None, m_1=None, m_2=None, a_i=None, f_orb_i=None, output_vars=['ecc', 'f_orb'], n_proc=1, avoid_merger=True, exact_t_merge=False, t_before=<Quantity 1. Myr>, t_merge=None)[source]

Evolve an array of eccentric binaries for t_evol time

This function use Peters & Mathews (1964) Eq. 5.11 and 5.13.

Note that all of {beta, m_1, m_2, ecc_i, a_i, f_orb_i} must have the same dimensions.

Parameters
ecc_ifloat/array

Initial eccentricity

t_evolfloat/array

Amount of time for which to evolve each binaries. Required if timesteps is None. Defaults to merger times.

n_stepsint

Number of timesteps to take between t=0 and t=``t_evol``. Required if timesteps is None. Defaults to 100.

timestepsfloat/array

Array of exact timesteps to take when evolving each binary. Must be monotonically increasing and start with t=0. Either supply a 1D array to use for every binary or a 2D array that has a different array of timesteps for each binary. timesteps is used in place of t_evol and n_steps and takes precedence over them.

betafloat/array

Constant defined in Peters and Mathews (1964) Eq. 5.9. See legwork.utils.beta() (if supplied m_1 and m_2 are ignored)

m_1float/array

Primary mass (required if beta is None or if output_vars contains a frequency)

m_2float/array

Secondary mass (required if beta is None or if output_vars contains a frequency)

a_ifloat/array

Initial semi-major axis (if supplied f_orb_i is ignored)

f_orb_ifloat/array

Initial orbital frequency (required if a_i is None)

output_varsarray

List of ordered output vars, choose from any of timesteps, ecc, a, f_orb and f_GW for which of timesteps, eccentricity, semi-major axis and orbital/GW frequency that you want. Default is [ecc, f_orb]

n_procint

Number of processors to split eccentricity evolution over, where the default is n_proc=1

avoid_mergerboolean

Whether to avoid integration around the merger of the binary. Warning: setting this to false will result in many LSODA errors to be outputted since the derivatives get so large.

exact_t_mergeboolean

Whether to calculate the merger time exactly or use a fit (only relevant when avoid_merger is set to True

t_beforefloat

How much time before the merger to cutoff the integration (default is 1 Myr - this will prevent all LSODA warnings for e < 0.95, you may need to increase this time if your sample is more eccentric than this)

t_mergefloat/array

Merger times for each source to be evolved. Only used when avoid_merger=True. If None then these will be automatically calculated either approximately or exactly based on the values of exact_t_merge.

Returns
evolutionarray

Array possibly containing eccentricity, semi-major axis, timesteps and frequency evolution. Content determined by output_vars