/
fit_model.py
626 lines (493 loc) · 21 KB
/
fit_model.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
"""
Module with functionalities for fitting atmospheric model spectra.
"""
import os
import math
import warnings
from multiprocessing import Pool, cpu_count
import emcee
import numpy as np
try:
import pymultinest
except:
warnings.warn('PyMultiNest could not be imported.')
from species.core import constants
from species.data import database
from species.read import read_model, read_object
from species.util import read_util
def lnprior(param,
bounds,
modelpar,
prior=None):
"""
Internal function for the prior probability.
Parameters
----------
param : numpy.ndarray
Parameter values.
bounds : dict
Parameter boundaries.
modelpar : list(str, )
Parameter names.
prior : tuple(str, float, float), None
Gaussian prior on one of the parameters. Currently only possible for the mass, e.g.
('mass', 13., 3.) for an expected mass of 13 Mjup with an uncertainty of 3 Mjup. Not
used if set to None.
Returns
-------
float
Log prior probability.
"""
if prior is not None:
modeldict = {}
for i, item in enumerate(modelpar):
modeldict[item] = param[i]
ln_prior = 0
for i, item in enumerate(modelpar):
if bounds[item][0] <= param[i] <= bounds[item][1]:
if prior is not None and prior[0] == 'mass' and item == 'logg':
mass = read_util.get_mass(modeldict)
ln_prior += -0.5*(mass-prior[1])**2/prior[2]**2
else:
ln_prior += 0.
else:
ln_prior = -np.inf
break
return ln_prior
def lnlike(param,
modelpar,
objphot,
distance,
spectrum,
modelphot,
modelspec):
"""
Internal function for the likelihood function.
Parameters
----------
param : numpy.ndarray
Parameter values.
modelpar : list(str, )
Parameter names.
objphot : list(tuple(float, float), )
List with the fluxes and uncertainties of the object.
distance : tuple(float, float)
Distance and uncertainty (pc).
spectrum : dict
Dictionary with the spectrum stored as wavelength (um), flux (W m-2 um-1),
and error (W m-2 um-1), and optionally the covariance matrix and the inverse of
the covariance matrix.
modelphot : list(species.read.read_model.ReadModel, )
List with the interpolated synthetic photometry.
modelspec : list(species.read.read_model.ReadModel, )
List with the interpolated synthetic spectra.
Returns
-------
float
Log likelihood probability.
"""
paramdict = {}
spec_scaling = {}
for i, item in enumerate(modelpar):
if item == 'radius':
radius = param[i]
elif item[:8] == 'scaling_' and item[8:] in spectrum:
spec_scaling[item[8:]] = param[i]
else:
paramdict[item] = param[i]
scaling = (radius*constants.R_JUP)**2 / (distance[0]*constants.PARSEC)**2
chisq = 0.
if objphot is not None:
for i, item in enumerate(objphot):
flux = scaling * modelphot[i].spectrum_interp(list(paramdict.values()))
chisq += (item[0]-flux)**2 / item[1]**2
if spectrum is not None:
for i, item in enumerate(spectrum.keys()):
flux = scaling * modelspec[i].spectrum_interp(list(paramdict.values()))[0, :]
if item in spec_scaling:
flux_obs = spec_scaling[item]*spectrum[item][0][:, 1]
else:
flux_obs = spectrum[item][0][:, 1]
if spectrum[item][2] is not None:
spec_diff = flux_obs - flux
chisq += np.dot(spec_diff, np.dot(spectrum[item][2], spec_diff))
else:
chisq += np.nansum((flux_obs-flux)**2 / spectrum[item][0][:, 2]**2)
return -0.5*chisq
def lnprob(param,
bounds,
modelpar,
objphot,
distance,
prior,
spectrum,
modelphot,
modelspec):
"""
Internal function for the posterior probability.
Parameters
----------
param : numpy.ndarray
Parameter values.
bounds : dict
Parameter boundaries.
modelpar : list(str, )
Parameter names.
objphot : list(tuple(float, float), )
List with the fluxes and uncertainties of the object.
distance : tuple(float, float)
Distance and uncertainty (pc).
prior : tuple(str, float, float)
Gaussian prior on one of the parameters. Currently only possible for the mass, e.g.
('mass', 13., 3.) for an expected mass of 13 Mjup with an uncertainty of 3 Mjup. Not
used if set to None.
spectrum : dict
Wavelength (um), apparent flux (W m-2 um-1), and flux error (W m-2 um-1).
modelphot : list(species.read.read_model.ReadModel, )
List with the interpolated synthetic fluxes.
modelspec : list(species.read.read_model.ReadModel, )
List with the interpolated synthetic spectra.
Returns
-------
float
Log posterior probability.
"""
ln_prior = lnprior(param, bounds, modelpar, prior)
if math.isinf(ln_prior):
ln_prob = -np.inf
else:
ln_prob = ln_prior + lnlike(param,
modelpar,
objphot,
distance,
spectrum,
modelphot,
modelspec)
if np.isnan(ln_prob):
ln_prob = -np.inf
return ln_prob
class FitModel:
"""
Class for fitting atmospheric model spectra to photometric and/or spectroscopic data.
"""
def __init__(self,
object_name,
filters,
model,
bounds=None,
inc_phot=True,
inc_spec=True):
"""
For each photometric point and spectrum, the model grid is linearly interpolated at the
required synthetic photometry and wavelength sampling before running the MCMC. Therefore,
the computation time of this initial interpolation depends on the wavelength range and
spectral resolution of the spectra that are stored in the database, and the prior
boundaries that are chosen with ``bounds``.
Parameters
----------
object_name : str
Object name in the database as created with
:func:`~species.data.database.Database.add_object` or
:func:`~species.data.database.Database.add_companion`.
filters : tuple(str, )
Filter names for which the photometry is selected. All available photometry of the
object is selected if set to None.
model : str
Atmospheric model (e.g. 'drift-phoenix', 'petitcode-cool-cloudy', or 'bt-settl').
bounds : dict, None
Parameter boundaries. The full range is used for each parameter if set to None. In that
case, the radius range is set to 0.8-1.5 Rjup. It is also possible to specify the bounds
for a subset of the parameters, for example, ``{'radius': (0.5, 10.)}``. Restricting
the boundaries will decrease the computation time with the interpolation prior to the
MCMC sampling. An additional scaling parameter can be fitted for each spectrum in which
case, the boundaries have to be provided with the database tag of the spectrum.
For example, ``{'sphere_ifs': (0.5, 2.)}`` if the spectrum is stored as ``sphere_ifs``
with :func:`~species.data.database.Database.add_object`.
inc_phot : bool
Include photometric data in the fit.
inc_spec : bool
Include spectroscopic data in the fit.
Returns
-------
NoneType
None
"""
self.object = read_object.ReadObject(object_name)
self.distance = self.object.get_distance()
self.model = model
self.bounds = bounds
if not inc_phot and not inc_spec:
raise ValueError('No photometric or spectral data has been selected.')
if self.bounds is not None:
readmodel = read_model.ReadModel(self.model)
bounds_grid = readmodel.get_bounds()
for item in bounds_grid:
if item not in self.bounds:
self.bounds[item] = bounds_grid[item]
else:
readmodel = read_model.ReadModel(self.model, None, None)
self.bounds = readmodel.get_bounds()
if 'radius' not in self.bounds:
self.bounds['radius'] = (0.8, 1.5)
print('Prior and interpolation boundaries:')
for key, value in self.bounds.items():
print(f' - {key} = {value}')
if inc_phot:
self.objphot = []
self.modelphot = []
if filters is None:
species_db = database.Database()
objectbox = species_db.get_object(object_name, None)
filters = objectbox.filters
for item in filters:
print(f'Interpolating {item}...', end='', flush=True)
readmodel = read_model.ReadModel(self.model, filter_name=item)
readmodel.interpolate_grid(self.bounds)
self.modelphot.append(readmodel)
print(f' [DONE]')
obj_phot = self.object.get_photometry(item)
self.objphot.append((obj_phot[2], obj_phot[3]))
else:
self.objphot = None
self.modelphot = None
if inc_spec:
self.spectrum = self.object.get_spectrum()
self.modelspec = []
for key, value in self.spectrum.items():
print(f'\rInterpolating {key}...', end='', flush=True)
wavel_range = (0.9*value[0][0, 0], 1.1*value[0][-1, 0])
readmodel = read_model.ReadModel(self.model, wavel_range=wavel_range)
readmodel.interpolate_grid(self.bounds,
wavel_resample=self.spectrum[key][0][:, 0],
smooth=True,
spec_res=self.spectrum[key][3])
self.modelspec.append(readmodel)
print(f' [DONE]')
else:
self.spectrum = None
self.modelspec = None
self.modelpar = readmodel.get_parameters()
self.modelpar.append('radius')
if self.spectrum is not None:
for item in self.spectrum:
if item in bounds:
self.modelpar.append(f'scaling_{item}')
self.bounds[f'scaling_{item}'] = (bounds[item][0], bounds[item][1])
del self.bounds[item]
def run_mcmc(self,
tag,
guess,
nwalkers=200,
nsteps=1000,
prior=None):
"""
Function to run the MCMC sampler of ``emcee``.
Parameters
----------
tag : str
Database tag where the samples will be stored.
guess : dict, None
Guess for the parameter values. Random values between the boundary values are used
if set to None.
nwalkers : int
Number of walkers.
nsteps : int
Number of steps per walker.
prior : tuple(str, float, float), None
Gaussian prior on one of the parameters. Currently only possible for the mass, e.g.
('mass', 13., 3.) for an expected mass of 13 Mjup with an uncertainty of 3 Mjup. Not
used if set to None.
Returns
-------
NoneType
None
"""
sigma = {'teff': 5., 'logg': 0.01, 'feh': 0.01, 'fsed': 0.01, 'co': 0.01, 'radius': 0.01}
if self.spectrum is not None:
for item in self.spectrum:
if f'scaling_{item}' in self.bounds:
sigma[f'scaling_{item}'] = 0.01
guess[f'scaling_{item}'] = guess[item]
del guess[item]
print('Running MCMC...')
ndim = len(self.bounds)
initial = np.zeros((nwalkers, ndim))
for i, item in enumerate(self.modelpar):
if guess[item] is not None:
initial[:, i] = guess[item] + np.random.normal(0, sigma[item], nwalkers)
else:
initial[:, i] = np.random.uniform(low=self.bounds[item][0],
high=self.bounds[item][1],
size=nwalkers)
with Pool(processes=cpu_count()):
ens_sampler = emcee.EnsembleSampler(nwalkers,
ndim,
lnprob,
args=([self.bounds,
self.modelpar,
self.objphot,
self.distance,
prior,
self.spectrum,
self.modelphot,
self.modelspec]))
ens_sampler.run_mcmc(initial, nsteps, progress=True)
species_db = database.Database()
spec_labels = []
if self.spectrum is not None:
for item in self.spectrum:
if f'scaling_{item}' in self.bounds:
spec_labels.append(f'scaling_{item}')
else:
spec_labels = None
species_db.add_samples(sampler='emcee',
samples=ens_sampler.chain,
ln_prob=ens_sampler.lnprobability,
mean_accept=np.mean(ens_sampler.acceptance_fraction),
spectrum=('model', self.model),
tag=tag,
modelpar=self.modelpar,
distance=self.distance[0],
spec_labels=spec_labels)
def run_multinest(self,
tag,
n_live_points=4000,
output='multinest/'):
"""
Function to run the ``PyMultiNest`` wrapper of the ``MultiNest`` sampler. While
``PyMultiNest`` can be installed with ``pip`` from the PyPI repository, ``MultiNest``
has to to be build manually. See the ``PyMultiNest`` documentation for details:
http://johannesbuchner.github.io/PyMultiNest/install.html. Note that the library path
of ``MultiNest`` should be set to the environmental variable ``LD_LIBRARY_PATH`` on a
Linux machine and ``DYLD_LIBRARY_PATH`` on a Mac. Alternatively, the variable can be
set before importing the ``species`` package, for example:
.. code-block:: python
>>> import os
>>> os.environ['DYLD_LIBRARY_PATH'] = '/path/to/MultiNest/lib'
>>> import species
Parameters
----------
tag : str
Database tag where the samples will be stored.
n_live_points : int
Number of live points.
output : str
Path that is used for the output files from MultiNest.
Returns
-------
NoneType
None
"""
print('Running nested sampling...')
# create the output folder if required
if not os.path.exists(output):
os.mkdir(output)
# create a dictionary with the cube indices of the parameters
cube_index = {}
for i, item in enumerate(self.modelpar):
cube_index[item] = i
def lnprior_multinest(cube, n_dim, n_param):
"""
Function to transform the unit cube into the parameter cube. It is not clear how to
pass additional arguments to the function, therefore it is placed here.
Parameters
----------
cube : pymultinest.run.LP_c_double
Unit cube.
Returns
-------
NoneType
None
"""
# Effective temperature (K)
cube[cube_index['teff']] = self.bounds['teff'][0] + \
(self.bounds['teff'][1]-self.bounds['teff'][0])*cube[cube_index['teff']]
# Surface gravity (dex)
cube[cube_index['logg']] = self.bounds['logg'][0] + \
(self.bounds['logg'][1]-self.bounds['logg'][0])*cube[cube_index['logg']]
# Radius (Rjup)
cube[cube_index['radius']] = self.bounds['radius'][0] + \
(self.bounds['radius'][1]-self.bounds['radius'][0])*cube[cube_index['radius']]
# Metallicity [Fe/H]
if 'feh' in self.bounds:
cube[cube_index['feh']] = self.bounds['feh'][0] + \
(self.bounds['feh'][1]-self.bounds['feh'][0])*cube[cube_index['feh']]
# C/O ratio
if 'co' in self.bounds:
cube[cube_index['co']] = self.bounds['co'][0] + \
(self.bounds['co'][1]-self.bounds['co'][0])*cube[cube_index['co']]
# Sedimentation parameter
if 'fsed' in self.bounds:
cube[cube_index['fsed']] = self.bounds['fsed'][0] + \
(self.bounds['fsed'][1]-self.bounds['fsed'][0])*cube[cube_index['fsed']]
# Spectrum scaling
if self.spectrum is not None:
for item in self.spectrum:
if f'scaling_{item}' in self.bounds:
cube[cube_index[f'scaling_{item}']] = self.bounds[f'scaling_{item}'][0] + \
(self.bounds[f'scaling_{item}'][1]-self.bounds[f'scaling_{item}'][0]) \
* cube[cube_index[f'scaling_{item}']]
def lnlike_multinest(cube, n_dim, n_param):
"""
Function for the logarithm of the likelihood, computed from the parameter cube.
Parameters
----------
cube : pymultinest.run.LP_c_double
Unit cube.
Returns
-------
float
The logarithm of the likelihood.
"""
paramdict = {}
spec_scaling = {}
for i, item in enumerate(self.modelpar):
if item == 'radius':
radius = cube[cube_index['radius']]
elif item[:8] == 'scaling_' and item[8:] in self.spectrum:
spec_scaling[item[8:]] = cube[cube_index[item]]
else:
paramdict[item] = cube[cube_index[item]]
scaling = (radius*constants.R_JUP)**2 / (self.distance[0]*constants.PARSEC)**2
chisq = 0.
if self.objphot is not None:
for i, item in enumerate(self.objphot):
flux = scaling * self.modelphot[i].spectrum_interp(list(paramdict.values()))
chisq += (item[0]-flux)**2 / item[1]**2
if self.spectrum is not None:
for i, item in enumerate(self.spectrum.keys()):
flux = scaling * self.modelspec[i].spectrum_interp(list(paramdict.values()))[0, :]
if item in spec_scaling:
flux_obs = spec_scaling[item]*self.spectrum[item][0][:, 1]
else:
flux_obs = self.spectrum[item][0][:, 1]
if self.spectrum[item][2] is not None:
spec_diff = flux_obs - flux
chisq += np.dot(spec_diff, np.dot(self.spectrum[item][2], spec_diff))
else:
chisq += np.nansum((flux_obs-flux)**2 / self.spectrum[item][0][:, 2]**2)
return -0.5*chisq
pymultinest.run(lnlike_multinest,
lnprior_multinest,
len(self.modelpar),
outputfiles_basename=output,
resume=False,
n_live_points=n_live_points)
samples = np.loadtxt(f'{output}/post_equal_weights.dat')
species_db = database.Database()
spec_labels = []
if self.spectrum is not None:
for item in self.spectrum:
if f'scaling_{item}' in self.bounds:
spec_labels.append(f'scaling_{item}')
else:
spec_labels = None
species_db.add_samples(sampler='multinest',
samples=samples[:, :-1],
ln_prob=samples[:, -1],
mean_accept=None,
spectrum=('model', self.model),
tag=tag,
modelpar=self.modelpar,
distance=self.distance[0],
spec_labels=spec_labels)