ModelGridInterpolator

In practice, interaction with the model grid and bolometric correction objects is easiest through a ModelGridInterpolator object, which brings the two together. This object is the replacement of the Isochrone object from previous generations of this package, though it has a slightly different API. It is mostly backward compatible, except for the removal of the .mag function dictionary for interpolating apparent magnitudes, this being replaced by the .interp_mag method.

Isochrones

An IsochroneInterpolator object takes [EEP, log(age), feh] as parameters.

[1]:
from isochrones.mist import MIST_Isochrone

mist = MIST_Isochrone()

pars = [353, 9.78, -1.24]  # eep, log(age), feh
mist.interp_value(pars, ['mass', 'radius', 'Teff'])
[1]:
array([7.93829519e-01, 7.91444054e-01, 6.30305932e+03])

To interpolate apparent magnitudes, add distance [pc] and \(A_V\) extinction as parameters.

[2]:
mist.interp_mag(pars + [200, 0.11], ['K', 'BP', 'RP'])  # Returns Teff, logg, feh, mags
[2]:
(6303.059322477636,
 4.540738764316164,
 -1.377262817643937,
 array([10.25117074, 11.73997159, 11.06529993]))

Evolution tracks

Note that you can do the same using an EvolutionTrackInterpolator rather than an isochrone grid, using [mass, EEP, feh] as parameters:

[3]:
from isochrones.mist import MIST_EvolutionTrack

mist_track = MIST_EvolutionTrack()

pars = [0.794, 353, -1.24]  # mass, eep, feh [matching above]
mist_track.interp_value(pars, ['mass', 'radius', 'Teff', 'age'])
[3]:
array([7.93843749e-01, 7.91818696e-01, 6.31006708e+03, 9.77929505e+00])
[4]:
mist_track.interp_mag(pars + [200, 0.11], ['K', 'BP', 'RP'])
[4]:
(6310.067080800683,
 4.54076772643659,
 -1.372925841944066,
 array([10.24893319, 11.73358578, 11.06056746]))

There are also convenience methods (for both isochrones and tracks) if you prefer (and for backward compatibility—note that the parameters must be unpacked, unlike the calls to .interp_value and .interp_mag), though it is slower to call multiple of these than to call .interp_value once with several desired outputs:

[5]:
mist_track.mass(*pars)
[5]:
array(0.79384375)

You can also get the dataframe of a single isochrone (interpolated to any age or metallicity) as follows:

[6]:
mist.isochrone(9.53, 0.1).head()  # just show first few rows
[6]:
eep age feh mass initial_mass radius density logTeff Teff logg ... H_mag K_mag G_mag BP_mag RP_mag W1_mag W2_mag W3_mag TESS_mag Kepler_mag
223 223.0 9.53 0.150280 0.143050 0.143050 0.174516 42.182044 3.477544 3003.536405 5.121475 ... 8.785652 8.559155 12.766111 14.751368 11.522764 8.398324 8.200245 8.032482 11.381237 12.864034
224 224.0 9.53 0.150322 0.147584 0.147584 0.178799 40.088758 3.479902 3019.769652 5.112821 ... 8.713187 8.487450 12.662468 14.612205 11.426131 8.327414 8.129809 7.964879 11.287794 12.755405
225 225.0 9.53 0.150371 0.152520 0.152521 0.183594 37.948464 3.482375 3036.910262 5.103613 ... 8.635963 8.411037 12.552453 14.464800 11.323512 8.251886 8.054820 7.892865 11.188540 12.640135
226 226.0 9.53 0.150419 0.157318 0.157319 0.184463 37.208965 3.480519 3024.116433 5.101786 ... 8.629300 8.403586 12.569050 14.507862 11.334820 8.243224 8.045057 7.881000 11.197325 12.660600
227 227.0 9.53 0.150468 0.161795 0.161796 0.189168 35.381629 3.482801 3040.176145 5.093340 ... 8.558717 8.333774 12.467864 14.371759 11.240553 8.174286 7.976668 7.815386 11.106209 12.554499

