Isochrone

class artpop.stars.Isochrone(mini, mact, mags, eep=None, log_L=None, log_Teff=None)[source]

Bases: object

Class for storing generic isochrones. It also has methods for calculating IMF-weighted properties of a simple stellar population with the given age and metallicity.

Parameters
minindarray

Initial stellar masses. The masses must be monotonically increasing.

mactndarray

Actual stellar masses after mass loss.

mag_tableTable, dict, or structured ndarray

The stellar magnitudes.

eepndarray, optional

The Primary Equivalent Evolutionary Points. Needed to identify phases of stellar evolution.

log_Lndarray, optional

Stellar luminosities.

log_Teffndarray, optional

Stellar effective temperatures.

Attributes Summary

filters

List of filters in the given photometric system(s).

m_max

The maximum mass of the isochrone.

m_min

The minimum mass of the isochrone.

Methods Summary

calculate_mag_limit(imf, bandpass[, ...])

Calculate the limiting faint magnitude to sample a given fraction of mass or number of a stellar population.

from_parsec(file_name[, log_age, zini, ...])

Read isochrone generated from the PARSEC website.

imf_weights(imf[, m_min_norm, m_max_norm, ...])

Calculate IMF weights.

interpolate(y_name, x_interp[, x_name, ...])

Interpolate isochrone.

mag_to_mass(mag, bandpass)

Interpolate isochrone to calculate initial stellar masses that have the given magnitude.

nearest_mini(m)

Find the nearest mass to m in the initial mass array.

ssp_color(blue, red[, imf])

Calculate IMF-weighted integrated color.

ssp_mag(bandpass[, imf, norm_type])

Calculate IMF-weighted magnitude.

ssp_mean_star_mass(imf)

Calculate IMF-weighted mean stellar mass, where the IMF is normalized by number to one star formed.

ssp_sbf_mag(bandpass[, imf])

Calculate IMF-weighted SBF magnitude.

ssp_surviving_mass(imf[, m_min, m_max, ...])

Calculate IMF-weighted mass normalized such that 1 solar mass is formed by the age of the population.

Attributes Documentation

filters

List of filters in the given photometric system(s).

m_max

The maximum mass of the isochrone.

m_min

The minimum mass of the isochrone.

Methods Documentation

calculate_mag_limit(imf, bandpass, frac_mass_sampled=None, frac_num_sampled=None, distance=<Quantity 10. pc>)[source]

Calculate the limiting faint magnitude to sample a given fraction of mass or number of a stellar population. This is used for when you only what to sample stars that are brighter than a magnitude limit. You must specify either frac_mass_sampled or frac_num_sampled.

Note

You must give frac_mass_sampled or frac_num_sampled.

Parameters
imfstr or dict

Which IMF to use, if str then must be one of pre-defined: ‘kroupa’, ‘scalo’ or ‘salpeter’. Can also specify custom (broken) power law as dict, which must contain either ‘a’ as a Float (describing the slope of a single power law) or ‘a’ (a list with 3 elements describing the slopes of a broken power law) and ‘b’ (a list with 2 elements describing the locations of the breaks).

bandpassstr

Observation bandpass of the limiting magnitude.

frac_mass_sampled: float, optional

Fraction of mass that will be sampled if only stars that are brighter than the magnitude limit are sampled.

frac_num_sampled: float, optional

Fraction of stars by number that will be sampled if only stars that are brighter than the magnitude limit are sampled.

distancefloat or Quantity, optional

Distance to source. If float is given, the units are assumed to be Mpc. Default distance is 10 pc (i.e., the mags are in absolute units).

Returns
mag_limitfloat

The desired limiting magnitude.

mass_limitfloat

Mass of stars (in solar masses) that have the limiting magnitude.

static from_parsec(file_name, log_age=None, zini=None, num_filters=None)[source]

Read isochrone generated from the PARSEC website. If more than one age and/or metallicity is included in the file, you must provide the log_age and/or zini parameters.

Parameters
file_namestr

Isochrone file name.

log_agefloat, optional

Log of age in years. You must provided this parameter if there is more than one age in the isochrone file. Note this function does not interpolate ages, so it must be included in the file.

zinifloat, optional

Initial metal fraction. You must provided this parameter if there is more than one metallicity in the isochrone file. Note this function does not interpolate metallicity, so it must be included in the file.

num_filtersint, optional

Number of filters included in the isochrone file. If None, will assume the last non-filter parameter is mbolmag.

Returns
isoartpop.stars.Isochrone

PARSEC Isochrone object.

imf_weights(imf, m_min_norm=None, m_max_norm=None, norm_type='mass', **kwargs)[source]

Calculate IMF weights.

Parameters
imfstr or dict

Which IMF to use, if str then must be one of pre-defined: ‘kroupa’, ‘scalo’ or ‘salpeter’. Can also specify custom (broken) power law as dict, which must contain either ‘a’ as a Float (describing the slope of a single power law) or ‘a’ (a list with 3 elements describing the slopes of a broken power law) and ‘b’ (a list with 2 elements describing the locations of the breaks). of the breaks.

