Passed
Push — master ( c3c30b...6ed46e )
by Stefan
04:58
created

HarmonicModelAmpPhase.compute_gradient()   A

Complexity

Conditions 1

Size

Total Lines 6
Code Lines 6

Duplication

Lines 6
Ratio 100 %

Code Coverage

Tests 1
CRAP Score 1.5786

Importance

Changes 0
Metric Value
cc 1
eloc 6
nop 2
dl 6
loc 6
ccs 1
cts 6
cp 0.1666
crap 1.5786
rs 10
c 0
b 0
f 0
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):
0 ignored issues
show
Duplication introduced by
This code seems to be duplicated in your project.
Loading history...
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):
0 ignored issues
show
Duplication introduced by
This code seems to be duplicated in your project.
Loading history...
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):
0 ignored issues
show
Duplication introduced by
This code seems to be duplicated in your project.
Loading history...
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.times = proxy_times
176 1
		self.dt = 1.
177 1
		self.values = proxy_vals - self.mean
178 1
		self.sza_intp = sza_intp
179 1
		self.fit_phase = fit_phase
180 1
		self.days_per_time_unit = days_per_time_unit
181 1
		self.omega = 2 * np.pi * days_per_time_unit / 365.25
182 1
		self.lifetime_prior = lifetime_prior
183 1
		self.lifetime_metric = lifetime_metric
184 1
		super(ProxyModel, self).__init__(*args, **kwargs)
185
186 1
	def _lt_corr(self, t, tau, tmax=60.):
187
		"""Lifetime corrected values
188
189
		Corrects for a finite lifetime by summing over the last `tmax`
190
		days with an exponential decay given of lifetime(s) `taus`.
191
		"""
192 1
		bs = np.arange(self.dt, tmax + self.dt, self.dt)
193 1
		yp = np.zeros_like(t)
194 1
		tauexp = np.exp(-self.dt / tau)
195 1
		taufac = np.ones_like(tau)
196 1
		for b in bs:
197 1
			taufac *= tauexp
198 1
			yp += np.interp(t - self.lag - b / self.days_per_time_unit,
199
					self.times, self.values, left=0., right=0.) * taufac
200 1
		return yp * self.dt
201
202 1
	def _lt_corr_grad(self, t, tau, tmax=60.):
203
		"""Lifetime corrected gradient
204
205
		Corrects for a finite lifetime by summing over the last `tmax`
206
		days with an exponential decay given of lifetime(s) `taus`.
207
		"""
208 1
		bs = np.arange(self.dt, tmax + self.dt, self.dt)
209 1
		ypg = np.zeros_like(t)
210 1
		tauexp = np.exp(-self.dt / tau)
211 1
		taufac = np.ones_like(tau)
212 1
		for b in bs:
213 1
			taufac *= tauexp
214 1
			ypg += np.interp(t - self.lag - b / self.days_per_time_unit,
215
					self.times, self.values, left=0., right=0.) * taufac * b
216 1
		return ypg * self.dt / tau**2
217
218 1
	def get_value(self, t):
219 1
		t = np.atleast_1d(t)
220 1
		proxy_val = np.interp(t - self.lag,
221
				self.times, self.values, left=0., right=0.)
222 1
		if self.ltscan == 0:
223
			# no lifetime, nothing else to do
224 1
			return self.amp * proxy_val
225
		# annual variation of the proxy lifetime
226 1
		if self.sza_intp is not None:
227
			# using the solar zenith angle
228
			tau_cs = (self.taucos1 * np.cos(np.radians(self.sza_intp(t)))
229
					+ self.tausin1 * np.sin(np.radians(self.sza_intp(t))))
230 1
		elif self.fit_phase:
231
			# using time (cos) and phase (sin)
232
			tau_cs = (self.taucos1 * np.cos(1 * self.omega * t + self.tausin1)
233
					+ self.taucos2 * np.cos(2 * self.omega * t + self.tausin2))
234
		else:
235
			# using time
236 1
			tau_cs = (self.taucos1 * np.cos(1 * self.omega * t)
237
					+ self.tausin1 * np.sin(1 * self.omega * t)
238
					+ self.taucos2 * np.cos(2 * self.omega * t)
239
					+ self.tausin2 * np.sin(2 * self.omega * t))