5 rows × 27 columns

Generating synthetic properties

Often one wants to use stellar model grids to generate synthetic properties of stars. This can be done in a couple different ways, depending on what information you are able to provide. If you happen to have EEP values, you can use the fact that a ModelGridInterpolator is callable. Note that it takes the same parameters as all the other interpolation calls, with distance and AV as optional keyword parameters.

[7]:
from isochrones.mist import MIST_EvolutionTrack

mist_track = MIST_EvolutionTrack()
mist_track([0.8, 0.9, 1.0], 350, 0.0, distance=100, AV=0.1)
[7]:
nu_max logg eep initial_mass radius logTeff mass density Mbol phase ... H_mag K_mag G_mag BP_mag RP_mag W1_mag W2_mag W3_mag TESS_mag Kepler_mag
0 4254.629601 4.548780 350.0 0.8 0.787407 3.707984 0.799894 2.309938 5.792554 0.0 ... 9.040105 8.972502 10.872154 11.328425 10.258543 8.945414 8.989254 8.921756 10.247984 10.773706
1 3622.320906 4.495440 350.0 0.9 0.888064 3.741043 0.899876 1.811405 5.200732 0.0 ... 8.667003 8.614974 10.224076 10.602874 9.678976 8.593946 8.622577 8.575349 9.671007 10.129692
2 3041.107996 4.432089 350.0 1.0 1.006928 3.766249 0.999860 1.380733 4.675907 0.0 ... 8.312159 8.270380 9.679997 10.005662 9.186910 8.253638 8.269467 8.238306 9.180275 9.590731

3 rows × 29 columns

Often, however, you will not know the EEP values at which you wish to simulate your synthetic population. In this case, you can use the .generate() method.

[8]:
mist_track.generate([0.81, 0.91, 1.01], 9.51, 0.01)
[8]:
nu_max logg eep initial_mass radius logTeff mass density Mbol phase ... H K G BP RP W1 W2 W3 TESS Kepler
0 4787.598310 4.595858 320.808 0.81 0.750611 3.699978 0.809963 2.703461 5.977047 0.0 ... 4.154396 4.088644 5.988091 6.444688 5.375415 4.066499 4.117992 4.047535 5.365712 5.887722
1 3986.671794 4.535170 332.280 0.91 0.853120 3.737424 0.909935 2.066995 5.324246 0.0 ... 3.747329 3.699594 5.264620 5.632088 4.731978 3.684034 3.718112 3.670736 4.725020 5.169229
2 3154.677953 4.447853 343.800 1.01 0.993830 3.766201 1.009887 1.451510 4.705019 0.0 ... 3.322241 3.286761 4.620132 4.925805 4.148936 3.276062 3.295002 3.266166 4.143362 4.531319

3 rows × 29 columns

Under the hood, .generate() uses an interpolation step to approximate the eep value(s) corresponding to the requested value(s) of mass, age, and metallicity:

[9]:
mist_track.get_eep(1.01, 9.51, 0.01)
[9]:
343.8

Because this is fast, it is pretty inexpensive to generate a population of stars with given properties:

[10]:
import numpy as np

N = 10000
mass = np.ones(N) * 1.01
age = np.ones(N) * 9.82
feh = np.ones(N) * 0.02
%timeit mist_track.generate(mass, age, feh)
10 loops, best of 3: 112 ms per loop

Note though, that this interpolation doesn’t do great for evolved stars (this is the fundamental reason why isochrones always fits with EEP as one of the parameters). However, if you do want to compute more precise EEP values for given physical properties, you can set the accurate keyword parameter, which performs a function minimization:

[11]:
mist_track.get_eep(1.01, 9.51, 0.01, accurate=True)
[11]:
343.1963539123535

This is more accurate, but slow because it is actually performing a function minimization:

