1
|
|
|
# -*- coding: utf-8 -*- |
2
|
|
|
# vim:fileencoding=utf-8 |
3
|
|
|
# |
4
|
|
|
# Copyright (c) 2017-2018 Stefan Bender |
5
|
|
|
# |
6
|
|
|
# This module is part of sciapy. |
7
|
|
|
# sciapy is free software: you can redistribute it or modify |
8
|
|
|
# it under the terms of the GNU General Public License as published |
9
|
|
|
# by the Free Software Foundation, version 2. |
10
|
|
|
# See accompanying LICENSE file or http://www.gnu.org/licenses/gpl-2.0.html. |
11
|
1 |
|
"""SCIAMACHY regression models (celerite version) |
12
|
|
|
|
13
|
|
|
Model classes for SCIAMACHY data regression fits using the |
14
|
|
|
:mod:`celerite` [#]_ modeling protocol. |
15
|
|
|
|
16
|
|
|
.. [#] https://celerite.readthedocs.io |
17
|
|
|
""" |
18
|
1 |
|
from __future__ import absolute_import, division, print_function |
19
|
|
|
|
20
|
1 |
|
import numpy as np |
21
|
1 |
|
from scipy.interpolate import interp1d |
22
|
|
|
|
23
|
1 |
|
from celerite.modeling import Model, ModelSet, ConstantModel |
24
|
|
|
|
25
|
1 |
|
__all__ = ["ConstantModel", |
26
|
|
|
"HarmonicModelCosineSine", "HarmonicModelAmpPhase", |
27
|
|
|
"ProxyModel", "TraceGasModelSet", |
28
|
|
|
"setup_proxy_model_with_bounds", "trace_gas_model"] |
29
|
|
|
|
30
|
|
|
|
31
|
1 |
View Code Duplication |
class HarmonicModelCosineSine(Model): |
|
|
|
|
32
|
|
|
"""Model for harmonic terms |
33
|
|
|
|
34
|
|
|
Models harmonic terms using a cosine and sine part. |
35
|
|
|
The total amplitude and phase can be inferred from that. |
36
|
|
|
|
37
|
|
|
Parameters |
38
|
|
|
---------- |
39
|
|
|
freq : float |
40
|
|
|
The frequency in years^-1 |
41
|
|
|
cos : float |
42
|
|
|
The amplitude of the cosine part |
43
|
|
|
sin : float |
44
|
|
|
The amplitude of the sine part |
45
|
|
|
""" |
46
|
1 |
|
parameter_names = ("freq", "cos", "sin") |
47
|
|
|
|
48
|
1 |
|
def get_value(self, t): |
49
|
1 |
|
t = np.atleast_1d(t) |
50
|
1 |
|
return (self.cos * np.cos(self.freq * 2 * np.pi * t) + |
51
|
|
|
self.sin * np.sin(self.freq * 2 * np.pi * t)) |
52
|
|
|
|
53
|
1 |
|
def get_amplitude(self): |
54
|
|
|
return np.sqrt(self.cos**2 + self.sin**2) |
55
|
|
|
|
56
|
1 |
|
def get_phase(self): |
57
|
|
|
return np.arctan2(self.sin, self.cos) |
58
|
|
|
|
59
|
1 |
|
def compute_gradient(self, t): |
60
|
1 |
|
t = np.atleast_1d(t) |
61
|
1 |
|
dcos = np.cos(self.freq * 2 * np.pi * t) |
62
|
1 |
|
dsin = np.sin(self.freq * 2 * np.pi * t) |
63
|
1 |
|
df = 2 * np.pi * t * (self.sin * dcos - self.cos * dsin) |
64
|
1 |
|
return np.array([df, dcos, dsin]) |
65
|
|
|
|
66
|
|
|
|
67
|
1 |
View Code Duplication |
class HarmonicModelAmpPhase(Model): |
|
|
|
|
68
|
|
|
"""Model for harmonic terms |
69
|
|
|
|
70
|
|
|
Models harmonic terms using a cosine and sine part. |
71
|
|
|
The total amplitude and phase can be inferred from that. |
72
|
|
|
|
73
|
|
|
Parameters |
74
|
|
|
---------- |
75
|
|
|
freq : float |
76
|
|
|
The frequency in years^-1 |
77
|
|
|
amp : float |
78
|
|
|
The amplitude of the harmonic term |
79
|
|
|
phase : float |
80
|
|
|
The phase of the harmonic part |
81
|
|
|
""" |
82
|
1 |
|
parameter_names = ("freq", "amp", "phase") |
83
|
|
|
|
84
|
1 |
|
def get_value(self, t): |
85
|
|
|
t = np.atleast_1d(t) |
86
|
|
|
return self.amp * np.cos(self.freq * 2 * np.pi * t + self.phase) |
87
|
|
|
|
88
|
1 |
|
def get_amplitude(self): |
89
|
|
|
return self.amp |
90
|
|
|
|
91
|
1 |
|
def get_phase(self): |
92
|
|
|
return self.phase |
93
|
|
|
|
94
|
1 |
|
def compute_gradient(self, t): |
95
|
|
|
t = np.atleast_1d(t) |
96
|
|
|
damp = np.cos(self.freq * 2 * np.pi * t + self.phase) |
97
|
|
|
dphi = -self.amp * np.sin(self.freq * 2 * np.pi * t + self.phase) |
98
|
|
|
df = 2 * np.pi * t * dphi |
99
|
|
|
return np.array([df, damp, dphi]) |
100
|
|
|
|
101
|
|
|
|
102
|
1 |
View Code Duplication |
class ProxyModel(Model): |
|
|
|
|
103
|
|
|
"""Model for proxy terms |
104
|
|
|
|
105
|
|
|
Models proxy terms with a finite and (semi-)annually varying life time. |
106
|
|
|
|
107
|
|
|
Parameters |
108
|
|
|
---------- |
109
|
|
|
proxy_times : (N,) array_like |
110
|
|
|
The data times of the proxy values |
111
|
|
|
proxy_vals : (N,) array_like |
112
|
|
|
The proxy values at `proxy_times` |
113
|
|
|
amp : float |
114
|
|
|
The amplitude of the proxy term |
115
|
|
|
lag : float |
116
|
|
|
The lag of the proxy value in years. |
117
|
|
|
tau0 : float |
118
|
|
|
The base life time of the proxy |
119
|
|
|
taucos1 : float |
120
|
|
|
The amplitude of the cosine part of the annual life time variation. |
121
|
|
|
tausin1 : float |
122
|
|
|
The amplitude of the sine part of the annual life time variation. |
123
|
|
|
taucos2 : float |
124
|
|
|
The amplitude of the cosine part of the semi-annual life time variation. |
125
|
|
|
tausin2 : float |
126
|
|
|
The amplitude of the sine part of the semi-annual life time variation. |
127
|
|
|
ltscan : float |
128
|
|
|
The number of days to sum the previous proxy values. If it is |
129
|
|
|
negative, the value will be set to three times the maximal lifetime. |
130
|
|
|
No lifetime adjustemets are calculated when set to zero. |
131
|
|
|
center : bool, optional |
132
|
|
|
Centers the proxy values by subtracting the overall mean. The mean is |
133
|
|
|
calculated from the whole `proxy_vals` array and is stored in the |
134
|
|
|
`mean` attribute. |
135
|
|
|
Default: False |
136
|
|
|
sza_intp : scipy.interpolate.interp1d() instance, optional |
137
|
|
|
When not `None`, cos(sza) and sin(sza) are used instead |
138
|
|
|
of the time to model the annual variation of the lifetime. |
139
|
|
|
Semi-annual variations are not used in that case. |
140
|
|
|
Default: None |
141
|
|
|
fit_phase : bool, optional |
142
|
|
|
Fit the phase shift directly instead of using sine and cosine |
143
|
|
|
terms for the (semi-)annual lifetime variations. If True, the fitted |
144
|
|
|
cosine parameter is the amplitude and the sine parameter the phase. |
145
|
|
|
Default: False (= fit sine and cosine terms) |
146
|
|
|
lifetime_prior : str, optional |
147
|
|
|
The prior probability density for each coefficient of the lifetime. |
148
|
|
|
Possible types are "flat" or `None` for a flat prior, "exp" for an |
149
|
|
|
exponential density ~ :math:`\\text{exp}(-|\\tau| / \\text{metric})`, |
150
|
|
|
and "normal" for a normal distribution |
151
|
|
|
~ :math:`\\text{exp}(-\\tau^2 / (2 * \\text{metric}^2))`. |
152
|
|
|
Default: None (= flat prior). |
153
|
|
|
lifetime_metric : float, optional |
154
|
|
|
The metric (scale) of the lifetime priors in days, see `prior`. |
155
|
|
|
Default 1. |
156
|
|
|
days_per_time_unit : float, optional |
157
|
|
|
The number of days per time unit, used to normalize the lifetime |
158
|
|
|
units. Use 365.25 if the times are in fractional years, or 1 if |
159
|
|
|
they are in days. |
160
|
|
|
Default: 365.25 |
161
|
|
|
""" |
162
|
1 |
|
parameter_names = ("amp", "lag", "tau0", |
163
|
|
|
"taucos1", "tausin1", "taucos2", "tausin2", |
164
|
|
|
"ltscan") |
165
|
|
|
|
166
|
1 |
|
def __init__(self, proxy_times, proxy_vals, |
167
|
|
|
center=False, |
168
|
|
|
sza_intp=None, fit_phase=False, |
169
|
|
|
lifetime_prior=None, lifetime_metric=1., |
170
|
|
|
days_per_time_unit=365.25, |
171
|
|
|
*args, **kwargs): |
172
|
1 |
|
self.mean = 0. |
173
|
1 |
|
if center: |
174
|
|
|
self.mean = np.nanmean(proxy_vals) |
175
|
1 |
|
self.intp = interp1d(proxy_times, proxy_vals - self.mean, |
176
|
|
|
bounds_error=False) |
177
|
1 |
|
self.sza_intp = sza_intp |
178
|
1 |
|
self.fit_phase = fit_phase |
179
|
1 |
|
self.days_per_time_unit = days_per_time_unit |
180
|
1 |
|
self.omega = 2 * np.pi * days_per_time_unit / 365.25 |
181
|
1 |
|
self.lifetime_prior = lifetime_prior |
182
|
1 |
|
self.lifetime_metric = lifetime_metric |
183
|
1 |
|
super(ProxyModel, self).__init__(*args, **kwargs) |
184
|
|
|
|
185
|
1 |
|
def get_value(self, t): |
186
|
1 |
|
t = np.atleast_1d(t) |
187
|
1 |
|
proxy_val = self.intp(t - self.lag) |
188
|
1 |
|
if self.ltscan == 0: |
189
|
|
|
# no lifetime, nothing else to do |
190
|
1 |
|
return self.amp * proxy_val |
191
|
|
|
# annual variation of the proxy lifetime |
192
|
1 |
|
if self.sza_intp is not None: |
193
|
|
|
# using the solar zenith angle |
194
|
|
|
tau_cs = (self.taucos1 * np.cos(np.radians(self.sza_intp(t))) |
195
|
|
|
+ self.tausin1 * np.sin(np.radians(self.sza_intp(t)))) |
196
|
1 |
|
elif self.fit_phase: |
197
|
|
|
# using time (cos) and phase (sin) |
198
|
|
|
tau_cs = (self.taucos1 * np.cos(1 * self.omega * t + self.tausin1) |
199
|
|
|
+ self.taucos2 * np.cos(2 * self.omega * t + self.tausin2)) |
200
|
|
|
else: |
201
|
|
|
# using time |
202
|
1 |
|
tau_cs = (self.taucos1 * np.cos(1 * self.omega * t) |
203
|
|
|
+ self.tausin1 * np.sin(1 * self.omega * t) |
204
|
|
|
+ self.taucos2 * np.cos(2 * self.omega * t) |
205
|
|
|
+ self.tausin2 * np.sin(2 * self.omega * t)) |
206
|
1 |
|
tau_cs[tau_cs < 0] = 0. # clip to zero |
207
|
1 |
|
tau = self.tau0 + tau_cs |
208
|
1 |
|
if self.ltscan > 0: |
209
|
1 |
|
_ltscn = int(np.floor(self.ltscan)) |
210
|
|
|
else: |
211
|
|
|
# infer the scan time from the maximal lifetime |
212
|
|
|
_ltscn = 3 * int(np.ceil(self.tau0 + |
213
|
|
|
np.sqrt(self.taucos1**2 + self.tausin1**2))) |
214
|
1 |
|
if np.all(tau > 0): |
215
|
1 |
|
bs = np.arange(1, _ltscn + 1, 1.)[None, :] |
216
|
1 |
|
taufacs = np.exp(-bs / tau[:, None]) |
217
|
1 |
|
proxy_val += np.sum( |
218
|
|
|
self.intp(t[:, None] - self.lag - |
219
|
|
|
bs / self.days_per_time_unit) * taufacs, |
220
|
|
|
axis=1) |
221
|
1 |
|
return self.amp * proxy_val |
222
|
|
|
|
223
|
1 |
|
def compute_gradient(self, t): |
224
|
1 |
|
t = np.atleast_1d(t) |
225
|
1 |
|
proxy_val = self.intp(t - self.lag) |
226
|
1 |
|
proxy_val_grad0 = self.intp(t - self.lag) |
227
|
|
|
# annual variation of the proxy lifetime |
228
|
1 |
|
if self.sza_intp is not None: |
229
|
|
|
# using the solar zenith angle |
230
|
|
|
dtau_cos1 = np.cos(np.radians(self.sza_intp(t))) |
231
|
|
|
dtau_sin1 = np.sin(np.radians(self.sza_intp(t))) |
232
|
|
|
dtau_cos2 = np.zeros_like(t) |
233
|
|
|
dtau_sin2 = np.zeros_like(t) |
234
|
|
|
tau_cs = self.taucos1 * dtau_cos1 + self.tausin1 * dtau_sin1 |
235
|
1 |
|
elif self.fit_phase: |
236
|
|
|
# using time (cos) and phase (sin) |
237
|
|
|
dtau_cos1 = np.cos(1 * self.omega * t + self.tausin1) |
238
|
|
|
dtau_sin1 = -self.taucos1 * np.sin(1 * self.omega * t + self.tausin1) |
239
|
|
|
dtau_cos2 = np.cos(2 * self.omega * t + self.tausin2) |
240
|
|
|
dtau_sin2 = -self.taucos2 * np.sin(2 * self.omega * t + self.tausin2) |
241
|
|
|
tau_cs = self.taucos1 * dtau_cos1 + self.taucos2 * dtau_cos2 |
242
|
|
|
else: |
243
|
|
|
# using time |
244
|
1 |
|
dtau_cos1 = np.cos(1 * self.omega * t) |
245
|
1 |
|
dtau_sin1 = np.sin(1 * self.omega * t) |
246
|
1 |
|
dtau_cos2 = np.cos(2 * self.omega * t) |
247
|
1 |
|
dtau_sin2 = np.sin(2 * self.omega * t) |
248
|
1 |
|
tau_cs = (self.taucos1 * dtau_cos1 + self.tausin1 * dtau_sin1 + |
249
|
|
|
self.taucos2 * dtau_cos2 + self.tausin2 * dtau_sin2) |
250
|
1 |
|
tau_cs[tau_cs < 0] = 0. # clip to zero |
251
|
1 |
|
tau = self.tau0 + tau_cs |
252
|
1 |
|
if self.ltscan > 0: |
253
|
1 |
|
_ltscn = int(np.floor(self.ltscan)) |
254
|
|
|
else: |
255
|
|
|
# infer the scan time from the maximal lifetime |
256
|
1 |
|
_ltscn = 3 * int(np.ceil(self.tau0 + |
257
|
|
|
np.sqrt(self.taucos1**2 + self.tausin1**2))) |
258
|
1 |
|
if np.all(tau > 0): |
259
|
1 |
|
bs = np.arange(1, _ltscn + 1, 1.)[None, :] |
260
|
1 |
|
taufacs = np.exp(-bs / tau[:, None]) |
261
|
1 |
|
proxy_ts = self.intp(t[:, None] - self.lag - |
262
|
|
|
bs / self.days_per_time_unit) * taufacs |
263
|
1 |
|
proxy_val += np.sum(proxy_ts, axis=1) |
264
|
1 |
|
proxy_val_grad0 += np.sum(proxy_ts * bs / tau[:, None]**2, axis=1) |
265
|
1 |
|
return np.array([proxy_val, |
266
|
|
|
# set the gradient wrt lag to zero for now |
267
|
|
|
np.zeros_like(t), |
268
|
|
|
self.amp * proxy_val_grad0, |
269
|
|
|
self.amp * proxy_val_grad0 * dtau_cos1, |
270
|
|
|
self.amp * proxy_val_grad0 * dtau_sin1, |
271
|
|
|
self.amp * proxy_val_grad0 * dtau_cos2, |
272
|
|
|
self.amp * proxy_val_grad0 * dtau_sin2, |
273
|
|
|
# set the gradient wrt lifetime scan to zero for now |
274
|
|
|
np.zeros_like(t)]) |
275
|
|
|
|
276
|
1 |
|
def _log_prior_normal(self): |
277
|
|
|
l_prior = super(ProxyModel, self).log_prior() |
278
|
|
|
if not np.isfinite(l_prior): |
279
|
|
|
return -np.inf |
280
|
|
|
for n, p in self.get_parameter_dict().items(): |
281
|
|
|
if n.startswith("tau"): |
282
|
|
|
# Gaussian prior for the lifetimes |
283
|
|
|
l_prior -= 0.5 * (p / self.lifetime_metric)**2 |
284
|
|
|
return l_prior |
285
|
|
|
|
286
|
1 |
|
def _log_prior_exp(self): |
287
|
1 |
|
l_prior = super(ProxyModel, self).log_prior() |
288
|
1 |
|
if not np.isfinite(l_prior): |
289
|
|
|
return -np.inf |
290
|
1 |
|
for n, p in self.get_parameter_dict().items(): |
291
|
1 |
|
if n.startswith("tau"): |
292
|
|
|
# exponential prior for the lifetimes |
293
|
1 |
|
l_prior -= np.abs(p / self.lifetime_metric) |
294
|
1 |
|
return l_prior |
295
|
|
|
|
296
|
1 |
|
def log_prior(self): |
297
|
1 |
|
_priors = {"exp": self._log_prior_exp, |
298
|
|
|
"normal": self._log_prior_normal} |
299
|
1 |
|
if self.lifetime_prior is None or self.lifetime_prior == "flat": |
300
|
1 |
|
return super(ProxyModel, self).log_prior() |
301
|
1 |
|
return _priors[self.lifetime_prior]() |
302
|
|
|
|
303
|
|
|
|
304
|
1 |
View Code Duplication |
class TraceGasModelSet(ModelSet): |
|
|
|
|
305
|
|
|
"""Combined model class for trace gases (and probably other data) |
306
|
|
|
|
307
|
|
|
Inherited from :class:`celerite.ModelSet`, provides `get_value()` |
308
|
|
|
and `compute_gradient()` methods. |
309
|
|
|
""" |
310
|
1 |
|
def get_value(self, t): |
311
|
1 |
|
v = np.zeros_like(t) |
312
|
1 |
|
for m in self.models.values(): |
313
|
1 |
|
v += m.get_value(t) |
314
|
1 |
|
return v |
315
|
|
|
|
316
|
1 |
|
def compute_gradient(self, t): |
317
|
1 |
|
grad = [] |
318
|
1 |
|
for m in self.models.values(): |
319
|
1 |
|
grad.extend(list(m.compute_gradient(t))) |
320
|
1 |
|
return np.array(grad) |
321
|
|
|
|
322
|
|
|
|
323
|
1 |
View Code Duplication |
def setup_proxy_model_with_bounds(times, values, |
|
|
|
|
324
|
|
|
max_amp=1e10, max_days=100, |
325
|
|
|
**kwargs): |
326
|
|
|
# extract setup from `kwargs` |
327
|
1 |
|
center = kwargs.get("center", False) |
328
|
1 |
|
fit_phase = kwargs.get("fit_phase", False) |
329
|
1 |
|
lag = kwargs.get("lag", 0.) |
330
|
1 |
|
lt_metric = kwargs.get("lifetime_metric", 1) |
331
|
1 |
|
lt_prior = kwargs.get("lifetime_prior", "exp") |
332
|
1 |
|
lt_scan = kwargs.get("lifetime_scan", 60) |
333
|
1 |
|
positive = kwargs.get("positive", False) |
334
|
1 |
|
sza_intp = kwargs.get("sza_intp", None) |
335
|
1 |
|
time_format = kwargs.get("time_format", "jyear") |
336
|
|
|
|
337
|
1 |
|
return ProxyModel(times, values, |
338
|
|
|
center=center, |
339
|
|
|
sza_intp=sza_intp, |
340
|
|
|
fit_phase=fit_phase, |
341
|
|
|
lifetime_prior=lt_prior, |
342
|
|
|
lifetime_metric=lt_metric, |
343
|
|
|
days_per_time_unit=1 if time_format.endswith("d") else 365.25, |
344
|
|
|
amp=0., |
345
|
|
|
lag=lag, |
346
|
|
|
tau0=0, |
347
|
|
|
taucos1=0, tausin1=0, |
348
|
|
|
taucos2=0, tausin2=0, |
349
|
|
|
ltscan=lt_scan, |
350
|
|
|
bounds=dict([ |
351
|
|
|
("amp", [0, max_amp] if positive else [-max_amp, max_amp]), |
352
|
|
|
("lag", [0, max_days]), |
353
|
|
|
("tau0", [0, max_days]), |
354
|
|
|
("taucos1", [0, max_days] if fit_phase else [-max_days, max_days]), |
355
|
|
|
("tausin1", [-np.pi, np.pi] if fit_phase else [-max_days, max_days]), |
356
|
|
|
# semi-annual cycles for the life time |
357
|
|
|
("taucos2", [0, max_days] if fit_phase else [-max_days, max_days]), |
358
|
|
|
("tausin2", [-np.pi, np.pi] if fit_phase else [-max_days, max_days]), |
359
|
|
|
("ltscan", [0, 200])]) |
360
|
|
|
) |
361
|
|
|
|
362
|
|
|
|
363
|
1 |
View Code Duplication |
def _default_proxy_config(tfmt="jyear"): |
|
|
|
|
364
|
|
|
from .load_data import load_dailymeanLya, load_dailymeanAE |
365
|
|
|
proxy_config = {} |
366
|
|
|
# Lyman-alpha |
367
|
|
|
plyat, plyadf = load_dailymeanLya(tfmt=tfmt) |
368
|
|
|
proxy_config.update({"Lya": { |
369
|
|
|
"times": plyat, |
370
|
|
|
"values": plyadf["Lya"], |
371
|
|
|
"center": False, |
372
|
|
|
"positive": False, |
373
|
|
|
"lifetime_scan": 0, |
374
|
|
|
}} |
375
|
|
|
) |
376
|
|
|
# AE index |
377
|
|
|
paet, paedf = load_dailymeanAE(name="GM", tfmt=tfmt) |
378
|
|
|
proxy_config.update({"GM": { |
379
|
|
|
"times": paet, |
380
|
|
|
"values": paedf["GM"], |
381
|
|
|
"center": False, |
382
|
|
|
"positive": True, |
383
|
|
|
"lifetime_scan": 60, |
384
|
|
|
}} |
385
|
|
|
) |
386
|
|
|
return proxy_config |
387
|
|
|
|
388
|
|
|
|
389
|
1 |
View Code Duplication |
def trace_gas_model(constant=True, freqs=None, proxy_config=None, **kwargs): |
|
|
|
|
390
|
|
|
"""Trace gas model setup |
391
|
|
|
|
392
|
|
|
Sets up the trace gas model for easy access. All parameters are optional, |
393
|
|
|
defaults to an offset, no harmonics, proxies are uncentered and unscaled |
394
|
|
|
Lyman-alpha and AE. AE with only positive amplitude and a seasonally |
395
|
|
|
varying lifetime. |
396
|
|
|
|
397
|
|
|
Parameters |
398
|
|
|
---------- |
399
|
|
|
constant : bool, optional |
400
|
|
|
Whether or not to include a constant (offset) term, default is True. |
401
|
|
|
freqs : list, optional |
402
|
|
|
Frequencies of the harmonic terms in 1 / a^-1 (inverse years). |
403
|
|
|
proxy_config : dict, optional |
404
|
|
|
Proxy configuration if different from the standard setup. |
405
|
|
|
**kwargs : optional |
406
|
|
|
Additional keyword arguments, all of them are also passed on to |
407
|
|
|
the proxy setup. For now, supported are the following which are |
408
|
|
|
also passed along to the proxy setup with |
409
|
|
|
`setup_proxy_model_with_bounds()`: |
410
|
|
|
|
411
|
|
|
* fit_phase : bool |
412
|
|
|
fit amplitude and phase instead of sine and cosine |
413
|
|
|
* scale : float |
414
|
|
|
the factor by which the data is scaled, used to constrain |
415
|
|
|
the maximum and minimum amplitudes to be fitted. |
416
|
|
|
* tfmt : string |
417
|
|
|
The `astropy.time.Time` format string to setup the time axis. |
418
|
|
|
* max_amp : float |
419
|
|
|
Maximum magnitude of the coefficients, used to constrain the |
420
|
|
|
parameter search. |
421
|
|
|
* max_days : float |
422
|
|
|
Maximum magnitude of the lifetimes, used to constrain the |
423
|
|
|
parameter search. |
424
|
|
|
|
425
|
|
|
Returns |
426
|
|
|
------- |
427
|
|
|
model : :class:`TraceGasModelSet` (extended :class:`celerite.ModelSet`) |
428
|
|
|
""" |
429
|
1 |
|
fit_phase = kwargs.get("fit_phase", False) |
430
|
1 |
|
scale = kwargs.get("scale", 1e-6) |
431
|
1 |
|
tfmt = kwargs.get("time_format", "jyear") |
432
|
|
|
|
433
|
1 |
|
max_amp = kwargs.pop("max_amp", 1e10 * scale) |
434
|
1 |
|
max_days = kwargs.pop("max_days", 100) |
435
|
|
|
|
436
|
1 |
|
offset_model = [] |
437
|
1 |
|
if constant: |
438
|
1 |
|
offset_model = [("offset", |
439
|
|
|
ConstantModel(value=0., |
440
|
|
|
bounds={"value": [-max_amp, max_amp]}))] |
441
|
|
|
|
442
|
1 |
|
freqs = freqs or [] |
443
|
1 |
|
harmonic_models = [] |
444
|
1 |
|
for freq in freqs: |
445
|
1 |
|
if not fit_phase: |
446
|
1 |
|
harm = HarmonicModelCosineSine(freq=freq, |
447
|
|
|
cos=0, sin=0, |
448
|
|
|
bounds=dict([ |
449
|
|
|
("cos", [-max_amp, max_amp]), |
450
|
|
|
("sin", [-max_amp, max_amp])]) |
451
|
|
|
) |
452
|
|
|
else: |
453
|
|
|
harm = HarmonicModelAmpPhase(freq=freq, |
454
|
|
|
amp=0, phase=0, |
455
|
|
|
bounds=dict([ |
456
|
|
|
("amp", [0, max_amp]), |
457
|
|
|
("phase", [-np.pi, np.pi])]) |
458
|
|
|
) |
459
|
1 |
|
harm.freeze_parameter("freq") |
460
|
1 |
|
harmonic_models.append(("f{0:.0f}".format(freq), harm)) |
461
|
|
|
|
462
|
1 |
|
proxy_config = proxy_config or _default_proxy_config(tfmt=tfmt) |
463
|
1 |
|
proxy_models = [] |
464
|
1 |
|
for pn, conf in proxy_config.items(): |
465
|
1 |
|
if "max_amp" not in conf: |
466
|
|
|
conf.update(dict(max_amp=max_amp)) |
467
|
1 |
|
if "max_days" not in conf: |
468
|
|
|
conf.update(dict(max_days=max_days)) |
469
|
1 |
|
kw = kwargs.copy() # don't mess with the passed arguments |
470
|
1 |
|
kw.update(conf) |
471
|
1 |
|
proxy_models.append( |
472
|
|
|
(pn, setup_proxy_model_with_bounds(**kw)) |
473
|
|
|
) |
474
|
|
|
|
475
|
|
|
return TraceGasModelSet(offset_model + harmonic_models + proxy_models) |
476
|
|
|
|