-
Notifications
You must be signed in to change notification settings - Fork 10
/
phot_util.py
289 lines (224 loc) · 10.1 KB
/
phot_util.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
"""
Utility functions for photometry.
"""
import math
import warnings
from typing import Optional, Union, Dict, List
import spectres
import numpy as np
from typeguard import typechecked
from species.core import box
from species.read import read_model, read_calibration, read_filter, read_planck
# from species.read import read_model, read_calibration, read_filter, read_planck, read_radtrans
def multi_photometry(datatype,
spectrum,
filters,
parameters):
"""
Parameters
----------
datatype : str
Data type ('model' or 'calibration').
spectrum : str
Spectrum name (e.g., 'drift-phoenix').
filters : tuple(str, )
Filter IDs.
parameters : dict
Parameters and values for the spectrum
Returns
-------
species.core.box.SynphotBox
Box with synthetic photometry.
"""
print('Calculating synthetic photometry...', end='', flush=True)
flux = {}
if datatype == 'model':
for item in filters:
if spectrum == 'planck':
readmodel = read_planck.ReadPlanck(filter_name=item)
else:
readmodel = read_model.ReadModel(spectrum, filter_name=item)
try:
flux[item] = readmodel.get_flux(parameters)[0]
except IndexError:
flux[item] = np.nan
warnings.warn(f'The wavelength range of the {item} filter does not match with '
f'the wavelength coverage of {spectrum}. The flux is set to NaN.')
elif datatype == 'calibration':
for item in filters:
readcalib = read_calibration.ReadCalibration(spectrum, filter_name=item)
flux[item] = readcalib.get_flux(parameters)[0]
print(' [DONE]')
return box.create_box('synphot', name='synphot', flux=flux)
def apparent_to_absolute(app_mag,
distance):
"""
Function for converting an apparent magnitude into an absolute magnitude. The uncertainty on
the distance is propagated into the uncertainty on the absolute magnitude.
Parameters
----------
app_mag : tuple(float, float), tuple(numpy.ndarray, numpy.ndarray)
Apparent magnitude and uncertainty (mag). The returned error on the absolute magnitude
is set to None if the error on the apparent magnitude is set to None, for example
``app_mag=(15., None)``.
distance : tuple(float, float), tuple(numpy.ndarray, numpy.ndarray)
Distance and uncertainty (pc). The error is not propagated into the error on the absolute
magnitude if set to None, for example ``distance=(20., None)``.
Returns
-------
float, numpy.ndarray
Absolute magnitude (mag).
float, numpy.ndarray, None
Uncertainty (mag).
"""
abs_mag = app_mag[0] - 5.*np.log10(distance[0]) + 5.
if app_mag[1] is not None and distance[1] is not None:
dist_err = distance[1] * (5./(distance[0]*math.log(10.)))
abs_err = math.sqrt(app_mag[1]**2 + dist_err**2)
elif app_mag[1] is not None and distance[1] is None:
abs_err = app_mag[1]
else:
abs_err = None
return abs_mag, abs_err
def absolute_to_apparent(abs_mag,
distance):
"""
Function for converting an absolute magnitude into an apparent magnitude.
Parameters
----------
abs_mag : tuple(float, float), tuple(numpy.ndarray, numpy.ndarray)
Absolute magnitude and uncertainty (mag). The same uncertainty is used for the
apparent magnitude.
distance : tuple(float, float), tuple(numpy.ndarray, numpy.ndarray)
Distance and uncertainty (pc).
Returns
-------
float, numpy.ndarray
Apparent magnitude (mag).
float, numpy.ndarray, None
Uncertainty (mag).
"""
app_mag = abs_mag[0] + 5.*np.log10(distance[0]) - 5.
return app_mag, abs_mag[1]
@typechecked
def get_residuals(datatype: str,
spectrum: str,
parameters: Dict[str, float],
objectbox: box.ObjectBox,
inc_phot: Union[bool, List[str]] = True,
inc_spec: Union[bool, List[str]] = True,
**kwargs_radtrans: Optional[dict]) -> box.ResidualsBox:
"""
Parameters
----------
datatype : str
Data type ('model' or 'calibration').
spectrum : str
Name of the atmospheric model or calibration spectrum.
parameters : dict
Parameters and values for the spectrum
objectbox : species.core.box.ObjectBox
Box with the photometry and/or spectra of an object. A scaling and/or error inflation of
the spectra should be applied with :func:`~species.util.read_util.update_spectra`
beforehand.
inc_phot : bool, list(str)
Include photometric data in the fit. If a boolean, either all (``True``) or none
(``False``) of the data are selected. If a list, a subset of filter names (as stored in
the database) can be provided.
inc_spec : bool, list(str)
Include spectroscopic data in the fit. If a boolean, either all (``True``) or none
(``False``) of the data are selected. If a list, a subset of spectrum names (as stored
in the database with :func:`~species.data.database.Database.add_object`) can be
provided.
Keyword arguments
-----------------
kwargs_radtrans : dict
Dictionary with the keyword arguments for the ``ReadRadtrans`` object, containing
``line_species``, ``cloud_species``, and ``scattering``.
Returns
-------
species.core.box.ResidualsBox
Box with the residuals.
"""
if 'filters' in kwargs_radtrans:
warnings.warn('The \'filters\' parameter has been deprecated. Please use the \'inc_phot\' '
'parameter instead. The \'filters\' parameter is ignored.')
if isinstance(inc_phot, bool) and inc_phot:
inc_phot = objectbox.filters
if inc_phot:
model_phot = multi_photometry(datatype=datatype,
spectrum=spectrum,
filters=inc_phot,
parameters=parameters)
res_phot = {}
for item in inc_phot:
transmission = read_filter.ReadFilter(item)
res_phot[item] = np.zeros(objectbox.flux[item].shape)
if objectbox.flux[item].ndim == 1:
res_phot[item][0] = transmission.mean_wavelength()
res_phot[item][1] = (objectbox.flux[item][0]-model_phot.flux[item]) / \
objectbox.flux[item][1]
elif objectbox.flux[item].ndim == 2:
for j in range(objectbox.flux[item].shape[1]):
res_phot[item][0, j] = transmission.mean_wavelength()
res_phot[item][1, j] = (objectbox.flux[item][0, j]-model_phot.flux[item]) / \
objectbox.flux[item][1, j]
else:
res_phot = None
if inc_spec:
res_spec = {}
readmodel = None
for key in objectbox.spectrum:
if isinstance(inc_spec, bool) or key in inc_spec:
wavel_range = (0.9*objectbox.spectrum[key][0][0, 0],
1.1*objectbox.spectrum[key][0][-1, 0])
wl_new = objectbox.spectrum[key][0][:, 0]
spec_res = objectbox.spectrum[key][3]
if spectrum == 'planck':
readmodel = read_planck.ReadPlanck(wavel_range=wavel_range)
model = readmodel.get_spectrum(model_param=parameters, spec_res=1000.)
flux_new = spectres.spectres(wl_new, model.wavelength, model.flux)
else:
if spectrum == 'petitradtrans':
# TODO change back
pass
# radtrans = read_radtrans.ReadRadtrans(line_species=kwargs_radtrans['line_species'],
# cloud_species=kwargs_radtrans['cloud_species'],
# scattering=kwargs_radtrans['scattering'],
# wavel_range=wavel_range)
#
# model = radtrans.get_model(parameters, spec_res=None)
#
# # separate resampling to the new wavelength points
#
# flux_new = spectres.spectres(wl_new, model.wavelength, model.flux)
else:
readmodel = read_model.ReadModel(spectrum, wavel_range=wavel_range)
# resampling to the new wavelength points is done in teh get_model function
model = readmodel.get_model(parameters,
spec_res=spec_res,
wavel_resample=wl_new,
smooth=True)
flux_new = model.flux
res_tmp = (objectbox.spectrum[key][0][:, 1]-flux_new)/objectbox.spectrum[key][0][:, 2]
res_spec[key] = np.column_stack([wl_new, res_tmp])
else:
res_spec = None
print('Calculating residuals... [DONE]')
print('Residuals (sigma):')
if res_phot is not None:
for item in inc_phot:
if res_phot[item].ndim == 1:
print(f' - {item}: {res_phot[item][1]:.2f}')
elif res_phot[item].ndim == 2:
for j in range(res_phot[item].shape[1]):
print(f' - {item}: {res_phot[item][1, j]:.2f}')
if res_spec is not None:
for key in objectbox.spectrum:
if isinstance(inc_spec, bool) or key in inc_spec:
print(f' - {key}: min: {np.amin(res_spec[key]):.2f}, '
f'max: {np.amax(res_spec[key]):.2f}')
return box.create_box(boxtype='residuals',
name=objectbox.name,
photometry=res_phot,
spectrum=res_spec)