240 1
		tau_cs[tau_cs < 0] = 0.  # clip to zero
241 1
		tau = self.tau0 + tau_cs
242 1
		if self.ltscan > 0:
243 1
			_ltscn = int(np.floor(self.ltscan))
244
		else:
245
			# infer the scan time from the maximal lifetime
246
			_ltscn = 3 * int(np.ceil(self.tau0 +
247
						np.sqrt(self.taucos1**2 + self.tausin1**2)))
248 1
		if np.all(tau > 0):
249 1
			proxy_val += self._lt_corr(t, tau, tmax=_ltscn)
250 1
		return self.amp * proxy_val
251
252 1
	def compute_gradient(self, t):
253 1
		t = np.atleast_1d(t)
254 1
		proxy_val = np.interp(t - self.lag,
255
				self.times, self.values, left=0., right=0.)
256 1
		proxy_val_grad0 = proxy_val.copy()
257
		# annual variation of the proxy lifetime
258 1
		if self.sza_intp is not None:
259
			# using the solar zenith angle
260
			dtau_cos1 = np.cos(np.radians(self.sza_intp(t)))
261
			dtau_sin1 = np.sin(np.radians(self.sza_intp(t)))
262
			dtau_cos2 = np.zeros_like(t)
263
			dtau_sin2 = np.zeros_like(t)
264
			tau_cs = self.taucos1 * dtau_cos1 + self.tausin1 * dtau_sin1
265 1
		elif self.fit_phase:
266
			# using time (cos) and phase (sin)
267
			dtau_cos1 = np.cos(1 * self.omega * t + self.tausin1)
268
			dtau_sin1 = -self.taucos1 * np.sin(1 * self.omega * t + self.tausin1)
269
			dtau_cos2 = np.cos(2 * self.omega * t + self.tausin2)
270
			dtau_sin2 = -self.taucos2 * np.sin(2 * self.omega * t + self.tausin2)
271
			tau_cs = self.taucos1 * dtau_cos1 + self.taucos2 * dtau_cos2
272
		else:
273
			# using time
274 1
			dtau_cos1 = np.cos(1 * self.omega * t)
275 1
			dtau_sin1 = np.sin(1 * self.omega * t)
276 1
			dtau_cos2 = np.cos(2 * self.omega * t)
277 1
			dtau_sin2 = np.sin(2 * self.omega * t)
278 1
			tau_cs = (self.taucos1 * dtau_cos1 + self.tausin1 * dtau_sin1 +
279
					self.taucos2 * dtau_cos2 + self.tausin2 * dtau_sin2)
280 1
		tau_cs[tau_cs < 0] = 0.  # clip to zero
281 1
		tau = self.tau0 + tau_cs
282 1
		if self.ltscan > 0:
283 1
			_ltscn = int(np.floor(self.ltscan))
284
		else:
285
			# infer the scan time from the maximal lifetime
286 1
			_ltscn = 3 * int(np.ceil(self.tau0 +
287
						np.sqrt(self.taucos1**2 + self.tausin1**2)))
288 1
		if np.all(tau > 0):
289 1
			proxy_val += self._lt_corr(t, tau, tmax=_ltscn)
290 1
			proxy_val_grad0 += self._lt_corr_grad(t, tau, tmax=_ltscn)
291 1
		return np.array([proxy_val,
292
				# set the gradient wrt lag to zero for now
293
				np.zeros_like(t),
294
				self.amp * proxy_val_grad0,
295
				self.amp * proxy_val_grad0 * dtau_cos1,
296
				self.amp * proxy_val_grad0 * dtau_sin1,
297
				self.amp * proxy_val_grad0 * dtau_cos2,
298
				self.amp * proxy_val_grad0 * dtau_sin2,
299
				# set the gradient wrt lifetime scan to zero for now
300
				np.zeros_like(t)])
301
302 1
	def _log_prior_normal(self):
303
		l_prior = super(ProxyModel, self).log_prior()
304
		if not np.isfinite(l_prior):
305
			return -np.inf
306
		for n, p in self.get_parameter_dict().items():
307
			if n.startswith("tau"):
308
				# Gaussian prior for the lifetimes
309
				l_prior -= 0.5 * (p / self.lifetime_metric)**2
310
		return l_prior
311
312 1
	def _log_prior_exp(self):
313 1
		l_prior = super(ProxyModel, self).log_prior()
314 1
		if not np.isfinite(l_prior):
315 1
			return -np.inf
316 1
		for n, p in self.get_parameter_dict().items():
317 1
			if n.startswith("tau"):
318
				# exponential prior for the lifetimes
319 1
				l_prior -= np.abs(p / self.lifetime_metric)
320 1
		return l_prior
321
322 1
	def log_prior(self):
323 1
		_priors = {"exp": self._log_prior_exp,
324
				"normal": self._log_prior_normal}
325 1
		if self.lifetime_prior is None or self.lifetime_prior == "flat":
326 1
			return super(ProxyModel, self).log_prior()
327 1
		return _priors[self.lifetime_prior]()
328
329
330 1 View Code Duplication
class TraceGasModelSet(ModelSet):
0 ignored issues
show
Duplication introduced by
This code seems to be duplicated in your project.
Loading history...
331
	"""Combined model class for trace gases (and probably other data)
332
333
	Inherited from :class:`celerite.ModelSet`, provides `get_value()`
334
	and `compute_gradient()` methods.
335
	"""
336 1
	def get_value(self, t):
337 1
		v = np.zeros_like(t)
338 1
		for m in self.models.values():
339 1
			v += m.get_value(t)
340 1
		return v
341
342 1
	def compute_gradient(self, t):
343 1
		grad = []
344 1
		for m in self.models.values():
345 1
			grad.extend(list(m.compute_gradient(t)))
346 1
		return np.array(grad)
347
348
349 1 View Code Duplication
def setup_proxy_model_with_bounds(times, values,
0 ignored issues
show
Duplication introduced by
This code seems to be duplicated in your project.
Loading history...
350
		max_amp=1e10, max_days=100,
351
		**kwargs):
352
	# extract setup from `kwargs`
353 1
	center = kwargs.get("center", False)
354 1
	fit_phase = kwargs.get("fit_phase", False)
355 1
	lag = kwargs.get("lag", 0.)
356 1
	lt_metric = kwargs.get("lifetime_metric", 1)
357 1
	lt_prior = kwargs.get("lifetime_prior", "exp")
358 1
	lt_scan = kwargs.get("lifetime_scan", 60)
359 1
	positive = kwargs.get("positive", False)
360 1
	sza_intp = kwargs.get("sza_intp", None)
361 1
	time_format = kwargs.get("time_format", "jyear")
362
363 1
	return ProxyModel(times, values,
364
			center=center,
365
			sza_intp=sza_intp,
366
			fit_phase=fit_phase,
367
			lifetime_prior=lt_prior,
368
			lifetime_metric=lt_metric,
369
			days_per_time_unit=1 if time_format.endswith("d") else 365.25,
370
			amp=0.,
371
			lag=lag,
372
			tau0=0,
373
			taucos1=0, tausin1=0,
374
			taucos2=0, tausin2=0,
375
			ltscan=lt_scan,
376
			bounds=dict([
377
				("amp", [0, max_amp] if positive else [-max_amp, max_amp]),
378
				("lag", [0, max_days]),
379
				("tau0", [0, max_days]),
380
				("taucos1", [0, max_days] if fit_phase else [-max_days, max_days]),
381
				("tausin1", [-np.pi, np.pi] if fit_phase else [-max_days, max_days]),
382
				# semi-annual cycles for the life time
383
				("taucos2", [0, max_days] if fit_phase else [-max_days, max_days]),
384
				("tausin2", [-np.pi, np.pi] if fit_phase else [-max_days, max_days]),
385
				("ltscan", [0, 200])])
386
			)
387
388
389 1 View Code Duplication
def _default_proxy_config(tfmt="jyear"):
0 ignored issues
show
Duplication introduced by
This code seems to be duplicated in your project.
Loading history...
390
	from .load_data import load_dailymeanLya, load_dailymeanAE
391
	proxy_config = {}
392
	# Lyman-alpha
393
	plyat, plyadf = load_dailymeanLya(tfmt=tfmt)
394
	proxy_config.update({"Lya": {
395
		"times": plyat,
396
		"values": plyadf["Lya"],
397
		"center": False,
398
		"positive": False,
399
		"lifetime_scan": 0,
400
		}}
401
	)