m_min_normNone or float, optional

Minimum mass for the normalization. Must be less than or equal to the mini mass of isochrone, which will be used if None is given.

m_max_normfloat, optional

Maximum mass for normalization.

norm_typestr, optional

Type of IMF normalization (mass or number).

Returns
wghtndarray

IMF weights calculated such that an integral over mass is simply given by SUM(m * wght).

interpolate(y_name, x_interp, x_name='mini', slice_interp=slice(None, None, None), **kwargs)[source]

Interpolate isochrone.

Parameters
y_namestr

Parameter name of interpolated variable.

x_interpndarray

New x values to interpolate over.

x_namestr

Parameter name of x values. Usually want this to be initial mass.

Returns
y_interpndarray

Interpolated y values.

mag_to_mass(mag, bandpass)[source]

Interpolate isochrone to calculate initial stellar masses that have the given magnitude. Since more than one mass can have the same magnitude (from stars being in different phases of stellar evolution), we step through each isochrone and check if the given magnitude falls between two bins, in which case we interpolate the associated mass.

Parameters
magfloat

Magnitude to be converted in to stellar mass(es).

bandpassstr

Photometric bandpass of the mag.

Returns
mass_interpndarray

Interpolated mass(es) associated with mag.

nearest_mini(m)[source]

Find the nearest mass to m in the initial mass array.

Parameters
mfloat:

The mass for which we want the nearest initial mass.

Returns
m_nearestfloat

Value of the nearest mass.

arg_nearestint

Argument of the nearest mass.

ssp_color(blue, red, imf='kroupa', **kwargs)[source]

Calculate IMF-weighted integrated color.

Parameters
bluestr

Blue bandpass.

redstr

Red bandpass.

imfstr or dict

Which IMF to use, if str then must be one of pre-defined: ‘kroupa’, ‘scalo’ or ‘salpeter’. Can also specify custom (broken) power law as dict, which must contain either ‘a’ as a Float (describing the slope of a single power law) or ‘a’ (a list with 3 elements describing the slopes of a broken power law) and ‘b’ (a list with 2 elements describing the locations of the breaks).

Returns
colorfloat

Integrated color of stellar population.

ssp_mag(bandpass, imf='kroupa', norm_type='mass', **kwargs)[source]

Calculate IMF-weighted magnitude.

Parameters
bandpassstr

Bandpass to of mean mag.

imfstr or dict

Which IMF to use, if str then must be one of pre-defined: ‘kroupa’, ‘scalo’ or ‘salpeter’. Can also specify custom (broken) power law as dict, which must contain either ‘a’ as a Float (describing the slope of a single power law) or ‘a’ (a list with 3 elements describing the slopes of a broken power law) and ‘b’ (a list with 2 elements describing the locations of the breaks).

norm_typestr, optional

Normalization type (mass or number)

Returns
mfloat

Integrated magnitude of stellar population.

ssp_mean_star_mass(imf)[source]

Calculate IMF-weighted mean stellar mass, where the IMF is normalized by number to one star formed. The mean stellar mass is equivalent to dividing a large of samples (distributed according to the IMF) by the number of stars.

Parameters
imfstr or dict

Which IMF to use, if str then must be one of pre-defined: ‘kroupa’, ‘scalo’ or ‘salpeter’. Can also specify custom (broken) power law as dict, which must contain either ‘a’ as a Float (describing the slope of a single power law) or ‘a’ (a list with 3 elements describing the slopes of a broken power law) and ‘b’ (a list with 2 elements describing the locations of the breaks).

Returns
——-
mean_massfloat

The mean stellar mass.

ssp_sbf_mag(bandpass, imf='kroupa', **kwargs)[source]

Calculate IMF-weighted SBF magnitude.

Parameters
bandpassstr

Bandpass to of SBF mag.

imfstr or dict

Which IMF to use, if str then must be one of pre-defined: ‘kroupa’, ‘scalo’ or ‘salpeter’. Can also specify custom (broken) power law as dict, which must contain either ‘a’ as a Float (describing the slope of a single power law) or ‘a’ (a list with 3 elements describing the slopes of a broken power law) and ‘b’ (a list with 2 elements describing the locations of the breaks).

Returns
sbffloat

SBF magnitude of stellar population.

ssp_surviving_mass(imf, m_min=None, m_max=120, add_remnants=True, mlim_bh=40.0, mlim_ns=8.5)[source]

Calculate IMF-weighted mass normalized such that 1 solar mass is formed by the age of the population.

The initial-mass-dependent remnant formulae are taken from Renzini & Ciotti 1993.

Parameters
imfstr or dict

Which IMF to use, if str then must be one of pre-defined: ‘kroupa’, ‘scalo’ or ‘salpeter’. Can also specify custom (broken) power law as dict, which must contain either ‘a’ as a Float (describing the slope of a single power law) or ‘a’ (a list with 3 elements describing the slopes of a broken power law) and ‘b’ (a list with 2 elements describing the locations of the breaks).

add_remnantsbool

If True, add stellar remnants to surviving mass.

Returns
massfloat

Total surviving mass.