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", "CeleriteModelSet"] |
28
|
|
|
|
29
|
1 |
View Code Duplication |
class HarmonicModelCosineSine(Model): |
|
|
|
|
30
|
|
|
"""Model for harmonic terms |
31
|
|
|
|
32
|
|
|
Models harmonic terms using a cosine and sine part. |
33
|
|
|
The total amplitude and phase can be inferred from that. |
34
|
|
|
|
35
|
|
|
Parameters |
36
|
|
|
---------- |
37
|
|
|
freq : float |
38
|
|
|
The frequency in years^-1 |
39
|
|
|
cos : float |
40
|
|
|
The amplitude of the cosine part |
41
|
|
|
sin : float |
42
|
|
|
The amplitude of the sine part |
43
|
|
|
""" |
44
|
1 |
|
parameter_names = ("freq", "cos", "sin") |
45
|
|
|
|
46
|
1 |
|
def get_value(self, t): |
47
|
|
|
return (self.cos * np.cos(self.freq * 2 * np.pi * t) + |
48
|
|
|
self.sin * np.sin(self.freq * 2 * np.pi * t)) |
49
|
|
|
|
50
|
1 |
|
def get_amplitude(self): |
51
|
|
|
return np.sqrt(self.cos**2 + self.sin**2) |
52
|
|
|
|
53
|
1 |
|
def get_phase(self): |
54
|
|
|
return np.arctan2(self.sin, self.cos) |
55
|
|
|
|
56
|
1 |
|
def compute_gradient(self, t): |
57
|
|
|
dcos = np.cos(self.freq * 2 * np.pi * t) |
58
|
|
|
dsin = np.sin(self.freq * 2 * np.pi * t) |
59
|
|
|
df = 2 * np.pi * t * (self.sin * dcos - self.cos * dsin) |
60
|
|
|
return np.array([df, dcos, dsin]) |
61
|
|
|
|
62
|
|
|
|
63
|
1 |
View Code Duplication |
class HarmonicModelAmpPhase(Model): |
|
|
|
|
64
|
|
|
"""Model for harmonic terms |
65
|
|
|
|
66
|
|
|
Models harmonic terms using a cosine and sine part. |
67
|
|
|
The total amplitude and phase can be inferred from that. |
68
|
|
|
|
69
|
|
|
Parameters |
70
|
|
|
---------- |
71
|
|
|
freq : float |
72
|
|
|
The frequency in years^-1 |
73
|
|
|
amp : float |
74
|
|
|
The amplitude of the harmonic term |
75
|
|
|
phase : float |
76
|
|
|
The phase of the harmonic part |
77
|
|
|
""" |
78
|
1 |
|
parameter_names = ("freq", "amp", "phase") |
79
|
|
|
|
80
|
1 |
|
def get_value(self, t): |
81
|
|
|
return self.amp * np.cos(self.freq * 2 * np.pi * t + self.phase) |
82
|
|
|
|
83
|
1 |
|
def get_amplitude(self): |
84
|
|
|
return self.amp |
85
|
|
|
|
86
|
1 |
|
def get_phase(self): |
87
|
|
|
return self.phase |
88
|
|
|
|
89
|
1 |
|
def compute_gradient(self, t): |
90
|
|
|
damp = np.cos(self.freq * 2 * np.pi * t + self.phase) |
91
|
|
|
dphi = -self.amp * np.sin(self.freq * 2 * np.pi * t + self.phase) |
92
|
|
|
df = 2 * np.pi * t * dphi |
93
|
|
|
return np.array([df, damp, dphi]) |
94
|
|
|
|
95
|
|
|
|
96
|
1 |
View Code Duplication |
class ProxyModel(Model): |
|
|
|
|
97
|
|
|
"""Model for proxy terms |
98
|
|
|
|
99
|
|
|
Models proxy terms with a finite and (semi-)annually varying life time. |
100
|
|
|
|
101
|
|
|
Parameters |
102
|
|
|
---------- |
103
|
|
|
proxy_times : (N,) array_like |
104
|
|
|
The data times of the proxy values |
105
|
|
|
proxy_vals : (N,) array_like |
106
|
|
|
The proxy values at `proxy_times` |
107
|
|
|
amp : float |
108
|
|
|
The amplitude of the proxy term |
109
|
|
|
lag : float |
110
|
|
|
The lag of the proxy value in years. |
111
|
|
|
tau0 : float |
112
|
|
|
The base life time of the proxy |
113
|
|
|
taucos1 : float |
114
|
|
|
The amplitude of the cosine part of the annual life time variation. |
115
|
|
|
tausin1 : float |
116
|
|
|
The amplitude of the sine part of the annual life time variation. |
117
|
|
|
taucos2 : float |
118
|
|
|
The amplitude of the cosine part of the semi-annual life time variation. |
119
|
|
|
tausin2 : float |
120
|
|
|
The amplitude of the sine part of the semi-annual life time variation. |
121
|
|
|
ltscan : float |
122
|
|
|
The number of days to sum the previous proxy values. If it is |
123
|
|
|
negative, the value will be set to three times the maximal lifetime. |
124
|
|
|
No lifetime adjustemets are calculated when set to zero. |
125
|
|
|
center : bool, optional |
126
|
|
|
Centers the proxy values by subtracting the overall mean. The mean is |
127
|
|
|
calculated from the whole `proxy_vals` array and is stored in the |
128
|
|
|
`mean` attribute. |
129
|
|
|
Default: False |
130
|
|
|
sza_intp : scipy.interpolate.interp1d() instance, optional |
131
|
|
|
When not `None`, cos(sza) and sin(sza) are used instead |
132
|
|
|
of the time to model the annual variation of the lifetime. |
133
|
|
|
Semi-annual variations are not used in that case. |
134
|
|
|
Default: None |
135
|
|
|
fit_phase : bool, optional |
136
|
|
|
Fit the phase shift directly instead of using sine and cosine |
137
|
|
|
terms for the (semi-)annual lifetime variations. If True, the fitted |
138
|
|
|
cosine parameter is the amplitude and the sine parameter the phase. |
139
|
|
|
Default: False (= fit sine and cosine terms) |
140
|
|
|
lifetime_prior : str, optional |
141
|
|
|
The prior probability density for each coefficient of the lifetime. |
142
|
|
|
Possible types are "flat" or `None` for a flat prior, "exp" for an |
143
|
|
|
exponential density ~ :math:`\\text{exp}(-|\\tau| / \\text{metric})`, |
144
|
|
|
and "normal" for a normal distribution |
145
|
|
|
~ :math:`\\text{exp}(-\\tau^2 / (2 * \\text{metric}^2))`. |
146
|
|
|
Default: None (= flat prior). |
147
|
|
|
lifetime_metric : float, optional |
148
|
|
|
The metric (scale) of the lifetime priors in days, see `prior`. |
149
|
|
|
Default 1. |
150
|
|
|
days_per_time_unit : float, optional |
151
|
|
|
The number of days per time unit, used to normalize the lifetime |
152
|
|
|
units. Use 365.25 if the times are in fractional years, or 1 if |
153
|
|
|
they are in days. |
154
|
|
|
Default: 365.25 |
155
|
|
|
""" |
156
|
1 |
|
parameter_names = ("amp", "lag", "tau0", |
157
|
|
|
"taucos1", "tausin1", "taucos2", "tausin2", |
158
|
|
|
"ltscan") |
159
|
|
|
|
160
|
1 |
|
def __init__(self, proxy_times, proxy_vals, |
161
|
|
|
center=False, |
162
|
|
|
sza_intp=None, fit_phase=False, |
163
|
|
|
lifetime_prior=None, lifetime_metric=1., |
164
|
|
|
days_per_time_unit=365.25, |
165
|
|
|
*args, **kwargs): |
166
|
|
|
self.mean = 0. |
167
|
|
|
if center: |
168
|
|
|
self.mean = np.nanmean(proxy_vals) |
169
|
|
|
self.intp = interp1d(proxy_times, proxy_vals - self.mean, |
170
|
|
|
bounds_error=False) |
171
|
|
|
self.sza_intp = sza_intp |
172
|
|
|
self.fit_phase = fit_phase |
173
|
|
|
self.days_per_time_unit = days_per_time_unit |
174
|
|
|
self.omega = 2 * np.pi * days_per_time_unit / 365.25 |
175
|
|
|
self.lifetime_prior = lifetime_prior |
176
|
|
|
self.lifetime_metric = lifetime_metric |
177
|
|
|
super(ProxyModel, self).__init__(*args, **kwargs) |
178
|
|
|
|
179
|
1 |
|
def get_value(self, t): |
180
|
|
|
proxy_val = self.intp(t - self.lag) |
181
|
|
|
if self.ltscan == 0: |
182
|
|
|
# no lifetime, nothing else to do |
183
|
|
|
return self.amp * proxy_val |
184
|
|
|
# annual variation of the proxy lifetime |
185
|
|
|
if self.sza_intp is not None: |
186
|
|
|
# using the solar zenith angle |
187
|
|
|
tau_cs = (self.taucos1 * np.cos(np.radians(self.sza_intp(t))) |
188
|
|
|
+ self.tausin1 * np.sin(np.radians(self.sza_intp(t)))) |
189
|
|
|
elif self.fit_phase: |
190
|
|
|
# using time (cos) and phase (sin) |
191
|
|
|
tau_cs = (self.taucos1 * np.cos(1 * self.omega * t + self.tausin1) |
192
|
|
|
+ self.taucos2 * np.cos(2 * self.omega * t + self.tausin2)) |
193
|
|
|
else: |
194
|
|
|
# using time |
195
|
|
|
tau_cs = (self.taucos1 * np.cos(1 * self.omega * t) |
196
|
|
|
+ self.tausin1 * np.sin(1 * self.omega * t) |
197
|
|
|
+ self.taucos2 * np.cos(2 * self.omega * t) |
198
|
|
|
+ self.tausin2 * np.sin(2 * self.omega * t)) |
199
|
|
|
tau_cs[tau_cs < 0] = 0. # clip to zero |
200
|
|
|
tau = self.tau0 + tau_cs |
201
|
|
|
if self.ltscan > 0: |
202
|
|
|
_ltscn = int(np.floor(self.ltscan)) |
203
|
|
|
else: |
204
|
|
|
# infer the scan time from the maximal lifetime |
205
|
|
|
_ltscn = 3 * int(np.ceil(self.tau0 + |
206
|
|
|
np.sqrt(self.taucos1**2 + self.tausin1**2))) |
207
|
|
|
if np.all(tau > 0): |
208
|
|
|
bs = np.arange(1, _ltscn + 1, 1.)[None, :] |
209
|
|
|
taufacs = np.exp(-bs / tau[:, None]) |
210
|
|
|
proxy_val += np.sum( |
211
|
|
|
self.intp(t[:, None] - self.lag - |
212
|
|
|
bs / self.days_per_time_unit) * taufacs, |
213
|
|
|
axis=1) |
214
|
|
|
return self.amp * proxy_val |
215
|
|
|
|
216
|
1 |
|
def compute_gradient(self, t): |
217
|
|
|
proxy_val = self.intp(t - self.lag) |
218
|
|
|
proxy_val_grad0 = self.intp(t - self.lag) |
219
|
|
|
# annual variation of the proxy lifetime |
220
|
|
|
if self.sza_intp is not None: |
221
|
|
|
# using the solar zenith angle |
222
|
|
|
dtau_cos1 = np.cos(np.radians(self.sza_intp(t))) |
223
|
|
|
dtau_sin1 = np.sin(np.radians(self.sza_intp(t))) |
224
|
|
|
dtau_cos2 = np.zeros_like(t) |
225
|
|
|
dtau_sin2 = np.zeros_like(t) |
226
|
|
|
tau_cs = self.taucos1 * dtau_cos1 + self.tausin1 * dtau_sin1 |
227
|
|
|
elif self.fit_phase: |
228
|
|
|
# using time (cos) and phase (sin) |
229
|
|
|
dtau_cos1 = np.cos(1 * self.omega * t + self.tausin1) |
230
|
|
|
dtau_sin1 = -self.taucos1 * np.sin(1 * self.omega * t + self.tausin1) |
231
|
|
|
dtau_cos2 = np.cos(2 * self.omega * t + self.tausin2) |
232
|
|
|
dtau_sin2 = -self.taucos2 * np.sin(2 * self.omega * t + self.tausin2) |
233
|
|
|
tau_cs = self.taucos1 * dtau_cos1 + self.taucos2 * dtau_cos2 |
234
|
|
|
else: |
235
|
|
|
# using time |
236
|
|
|
dtau_cos1 = np.cos(1 * self.omega * t) |
237
|
|
|
dtau_sin1 = np.sin(1 * self.omega * t) |
238
|
|
|
dtau_cos2 = np.cos(2 * self.omega * t) |
239
|
|
|
dtau_sin2 = np.sin(2 * self.omega * t) |
240
|
|
|
tau_cs = (self.taucos1 * dtau_cos1 + self.tausin1 * dtau_sin1 + |
241
|
|
|
self.taucos2 * dtau_cos2 + self.tausin2 * dtau_sin2) |
242
|
|
|
tau_cs[tau_cs < 0] = 0. # clip to zero |
243
|
|
|
tau = self.tau0 + tau_cs |
244
|
|
|
if self.ltscan > 0: |
245
|
|
|
_ltscn = int(np.floor(self.ltscan)) |
246
|
|
|
else: |
247
|
|
|
# infer the scan time from the maximal lifetime |
248
|
|
|
_ltscn = 3 * int(np.ceil(self.tau0 + |
249
|
|
|
np.sqrt(self.taucos1**2 + self.tausin1**2))) |
250
|
|
|
if np.all(tau > 0): |
251
|
|
|
bs = np.arange(1, _ltscn + 1, 1.)[None, :] |
252
|
|
|
taufacs = np.exp(-bs / tau[:, None]) |
253
|
|
|
proxy_ts = self.intp(t[:, None] - self.lag - |
254
|
|
|
bs / self.days_per_time_unit) * taufacs |
255
|
|
|
proxy_val += np.sum(proxy_ts, axis=1) |
256
|
|
|
proxy_val_grad0 += np.sum(proxy_ts * bs / tau[:, None]**2, axis=1) |
257
|
|
|
return np.array([proxy_val, |
258
|
|
|
# set the gradient wrt lag to zero for now |
259
|
|
|
np.zeros_like(t), |
260
|
|
|
self.amp * proxy_val_grad0, |
261
|
|
|
self.amp * proxy_val_grad0 * dtau_cos1, |
262
|
|
|
self.amp * proxy_val_grad0 * dtau_sin1, |
263
|
|
|
self.amp * proxy_val_grad0 * dtau_cos2, |
264
|
|
|
self.amp * proxy_val_grad0 * dtau_sin2, |
265
|
|
|
# set the gradient wrt lifetime scan to zero for now |
266
|
|
|
np.zeros_like(t)]) |
267
|
|
|
|
268
|
1 |
|
def _log_prior_normal(self): |
269
|
|
|
l_prior = super(ProxyModel, self).log_prior() |
270
|
|
|
if not np.isfinite(l_prior): |
271
|
|
|
return -np.inf |
272
|
|
|
for n, p in self.get_parameter_dict().items(): |
273
|
|
|
if n.startswith("tau"): |
274
|
|
|
# Gaussian prior for the lifetimes |
275
|
|
|
l_prior -= 0.5 * (p / self.lifetime_metric)**2 |
276
|
|
|
return l_prior |
277
|
|
|
|
278
|
1 |
|
def _log_prior_exp(self): |
279
|
|
|
l_prior = super(ProxyModel, self).log_prior() |
280
|
|
|
if not np.isfinite(l_prior): |
281
|
|
|
return -np.inf |
282
|
|
|
for n, p in self.get_parameter_dict().items(): |
283
|
|
|
if n.startswith("tau"): |
284
|
|
|
# exponential prior for the lifetimes |
285
|
|
|
l_prior -= np.abs(p / self.lifetime_metric) |
286
|
|
|
return l_prior |
287
|
|
|
|
288
|
1 |
|
def log_prior(self): |
289
|
|
|
_priors = {"exp": self._log_prior_exp, |
290
|
|
|
"normal": self._log_prior_normal} |
291
|
|
|
if self.lifetime_prior is None or self.lifetime_prior == "flat": |
292
|
|
|
return super(ProxyModel, self).log_prior() |
293
|
|
|
return _priors[self.lifetime_prior]() |
294
|
|
|
|
295
|
|
|
|
296
|
1 |
View Code Duplication |
class CeleriteModelSet(ModelSet): |
|
|
|
|
297
|
|
|
|
298
|
1 |
|
def get_value(self, t): |
299
|
|
|
v = np.zeros_like(t) |
300
|
|
|
for m in self.models.values(): |
301
|
|
|
v += m.get_value(t) |
302
|
|
|
return v |
303
|
|
|
|
304
|
1 |
|
def compute_gradient(self, t): |
305
|
|
|
grad = [] |
306
|
|
|
for m in self.models.values(): |
307
|
|
|
grad.extend(list(m.compute_gradient(t))) |
308
|
|
|
return np.array(grad) |
309
|
|
|
|
310
|
|
|
|
311
|
1 |
View Code Duplication |
def _setup_proxy_model_with_bounds(times, values, |
|
|
|
|
312
|
|
|
max_amp=1e10, max_days=100, |
313
|
|
|
**kwargs): |
314
|
|
|
# extract setup from `kwargs` |
315
|
|
|
center = kwargs.get("center", False) |
316
|
|
|
fit_phase = kwargs.get("fit_phase", False) |
317
|
|
|
lag = kwargs.get("lag", 0.) |
318
|
|
|
lt_metric = kwargs.get("lifetime_metric", 1) |
319
|
|
|
lt_prior = kwargs.get("lifetime_prior", "exp") |
320
|
|
|
lt_scan = kwargs.get("lifetime_scan", 60) |
321
|
|
|
positive = kwargs.get("positive", False) |
322
|
|
|
sza_intp = kwargs.get("sza_intp", None) |
323
|
|
|
time_format = kwargs.get("time_format", "jyear") |
324
|
|
|
|
325
|
|
|
return ProxyModel(times, values, |
326
|
|
|
center=center, |
327
|
|
|
sza_intp=sza_intp, |
328
|
|
|
fit_phase=fit_phase, |
329
|
|
|
lifetime_prior=lt_prior, |
330
|
|
|
lifetime_metric=lt_metric, |
331
|
|
|
days_per_time_unit=1 if time_format.endswith("d") else 365.25, |
332
|
|
|
amp=0., |
333
|
|
|
lag=lag, |
334
|
|
|
tau0=0, |
335
|
|
|
taucos1=0, tausin1=0, |
336
|
|
|
taucos2=0, tausin2=0, |
337
|
|
|
ltscan=lt_scan, |
338
|
|
|
bounds=dict([ |
339
|
|
|
("amp", [0, max_amp] if positive else [-max_amp, max_amp]), |
340
|
|
|
("lag", [0, max_days]), |
341
|
|
|
("tau0", [0, max_days]), |
342
|
|
|
("taucos1", [0, max_days] if fit_phase else [-max_days, max_days]), |
343
|
|
|
("tausin1", [-np.pi, np.pi] if fit_phase else [-max_days, max_days]), |
344
|
|
|
# semi-annual cycles for the life time |
345
|
|
|
("taucos2", [0, max_days] if fit_phase else [-max_days, max_days]), |
346
|
|
|
("tausin2", [-np.pi, np.pi] if fit_phase else [-max_days, max_days]), |
347
|
|
|
("ltscan", [0, 200])]) |
348
|
|
|
) |
349
|
|
|
|