[12]:
%timeit mist_track.get_eep(1.01, 9.51, 0.01, accurate=True)
%timeit mist_track.get_eep(1.01, 9.51, 0.01)
100 loops, best of 3: 4.56 ms per loop
The slowest run took 4.98 times longer than the fastest. This could mean that an intermediate result is being cached.
100000 loops, best of 3: 4.26 µs per loop

Here we can see the effect of accuracy by plugging back in the estimated EEP into the interpolation:

[13]:
[mist_track.interp_value([1.01, e, 0.01], ['age']) for e in [343.8, 343.1963539123535]]
[13]:
[array([9.51806019]), array([9.50999994])]

So if accuracy is required, definitely use accurate=True, but for most purposes, the default should be fine. You can request that .generate() run in “accurate” mode, which uses this more expensive EEP computation (it will be correspondingly slower).

[14]:
mist_track.generate([0.81, 0.91, 1.01], 9.51, 0.01, accurate=True)
[14]:
nu_max logg eep initial_mass radius logTeff mass density Mbol phase ... H K G BP RP W1 W2 W3 TESS Kepler
0 4794.035436 4.596385 320.219650 0.81 0.750156 3.699863 0.809963 2.708365 5.979507 0.0 ... 4.156117 4.090301 5.990784 6.447681 5.377849 4.068141 4.119700 4.049167 5.368138 5.890400
1 3995.692509 4.536089 331.721363 0.91 0.852218 3.737300 0.909936 2.073560 5.327785 0.0 ... 3.750018 3.702214 5.268320 5.636100 4.735394 3.686635 3.720795 3.673334 4.728428 5.172899
2 3168.148566 4.449647 343.196354 1.01 0.991781 3.766083 1.009890 1.460523 4.710671 0.0 ... 3.327067 3.291533 4.625783 4.931724 4.154311 3.280826 3.299859 3.270940 4.148735 4.536929

3 rows × 29 columns

Just for curiosity, let’s look at the difference in the predictions:

[15]:
df0 = mist_track.generate([0.81, 0.91, 1.01], 9.51, 0.01, accurate=True)
df1 = mist_track.generate([0.81, 0.91, 1.01], 9.51, 0.01)
((df1 - df0) / df0).mean()
[15]:
nu_max         -0.002617
logg           -0.000240
eep             0.001760
initial_mass    0.000000
radius          0.001243
logTeff         0.000032
mass           -0.000002
density        -0.003716
Mbol           -0.000759
phase                NaN
feh            -0.057173
Teff            0.000273
logL            0.061576
delta_nu       -0.001803
interpolated         NaN
star_age        0.018487
age             0.000837
dt_deep        -0.007171
J              -0.000848
H              -0.000861
K              -0.000854
G              -0.000791
BP             -0.000792
RP             -0.000823
W1             -0.000854
W2             -0.000869
W3             -0.000857
TESS           -0.000823
Kepler         -0.000800
dtype: float64

Not too bad, for this example!

Demo: Visualize

Now let’s make sure that interpolated isochrones fall nicely between ones that are actually part of the grid. In order to execute this code, you will need to

conda install -c pyviz pyviz

and to execute in JupyterLab, you will need to

jupyter labextension install @pyviz/jupyterlab_pyviz
[16]:
import hvplot.pandas

iso1 = mist.model_grid.df.xs((9.5, 0.0), level=(0, 1))   # extract subgrid at log_age=9.5, feh=0.0
iso2 = mist.model_grid.df.xs((9.5, 0.25), level=(0, 1))  # extract subgrid at log_age=9.5, feh=0.25
iso3 = mist.isochrone(9.5, 0.12)  # should be between the other two

plot1 = iso1.hvplot.line('logTeff', 'logL', label='[Fe/H] = 0.0')
plot2 = iso2.hvplot.line('logTeff', 'logL', label='[Fe/H] = 0.25')
plot3 = iso3.hvplot.line('logTeff', 'logL', label='[Fe/H] = 0.12')

(plot1 * plot2 * plot3).options(invert_xaxis=True, legend_position='bottom_left', width=600)
[16]: