-
Notifications
You must be signed in to change notification settings - Fork 5
/
ssp.py
281 lines (245 loc) · 9.35 KB
/
ssp.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
""" MCSED - ssp.py
Single Stellar Population module for loading models
.. moduleauthor:: Greg Zeimann <[email protected]>
"""
import sfh
import sys
import numpy as np
import os.path as op
import scipy.interpolate as scint
from astropy.convolution import Gaussian1DKernel, convolve
def get_coarser_wavelength_fsps(wave, spec):
sel = np.where((wave > 500) * (wave < 1e5))[0]
spec = spec[sel, :]
wave = wave[sel]
G = Gaussian1DKernel(25)
nsel = np.where(np.abs(np.diff(wave)-0.9) < 0.5)[0]
for i in np.arange(spec.shape[1]):
spec[nsel, i] = convolve(spec[nsel, i], G)
ndw = 12.
nw = np.arange(wave[nsel[0]], wave[nsel[-1]+1], ndw)
nwb = np.hstack([nw, wave[nsel[-1]+1]])
nwave = np.hstack([wave[:nsel[0]], nw, wave[(nsel[-1]+1):]])
nspec = np.zeros((len(nwave), spec.shape[1]))
for i, sp in enumerate(spec.swapaxes(0, 1)):
hsp = (np.histogram(wave[nsel], nwb, weights=sp[nsel])[0] /
np.histogram(wave[nsel], nwb)[0])
nspec[:, i] = np.hstack([sp[:nsel[0]], hsp, sp[(nsel[-1]+1):]])
return nwave, nspec
def bin_ages_fsps(args, ages, spec):
sfh_class = getattr(sfh, args.sfh)()
sel = ages >= 6.
ages, spec = (ages[sel], spec[:, sel])
weight = np.diff(np.hstack([0., 10**ages]))
agebin = np.hstack([0., sfh_class.ages])
nspec = np.zeros((spec.shape[0], len(sfh_class.ages)))
for i in np.arange(len(sfh_class.ages)):
sel = np.where((ages >= agebin[i]) * (ages < agebin[i+1]))[0]
nspec[:, i] = np.dot(spec[:, sel], weight[sel]) / weight[sel].sum()
return 10**(np.array(sfh_class.ages)-9.), nspec
def read_fsps_neb(filename):
cnt = 0
Z, Age, logU, spec = [], [], [], []
with open(filename) as f:
for lines in f:
if cnt == 1:
wave = np.array(lines.split(), dtype=float)
if cnt > 1:
l = lines.split()
if len(l) == 3:
Z.append(float(l[0]))
Age.append(float(l[1]))
logU.append(float(l[2]))
else:
spec.append(np.array(l, dtype=float))
cnt += 1
return Z, Age, logU, spec, wave
def read_fsps_ssp(filename):
cnt = 0
ages, masses, spec = [], [], []
with open(filename) as f:
for lines in f:
if cnt == 9:
wave = np.array(lines.split(), dtype=float)
if cnt > 9:
l = lines.split()
if len(l) == 4:
ages.append(float(l[0]))
masses.append(float(l[1]))
else:
spec.append(np.array(l, dtype=float))
cnt += 1
return ages, masses, spec, wave
def read_fsps(args, metallicity):
'''Read in the stellar population models from fsps for a given isochrone
and metallicity.
Parameters
----------
args : class
The args class from mcsed.parse_args()
Returns
-------
ages : numpy array (1 dim)
ages of each SPS (Gyr)
masses : numpy array (1 dim)
Remnant mass at a given age (solar masses)
wave : numpy array (1 dim)
wavelength for each spectrum
spec : numpy array (2 dim)
Spectra in f_nu (ergs/s/cm^2/Hz) at 10pc
'''
pc10 = 10. * 3.08567758e18
solar_microjansky = 3.826e33 * 1e29 / (4. * np.pi * pc10**2)
filename = op.join('SSP', 'fsps_%s_%0.4f.spec' % (args.isochrone,
metallicity))
if not op.exists(filename):
print('Tried to open %s' % filename)
print('Metallicity entered, %0.4f, does not match any of the %s '
'isochrones of the %s models' % (args.metallicity,
args.isochrone, args.ssp))
print('Metallicity options [')
for met in args.metallicity_dict[args.isochrone]:
print('%0.4f ' % met)
print(']')
sys.exit(1)
ages, masses, spec, wave = read_fsps_ssp(filename)
spec = np.array(spec).swapaxes(0, 1) * solar_microjansky
ages, masses = (np.array(ages), np.array(masses))
# Total mass including remnants, so set to 1.
sel = (ages <= 9.5) * (ages >= 6.)
return 10**(ages[sel]-9), np.ones(ages[sel].shape), wave, spec[:, sel]
def add_nebular_emission(ages, wave, spec, logU, metallicity,
filename='nebular/ZAU_ND_pdva',
sollum=3.826e33):
cont_file = filename + '.cont'
lines_file = filename + '.lines'
cont_res = [np.array(x) for x in read_fsps_neb(cont_file)]
lines_res = [np.array(x) for x in read_fsps_neb(lines_file)]
# Make array of Z, age, U
V = np.array([10**cont_res[0]*0.019, cont_res[1]/1e6,
cont_res[2]]).swapaxes(0, 1)
C = scint.LinearNDInterpolator(V, cont_res[3]*1e48)
L = scint.LinearNDInterpolator(V, lines_res[3]*1e48)
garray = make_gaussian_emission(wave, lines_res[4])
nspec = spec * 1.
for i, age in enumerate(ages):
if age <= 1e-2:
cont = C(metallicity, age*1e3, logU)
lines = L(metallicity, age*1e3, logU)
qq = number_ionizing_photons(wave, spec[:, i]) / 1e48 * sollum
nspec[:, i] = (nspec[:, i] + np.interp(wave, cont_res[4], cont*qq)
+ (garray * lines * qq).sum(axis=1))
return nspec
def make_gaussian_emission(wavebig, wave, stddev=1., clight=2.99792e18):
gspec = np.zeros((len(wavebig), len(wave)))
G = Gaussian1DKernel(stddev).array
mid = len(G) / 2
dw = np.diff(wavebig)
for i, w in enumerate(wave):
xl = np.argmin(np.abs(wavebig - w))
if (xl > mid) and ((xl+mid) < len(wavebig)):
gspec[(xl-mid):(xl+mid+1), i] = G / clight * w**2 / dw[xl]
return gspec
def number_ionizing_photons(wave, spectrum, clight=2.99792e18,
hplanck=6.626e-27):
nu = clight / wave
xlim = np.searchsorted(wave, 912., side='right')
x1 = np.abs(np.diff(nu[:xlim]))
x2 = spectrum[:xlim] / nu[:xlim]
x3 = (x2[:-1] + x2[1:]) / 2.
return np.sum(x1 * x3) / hplanck
def read_ssp(args):
''' Read in SPS model and return ages, masses, wavelength, and spectra
Parameters
----------
args : class
The args class from mcsed.parse_args()
'''
import matplotlib.pyplot as plt
s, m = ([], [])
for met in args.metallicity_dict[args.isochrone]:
if args.ssp.lower() == 'fsps':
ages, masses, wave, spec = read_fsps(args, met)
if args.add_nebular:
spec = add_nebular_emission(ages, wave, spec, args.logU,
met)
if args.sfh == 'empirical' or args.sfh == 'empirical_direct':
ages, spec = bin_ages_fsps(args, np.log10(ages)+9., spec)
masses = np.ones(ages.shape)
wave, spec = get_coarser_wavelength_fsps(wave, spec)
wei = (np.diff(np.hstack([0., ages])) *
getattr(sfh, args.sfh)().evaluate(ages))
s.append(spec)
m.append(np.dot(spec, wei)/spec.shape[1])
fig = plt.figure(figsize=(8, 8))
import seaborn as sns
colors = sns.color_palette("coolwarm", s[-5].shape[1])
wei = np.diff(np.hstack([0., ages]))
for i in np.arange(s[-5].shape[1]):
plt.plot(wave, s[-5][:, i] * wei[i] / 1e8, color=colors[i])
plt.xlim([900., 40000.])
plt.xscale('log')
plt.yscale('log')
plt.ylim([1e-5, 10])
plt.savefig('template_spectra_plot.png')
plt.close(fig)
spec = np.moveaxis(np.array(s), 0, 2)
metallicities = args.metallicity_dict[args.isochrone]
return ages, masses, wave, spec, np.array(metallicities)
class fsps_freeparams:
''' Allowing metallicity to be free -0.39'''
def __init__(self, met=-0.7, met_lims=[-1.98, 0.2], met_delta=0.3,
fix_met=False):
''' Initialize this class
Parameters
----------
TODO Fill these in
'''
self.met = met
if fix_met:
self.nparams = 0
else:
self.nparams = 1
self.met_lims = met_lims
self.met_delta = met_delta
self.fix_met = fix_met
def get_params(self):
''' Return current parameters '''
l = []
if not self.fix_met:
l.append(self.met)
return l
def get_param_lims(self):
''' Return current parameters limits '''
l = []
if not self.fix_met:
l.append(self.met_lims)
return l
def get_param_deltas(self):
''' Return current parameter deltas '''
l = []
if not self.fix_met:
l.append(self.met_delta)
return l
def get_names(self):
''' Return names of each parameter '''
l = []
if not self.fix_met:
l.append('Log Z')
return l
def prior(self):
''' Uniform prior based on boundaries '''
flag = (self.met >= self.met_lims[0]) * (self.met <= self.met_lims[1])
return flag
def set_parameters_from_list(self, input_list, start_value):
''' Set parameters from a list and a start_value
Parameters
----------
input_list : list
list of input parameters (could be much larger than number of
parameters to be set)
start_value : int
initial index from list to read out parameters
'''
if not self.fix_met:
self.met = input_list[start_value]