posydon.binary_evol.DT
posydon.binary_evol.DT.double_CO
Detached evolution for double compact-object binaries.
- class posydon.binary_evol.DT.double_CO.DoubleCO(n_o_steps_interval=None)[source]
Bases:
object
The double compact-object step class.
Initialize a DoubleCO instance.
posydon.binary_evol.DT.key_library
Dictionary containing data column name (key) translations between POSYDON h5 file PSyGrid data names (items) and MESA data names (keys).
- posydon.binary_evol.DT.key_library.DEFAULT_PROFILE_KEYS = ('radius', 'mass', 'logRho', 'energy', 'x_mass_fraction_H', 'y_mass_fraction_He', 'z_mass_fraction_metals', 'neutral_fraction_H', 'neutral_fraction_He', 'avg_charge_He')
# states to ID an HMS star # TODO: build these from the other lists to (hopefully) # ensure consistency STAR_STATES_FOR_HMS_MATCHING = [“H-rich_Core_H_burning”,
“accreted_He_Core_H_burning”]
# states to ID a postMS star STAR_STATES_FOR_postMS_MATCHING = [
“H-rich_Shell_H_burning”, “H-rich_Core_He_burning”, “H-rich_Central_He_depleted”, “H-rich_Core_C_burning”, “H-rich_Central_C_depletion”, “H-rich_non_burning”, “accreted_He_non_burning”]
# states to ID an He star STAR_STATES_FOR_strippedHe_MATCHING = [
‘accreted_He_Core_He_burning’, ‘stripped_He_Core_He_burning’, ‘stripped_He_Shell_He_burning’, # includes stars burning C in core, does this exist? ‘stripped_He_Central_He_depleted’, # includes stars burning C in core ‘stripped_He_Central_C_depletion’, ‘stripped_He_non_burning’ ]
posydon.binary_evol.DT.step_detached
Detached evolution step.
- class posydon.binary_evol.DT.step_detached.detached_evolution(primary=None, secondary=None, do_wind_loss=True, do_tides=True, do_magnetic_braking=True, magnetic_braking_mode='RVJ83', do_stellar_evolution_and_spin_from_winds=True, do_gravitational_radiation=True, verbose=False)[source]
Bases:
object
- ev_rel_rlo1(t, y, primary, secondary)[source]
Relative difference between radius and Roche lobe. Used to
check if there is RLOF mass transfer during the detached binary evolution interpolation. Calculated for the secondary.
- Parameters:
t (float) – Time of the evolution, in years.
y (tuple(float)) – [separation, eccentricity] at that time. Separation should be in solar radii.
primary (SingleStar object) – A single star object, representing the primary (more evolved) star in the binary and containing its properties.
secondary (SingleStar object) – A single star object, representing the secondary (less evolved) star in the binary and containing its properties.
- Returns:
RL_rel_diff – Relative difference between stellar radius and Roche lobe radius.
- Return type:
- ev_rel_rlo2(t, y, primary, secondary)[source]
Relative difference between radius and Roche lobe. Used to
check if there is RLOF mass transfer during the detached binary evolution interpolation. Calculated for the primary.
- Parameters:
t (float) – Time of the evolution, in years.
y (tuple(float)) – [separation, eccentricity] at that time. Separation should be in solar radii.
primary (SingleStar object) – A single star object, representing the primary (more evolved) star in the binary and containing its properties.
secondary (SingleStar object) – A single star object, representing the secondary (less evolved) star in the binary and containing its properties.
- Returns:
RL_rel_diff – Relative difference between stellar radius and Roche lobe radius.
- Return type:
- ev_rlo1(t, y, primary, secondary)[source]
Difference between radius and Roche lobe at a given time. Used
to check if there is RLOF mass transfer during the detached binary evolution interpolation. Calculated for the secondary.
- Parameters:
t (float) – Time of the evolution, in years.
y (tuple(float)) – [separation, eccentricity] at that time. Separation should be in solar radii.
primary (SingleStar object) – A single star object, representing the primary (more evolved) star in the binary and containing its properties.
secondary (SingleStar object) – A single star object, representing the secondary (less evolved) star in the binary and containing its properties.
- Returns:
RL_diff – Difference between stellar radius and 95% of the Roche lobe radius in solar radii.
- Return type:
- ev_rlo2(t, y, primary, secondary)[source]
Difference between radius and Roche lobe at a given time. Used
to check if there is RLOF mass transfer during the detached binary evolution interpolation. Calculated for the primary.
- Parameters:
t (float) – Time of the evolution, in years
y (tuple(float)) – [separation, eccentricity] at that time. Separation should be in solar radii.
primary (SingleStar object) – A single star object, representing the primary (more evolved) star in the binary and containing its properties.
secondary (SingleStar object) – A single star object, representing the secondary (less evolved) star in the binary and containing its properties.
- Returns:
RL_diff – Difference between stellar radius and 95% of the Roche lobe radius in solar radii.
- Return type:
- class posydon.binary_evol.DT.step_detached.detached_step(dt=None, n_o_steps_history=None, do_wind_loss=True, do_tides=True, do_gravitational_radiation=True, do_magnetic_braking=True, magnetic_braking_mode='RVJ83', do_stellar_evolution_and_spin_from_winds=True, RLO_orbit_at_orbit_with_same_am=False, record_matching=False, verbose=False, grid_name_Hrich=None, grid_name_strippedHe=None, metallicity=None, path=PATH_TO_POSYDON_DATA, matching_method='minimize', list_for_matching_HMS=None, list_for_matching_postMS=None, list_for_matching_HeStar=None)[source]
Bases:
object
Evolve a detached binary.
The binary will be evolved until Roche-lobe overflow, core-collapse or maximum simulation time, using the standard equations that govern the orbital evolution.
- Parameters:
path (str) – Path to the directory that contains POSYDON data HDF5 files. Defaults to the PATH_TO_POSYDON_DATA environment variable. Used for track matching.
metallicity (float) –
The metallicity of the grid. This should be one of the eight supported metallicities:
[2e+00, 1e+00, 4.5e-01, 2e-01, 1e-01, 1e-02, 1e-03, 1e-04]
and this will be converted to a corresponding string (e.g., 1e+00 –> “1e+00_Zsun”). Used for track matching.
matching_method (str) –
Method to find the best match between a star from a previous step and a point in a single star evolution track. Options:
- ”root”: Tries to find a root of two matching quantities. It is
possible to not find one, causing the evolution to fail.
- ”minimize”: Minimizes the sum of squares of differences of
various quantities between the previous evolution step and a stellar evolution track.
Used for track matching.
grid_name_Hrich (str) –
Name of the single star H-rich grid h5 file, including its parent directory. This is set to (for example):
grid_name_Hrich = ‘single_HMS/1e+00_Zsun.h5’
by default if not specified. Used for track matching.
grid_name_strippedHe (str) –
Name of the single star He-rich grid h5 file. This is set to (for example):
grid_name_strippedHe = ‘single_HeMS/1e+00_Zsun.h5’
by default if not specified. Used for track matching.
list_for_matching_HMS (list) – A list of mixed type that specifies properties of the matching process for HMS stars. Used for track matching.
list_for_matching_postMS (list) – A list of mixed type that specifies properties of the matching process for postMS stars. Used for track matching.
list_for_matching_HeStar (list) – A list of mixed type that specifies properties of the matching process for He stars. Used for track matching.
record_matching (bool) – Whether properties of the matched star(s) should be recorded in the binary evolution history. Used for track matching.
- KEYS
Contains keywords corresponding to MESA data column names which are used to extract quantities from the single star evolution grids.
- dt
The timestep size, in years, to be appended to the history of the binary. None means only the final step. Note: do not select very small timesteps because it may mess with the solving of the ODE.
- Type:
- n_o_steps_history
Alternatively, we can define the number of timesteps to be appended to the history of the binary. None means only the final step. If both dt and n_o_steps_history are different than None, dt has priority.
- Type:
- do_wind_loss
If True, take into account change of separation due to mass loss from the star.
- Type:
- do_tides
If True, take into account change of separation, eccentricity and star spin due to tidal forces.
- Type:
- do_gravitational_radiation
If True, take into account change of separation and eccentricity due to gravitational wave radiation.
- Type:
- do_magnetic_braking
If True, take into account change of star spin due to magnetic braking.
- Type:
- magnetic_braking_mode
- A string corresponding to the desired magnetic braking prescription.
– RVJ83: Rappaport, Verbunt, & Joss 1983 (Default) – M15: Matt et al. 2015 – G18: Garraffo et al. 2018 – CARB: Van & Ivanova 2019
- Type:
- do_stellar_evolution_and_spin_from_winds
If True, take into account change of star spin due to change of its moment of inertia during its evolution and due to spin angular momentum loss due to winds.
- Type:
- RLO_orbit_at_orbit_with_same_am
Binaries are circularized instaneously when RLO occurs and this option dictates how that is handled. If False (default), place the binary in an orbit with separation equal to the binary’s separation at periastron. If True, circularize the orbit assuming that angular momentum is conserved w.r.t. the previously (possibly) eccentric orbit. In the latter case, the star may no longer fill its Roche lobe after circularization, and may be further evolved until RLO commences once again, but without changing the orbit.
- Type:
- translate
Dictionary containing data column name (key) translations between POSYDON h5 file PSyGrid data names (items) and MESA data names (keys).
- Type:
- track_matcher
The TrackMatcher object performs functions related to matching binary stellar evolution components to single star evolution models.
- Type:
TrackMatcher object
Initialize the step. See class documentation for details.
- get_time_after_evo(res, binary)[source]
After detached evolution, this uses the ODESolver result
to determine what the current time is.
- Parameters:
res (ODESolver object) – This is the ODESolver object produced by SciPy’s solve_ivp function that contains calculated values of the stars evolution through the detached step.
binary (BinaryStar object) – A binary star object, containing the binary system’s properties.
- Returns:
t – This is the time elapsed as a result of detached evolution in years. This is a float unless the user specifies a timestep (see n_o_steps_history or dt) to use via the simulation properties ini file, in which case it is an array.
- Return type:
- init_evo_kwargs()[source]
Store keyword args required to initialize detached evolution, based on step’s kwargs.
- update_after_evo(res, t, binary, primary, secondary)[source]
Update star and binary properties and interpolators with
ODESolver result from detached evolution. This update gives the binary/stars their appropriate values, according to the interpolation after detached evolution.
- Parameters:
res (ODESolver object) – This is the ODESolver object produced by SciPy’s solve_ivp function that contains calculated values of the stars evolution through the detached step.
t (float or array[float]) – This is the time elapsed as a result of detached evolution in years. This is a float unless the user specifies a timestep to use via the simulation properties ini file, in which case it is an array.
binary (BinaryStar object) – A binary star object, containing the binary system’s properties.
primary (SingleStar object) – A single star object, representing the primary (more evolved) star in the binary and containing its properties.
secondary (SingleStar object) – A single star object, representing the secondary (less evolved) star in the binary and containing its properties.
- Warns:
InappropriateValueWarning – If trying to compute log angular momentum for object with no spin.
posydon.binary_evol.DT.step_disrupted
Merging and isolated evolution step.
- class posydon.binary_evol.DT.step_disrupted.DisruptedStep(grid_name_Hrich=None, grid_name_strippedHe=None, path=PATH_TO_POSYDON_DATA, *args, **kwargs)[source]
Bases:
IsolatedStep
Prepare a runaway star to do an an isolated_step)
Initialize the step. See class documentation for details.
posydon.binary_evol.DT.step_initially_single
Merging and isolated evolution step.
- class posydon.binary_evol.DT.step_initially_single.InitiallySingleStep(grid_name_Hrich=None, grid_name_strippedHe=None, path=PATH_TO_POSYDON_DATA, *args, **kwargs)[source]
Bases:
IsolatedStep
Prepare a runaway star to do an an isolated_step)
Initialize the step. See class documentation for details.
posydon.binary_evol.DT.step_isolated
Isolated evolution step.
- class posydon.binary_evol.DT.step_isolated.IsolatedStep(grid_name_Hrich=None, grid_name_strippedHe=None, path=PATH_TO_POSYDON_DATA, do_wind_loss=False, do_tides=False, do_gravitational_radiation=False, do_magnetic_braking=False, *args, **kwargs)[source]
Bases:
detached_step
Evolve an isolated star (a single star, a merger product, a runaway star, etc.)
The star will be matched in the beginning of the step and will be evolved until core-collapse or maximum simulation time, based on a grid of single star HDF5 grid.
Initialize the step. See class documentation for details.
posydon.binary_evol.DT.step_merged
Merging and isolated evolution step.
- class posydon.binary_evol.DT.step_merged.MergedStep(grid_name_Hrich=None, grid_name_strippedHe=None, path=PATH_TO_POSYDON_DATA, merger_critical_rot=0.4, rel_mass_lost_HMS_HMS=0.1, list_for_matching_HMS=[['mass', 'center_h1', 'he_core_mass'], [20.0, 1.0, 10.0], ['log_min_max', 'min_max', 'min_max'], [None, None], [0, None]], list_for_matching_postMS=[['mass', 'center_he4', 'he_core_mass'], [20.0, 1.0, 10.0], ['log_min_max', 'min_max', 'min_max'], [None, None], [0, None]], list_for_matching_HeStar=[['he_core_mass', 'center_he4'], [10.0, 1.0], ['min_max', 'min_max'], [None, None], [0, None]], *args, **kwargs)[source]
Bases:
IsolatedStep
Prepare a merging star to do an an IsolatedStep
Initialize the step. See class documentation for details.
- merged_star_properties(star_base, comp)[source]
Make assumptions about the core/total mass, and abundances of the star of a merged product.
Similar to the table of merging in BSE
- star_base: Single Star
is our base star that engulfs its companion. The merged star will have this star as a base
- comp: Single Star
is the star that is engulfed
posydon.binary_evol.DT.track_match
- class posydon.binary_evol.DT.track_match.TrackMatcher(grid_name_Hrich, grid_name_strippedHe, path=PATH_TO_POSYDON_DATA, metallicity=None, matching_method='minimize', list_for_matching_HMS=None, list_for_matching_postMS=None, list_for_matching_HeStar=None, record_matching=False, verbose=False)[source]
Bases:
object
This class contains the functionality to match binary star components
to single star models for detached evolution. Typically, this is so that the binary star can continue evolving in a detached (non-interacting) state.
Several matching methods may be used, and metrics (physical
quantities) that are used to determine the quality of the match may be customized. By default, the metrics are as follows:
- calc_omega(star)[source]
Calculate the spin of a star from its (pre-match) moment
of inertia and angular momentum (or rotation rates). This is required because we match a rotating model (from the binary grids) to a non-rotating model (from the single star grids).
- Parameters:
star (SingleStar object) – Star object containing the star properties.
- Returns:
omega_in_rad_per_year – The rotation rate of the star in radians per year, calculated from the star’s (pre-match) angular momentum and moment of inertia.
- Return type:
- Warns:
InappropriateValueWarning – If the pre-match rotation quantities are NaN or None and can not calculate a the post-match rotation rate, we setting the post-match rotation rate to zero.
- determine_star_states(binary)[source]
Determines which star is primary (further evolved) and which is
secondary (less evolved). Determines whether stars should be matched to the H- or He-rich grid, whether they exist, or if they are compact objects/massless remnants. This is used to determine how to match the stars.
- Parameters:
binary (BinaryStar object) – A binary star object, containing the binary system’s properties.
- Returns:
primary (SingleStar object) – A single star object, representing the primary (more evolved) star in the binary and containing its properties.
secondary (SingleStar object) – A single star object, representing the secondary (less evolved) star in the binary and containing its properties.
only_CO (bool) – A boolean indicating whether the binary system contains only a single compact object (True) or not (False). As in the detached step, it may be desirable to use this flag to exit an evolution step, as a single compact obejct (point mass) can not be evolved further.
- Raises:
POSYDONError – If there is no star to evolve (both are massless remnants), then evolution can no continue and detached step should not have been called.
ValueError – If the State of star 1 or 2 is not recognized.
POSYDONError – If the non_existent_companion of the binary is determined to be not equal to 0 (both stars exist), 1 (only star 2 exists), or 2 (only star 1 exists), something has gone wrong.
- do_matching(binary, step_name='step_match')[source]
Perform binary to single star grid matching. This is currently
used when transitioning to detached star evolution from binary but may be used in other steps. This performs several actions:
Determines which star is primary/secondary in the evolution and their evolutionary states. If evolvable, matching will proceed.
Match either one or both stars to a (non-rotating) single star evolution track.
Calculate the matched star’s rotation using the pre-match step’s angular momentum-related properties.
Returns the primary/secondary stars with interpolator objects that may be used to calculate quantities along the time series of each star for further evolution.
- Parameters:
binary (BinaryStar object) – A binary star object, containing the binary system’s properties.
step_name (str) – If self.record_matching is True, then the matched quantities of the star will be appended to the history. This is a string that can be used as a custom label in the BinaryStar object’s history, meant to indicate the relevant evolution step’s name. This should normally match the name of the step in which the matching was made, e.g., “step_detached”.
- Returns:
primary_out (tuple(SingleStar, PchipInterpolator, float)) – The first element is the SingleStar object of the primary (more evolved) star. The second element is the primary’s PchipInterpolator object, used to interpolate values along this star’s time series. The third is the primary’s calculated (post-match) rotation rate, using angular momentum-related quantites from the pre-match step.
secondary_out (tuple(SingleStar, PchipInterpolator, float)) – The first element is the SingleStar object of the secondary (less evolved) star. The second element is the secondary’s PchipInterpolator object, used to interpolate values along this star’s time series. The third is the secondary’s calculated (post-match) rotation rate, using angular momentum-related quantites from the pre-match step.
only_CO (bool) – A boolean indicating whether the binary system contains only a single compact object (True) or not (False). As in the detached step, it may be desirable to use this flag to exit an evolution step, as a single compact obejct (point mass) can not be evolved further.
- Raises:
ValueError – If the primary_normal and primary_not_normal flags are both determined to be False. One or the other should be True.
MatchingError – If the stellar matching to a single star model fails, or the PchipInterpolator object returned from matching is None.
- get_root0(attr_names, attr_vals, htrack, rescale_facs=None)[source]
Get the stellar evolution track in the single star grid with values
closest to the requested ones. This calculates the difference in stellar properties between a given track and those in the single star grids. It then returns the mass of and age along the track where the minimum difference occurs.
- Parameters:
attr_names (list[str]) – Contains the keys of the requested specific quantities that will be matched in the single star track.
attr_vals (list[float]) – Contains the latest values (from a previous POSYDON step) of the quantities of “keys” in the POSYDON SingleStar object. This should be the same length as “keys”.
htrack (bool) – Set True to search the single star H-rich grids, or False to search the He-rich grids.
rescale_facs (list[float]) – Contains normalization factors to be divided for rescaling attribute values. This should be the same length as attr_names.
- Returns:
m0 (float) – Initial mass (in solar units) of the minimum diff model for initial guess.
t0 (float) – Age (in years) of the minimum diff model for initial guess.
- get_star_final_values(star)[source]
This updates the final values of a SingleStar object,
given an initial stellar mass m0, typically found from matching to a single star track.
- Parameters:
star (SingleStar object) – A single star object that contains the star’s properties.
htrack (bool) – A boolean that specifies whether the star would be found in the hydrogen rich single star grid or not (in which case it is matched to the helium rich single star grid).
m0 (float) – Initial stellar mass (in solar units) of the single star track that we will grab values from and update star with.
- get_star_match_data(binary, star, copy_prev_m0=None, copy_prev_t0=None)[source]
Match a given component of a binary (i.e., a star) to a
single star model. This then creates and returns interpolator objects that may be used to calculate properties of the star as a function of time.
In the case of a compact object, radius, mdot, and Idot are
set to zero. One may use another star, e.g., the companion of the compact object to provide a mass and age.
- Parameters:
binary (BinaryStar object) – A binary star object, containing the binary system’s properties.
star (SingleStar object) – A single star object that contains the star’s properties.
copy_prev_m0 (float) – A mass value that may be copied from another star in the case where the target star is a compact object
copy_prev_t0 (float) – An age value that may be copied from another star in the case where the target star is a compact object
- Returns:
match_m0 (float) – Initial mass (in solar units) of the matched model.
match_t0 (float) – Age (in years) of the matched model.
- get_star_profile(star)[source]
This updates the stellar profile of a SingleStar object,
given an initial stellar mass m0, typically found from matching to a single star track. The profile of the SingleStar object is updated to become the profile of the (matched) single star track.
- Parameters:
star (SingleStar object) – A single star object that contains the star’s properties.
htrack (bool) – A boolean that specifies whether the star would be found in the hydrogen rich single star grid or not (in which case it is matched to the helium rich single star grid).
m0 (float) – Initial stellar mass (in solar units) of the single star track that we will grab values from and update star with.
- get_track_val(key, htrack, m0, t)[source]
Return a single value of a stellar property from the
interpolated time-series along a requested stellar track of mass m0 at an age of t.
- match_through_minimize(star, get_root0, get_track_val)[source]
Match the star through a minimization method. An initial
match is found using the get_root0() function. Then, those initial values are used in SciPy’s minimize method to find the closest matching point in evolution, based on several physical matching metrics.
- Parameters:
star (SingleStar object) – This is a SingleStar object, typically representing a member of a binary star system that we are trying to match to a single star track.
get_root0 (function) – Function that is used to get an initial guess for a closely matching stellar track.
get_track_val (function) – Function that returns a single value of a stellar property from the interpolated time-series along a requested stellar track of mass m0 at an age of t.
- Returns:
match_vals (array) – Mass and age of the star found through matching. These serve as starting points for subsequent (post-match) evolution.
best_sol (OptimizeResult object) – The OptimizeResult object that contains attributes of the closest matching stellar track. Produced by SciPy’s minimize().
- match_through_root(star, get_root0, get_track_val)[source]
Match the star through a root method.
- Parameters:
star (SingleStar object) – This is a SingleStar object, typically representing a member of a binary star system that we are trying to match to a single star track.
get_root0 (function) – Function that is used to get an initial guess for a closely matching stellar track.
get_track_val (function) – Function that returns a single value of a stellar property from the interpolated time-series along a requested stellar track of mass m0 at an age of t.
- Returns:
match_vals (array) – Mass and age of the star found through matching. These serve as starting points for subsequent (post-match) evolution.
best_sol (OptimizeResult object) – The OptimizeResult object that contains attributes of the closest matching stellar track. Produced by SciPy’s root().
- match_to_single_star(star)[source]
Get the track in the grid that matches the time and mass of a star,
that has typically undergone prior binary star evolution. A match is made according to a given algorithm that minimizes the difference amongst several physical properties, e.g., mass, central hydrogen abundance, radius, and core helium mass, depending on the type of star being matched. However, these properties may also be customized by the user.
- Parameters:
star (SingleStar object) – A single star object that contains the star’s properties.
- Returns:
match_m0 (float) – Initial mass (in solar units) of the matched model.
match_t0 (float) – Age (in years) of the matched model.
- Warns:
EvolutionWarning – If attempting to match an He-star with an H-rich grid or post-MS star with a stripped-He grid. This can happen if an initial matching attempt fails; alternative grids are checked for a match in such cases.
- Raises:
AttributeError – If a matching parameter does not exist in the single star grid options
NumericalError – If SciPy numerical differentiation occured outside boundary while matching to single star track.
MatchingError – If the initial mass to be matched is out of the single star grid range, thereby preventing a possible match.
- scale(attr_name, htrack, scaler_method)[source]
Normalize quantities in the single star grids to (0,1).
- Parameters:
attr_name (str) – Keyword of the requested quantity.
htrack (bool) – A boolean that specifies whether the star would be found in the hydrogen rich single star grid or not (in which case it is matched to the helium rich single star grid).
scaler_method (str) – Scaling method in the DataScaler class. See posydon.interpolation.data_scaling.DataScaler().fit() for more details.
- Returns:
scaler – This is a DataScaler object, trained to rescale the requested attribute to the range (0, 1).
- Return type:
DataScaler object
- update_rotation_info(primary, secondary)[source]
Once we have matched to a non-rotating single star, we need to
calculate what the single star’s spin should be, based the star’s rotation before matching. This calculates a new rotation rate for the matched star from the previous moment of inertia and angular momentum (or rotation rate in lieu of those). This also updates the stars in the binary with the newly calculated values to reflect this.
- Parameters:
primary (SingleStar object) – A single star object, representing the primary (more evolved) star in the binary and containing its properties.
secondary (SingleStar object) – A single star object, representing the secondary (less evolved) star in the binary and containing its properties.
- Returns:
omega0_pri (float) – The rotation rate of the primary star in radians per year, calculated from the star’s (pre-match) angular momentum and moment of inertia.
omega0_sec (float) – The rotation rate of the secondary star in radians per year, calculated from the star’s (pre-match) angular momentum and moment of inertia.
- update_star_properties(star, htrack)[source]
This updates a SingleStar object (star) with the
values from a single star track that has initial mass m0 and age t0. This can be used after matching finds the closest m0, t0 to update the SingleStar object with the values of the best matching single star track. This will not update the log_total_angular_momentum or surf_avg_omega because the single star track is non-rotating and those quantities are poorly defined.
- Parameters:
star (SingleStar object) – A single star object that contains the star’s properties.
htrack (bool) – A boolean that specifies whether the star would be found in the hydrogen rich single star grid or not (in which case it is matched to the helium rich single star grid).