402
	# AE index
403
	paet, paedf = load_dailymeanAE(name="GM", tfmt=tfmt)
404
	proxy_config.update({"GM": {
405
		"times": paet,
406
		"values": paedf["GM"],
407
		"center": False,
408
		"positive": True,
409
		"lifetime_scan": 60,
410
		}}
411
	)
412
	return proxy_config
413
414
415 1 View Code Duplication
def trace_gas_model(constant=True, freqs=None, proxy_config=None, **kwargs):
0 ignored issues
show
Duplication introduced by
This code seems to be duplicated in your project.
Loading history...
416
	"""Trace gas model setup
417
418
	Sets up the trace gas model for easy access. All parameters are optional,
419
	defaults to an offset, no harmonics, proxies are uncentered and unscaled
420
	Lyman-alpha and AE. AE with only positive amplitude and a seasonally
421
	varying lifetime.
422
423
	Parameters
424
	----------
425
	constant : bool, optional
426
		Whether or not to include a constant (offset) term, default is True.
427
	freqs : list, optional
428
		Frequencies of the harmonic terms in 1 / a^-1 (inverse years).
429
	proxy_config : dict, optional
430
		Proxy configuration if different from the standard setup.
431
	**kwargs : optional
432
		Additional keyword arguments, all of them are also passed on to
433
		the proxy setup. For now, supported are the following which are
434
		also passed along to the proxy setup with
435
		`setup_proxy_model_with_bounds()`:
436
437
		* fit_phase : bool
438
			fit amplitude and phase instead of sine and cosine
439
		* scale : float
440
			the factor by which the data is scaled, used to constrain
441
			the maximum and minimum amplitudes to be fitted.
442
		* tfmt : string
443
			The `astropy.time.Time` format string to setup the time axis.
444
		* max_amp : float
445
			Maximum magnitude of the coefficients, used to constrain the
446
			parameter search.
447
		* max_days : float
448
			Maximum magnitude of the lifetimes, used to constrain the
449
			parameter search.
450
451
	Returns
452
	-------
453
	model : :class:`TraceGasModelSet` (extended :class:`celerite.ModelSet`)
454
	"""
455 1
	fit_phase = kwargs.get("fit_phase", False)
456 1
	scale = kwargs.get("scale", 1e-6)
457 1
	tfmt = kwargs.get("time_format", "jyear")
458
459 1
	max_amp = kwargs.pop("max_amp", 1e10 * scale)
460 1
	max_days = kwargs.pop("max_days", 100)
461
462 1
	offset_model = []
463 1
	if constant:
464 1
		offset_model = [("offset",
465
				ConstantModel(value=0.,
466
						bounds={"value": [-max_amp, max_amp]}))]
467
468 1
	freqs = freqs or []
469 1
	harmonic_models = []
470 1
	for freq in freqs:
471 1
		if not fit_phase:
472 1
			harm = HarmonicModelCosineSine(freq=freq,
473
					cos=0, sin=0,
474
					bounds=dict([
475
						("cos", [-max_amp, max_amp]),
476
						("sin", [-max_amp, max_amp])])
477
			)
478
		else:
479
			harm = HarmonicModelAmpPhase(freq=freq,
480
					amp=0, phase=0,
481
					bounds=dict([
482
						("amp", [0, max_amp]),
483
						("phase", [-np.pi, np.pi])])
484
			)
485 1
		harm.freeze_parameter("freq")
486 1
		harmonic_models.append(("f{0:.0f}".format(freq), harm))
487
488 1
	proxy_config = proxy_config or _default_proxy_config(tfmt=tfmt)
489 1
	proxy_models = []
490 1
	for pn, conf in proxy_config.items():
491 1
		if "max_amp" not in conf:
492
			conf.update(dict(max_amp=max_amp))
493 1
		if "max_days" not in conf:
494
			conf.update(dict(max_days=max_days))
495 1
		kw = kwargs.copy()  # don't mess with the passed arguments
496 1
		kw.update(conf)
497 1
		proxy_models.append(
498
			(pn, setup_proxy_model_with_bounds(**kw))
499
		)
500
501
	return TraceGasModelSet(offset_model + harmonic_models + proxy_models)
502