-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathtemporal_test.py
449 lines (351 loc) · 13.4 KB
/
temporal_test.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
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
# Licensed under a 3-clause BSD style license - see LICENSE.rst
"""Time-dependent models."""
import numpy as np
import scipy.interpolate
from astropy import units as u
from astropy.table import Table
from astropy.time import Time
from astropy.utils import lazyproperty
from gammapy.modeling import Parameter
from gammapy.utils.random import InverseCDFSampler, get_random_state
from gammapy.utils.scripts import make_path
from gammapy.utils.time import time_ref_from_dict
from .core import Model
# TODO: make this a small ABC to define a uniform interface.
class TemporalModel(Model):
"""Temporal model base class.
evaluates on astropy.time.Time objects"""
def __call__(self, time):
"""Evaluate model
Parameters
----------
time : `~astropy.time.Time`
Time object
"""
kwargs = {par.name: par.quantity for par in self.parameters}
time = u.Quantity(time.mjd, "day")
return self.evaluate(time, **kwargs)
@staticmethod
def time_sum(t_min, t_max):
"""
Total time between t_min and t_max
Parameters
----------
t_min, t_max: `~astropy.time.Time`
Lower and upper bound of integration range
Returns
-------
time_sum : `~astropy.time.TimeDelta`
Summed time in the intervals.
"""
return np.sum(t_max - t_min)
class ConstantTemporalModel(TemporalModel):
"""Constant temporal model."""
tag = "ConstantTemporalModel"
@staticmethod
def evaluate(time):
"""Evaluate at given times."""
return np.ones(time.shape)
def integral(self, t_min, t_max):
"""Evaluate the integrated flux within the given time intervals
Parameters
----------
t_min: `~astropy.time.Time`
Start times of observation
t_max: `~astropy.time.Time`
Stop times of observation
Returns
-------
norm : `~astropy.units.Quantity`
Integrated flux norm on the given time intervals
"""
return (t_max - t_min) / self.time_sum(t_min, t_max)
@staticmethod
def sample_time(n_events, t_min, t_max, random_state=0):
"""Sample arrival times of events.
Parameters
----------
n_events : int
Number of events to sample.
t_min : `~astropy.time.Time`
Start time of the sampling.
t_max : `~astropy.time.Time`
Stop time of the sampling.
random_state : {int, 'random-seed', 'global-rng', `~numpy.random.RandomState`}
Defines random number generator initialisation.
Passed to `~gammapy.utils.random.get_random_state`.
Returns
-------
time : `~astropy.units.Quantity`
Array with times of the sampled events.
"""
random_state = get_random_state(random_state)
t_min = Time(t_min)
t_max = Time(t_max)
t_stop = (t_max - t_min).sec
time_delta = random_state.uniform(high=t_stop, size=n_events) * u.s
return t_min + time_delta
class ExpDecayTemporalModel(TemporalModel):
"""Temporal model with an exponential decay.
..math::
F(t) = exp(t - t_ref)/t0
Parameters
----------
t0 : `~astropy.units.Quantity`
Decay time scale
t_ref: `~astropy.units.Quantity`
The reference time in mjd
"""
tag = "ExponentialDecayTemporalModel"
t0 = Parameter("t0", "1 d", frozen=False)
_t_ref_default = Time("2000-01-01")
t_ref = Parameter("t_ref", _t_ref_default.mjd, unit="day", frozen=True)
@staticmethod
def evaluate(time, t0, t_ref):
"""Evaluate at given times"""
return np.exp(-(time - t_ref) / t0)
def integral(self, t_min, t_max):
"""Evaluate the integrated flux within the given time intervals
Parameters
----------
t_min: `~astropy.time.Time`
Start times of observation
t_max: `~astropy.time.Time`
Stop times of observation
Returns
-------
norm : float
Integrated flux norm on the given time intervals
"""
pars = self.parameters
t0 = pars["t0"].quantity
t_ref = Time(pars["t_ref"].quantity, format="mjd")
value = self.evaluate(t_max, t0, t_ref) - self.evaluate(t_min, t0, t_ref)
return -t0 * value / self.time_sum(t_min, t_max)
class GaussianTemporalModel(TemporalModel):
"""A Gaussian temporal profile
Parameters
----------
t_ref: `~astropy.units.Quantity`
The reference time in mjd at the peak.
sigma : `~astropy.units.Quantity`
Width of the gaussian profile.
"""
tag = "GaussianTemporalModel"
_t_ref_default = Time("2000-01-01")
t_ref = Parameter("t_ref", _t_ref_default.mjd, unit="day", frozen=False)
sigma = Parameter("sigma", "1 d", frozen=False)
@staticmethod
def evaluate(time, t_ref, sigma):
return np.exp(-((time - t_ref) ** 2) / (2 * sigma ** 2))
def integral(self, t_min, t_max, **kwargs):
"""Evaluate the integrated flux within the given time intervals
Parameters
----------
t_min: `~astropy.time.Time`
Start times of observation
t_max: `~astropy.time.Time`
Stop times of observation
Returns
-------
norm : float
Integrated flux norm on the given time intervals
"""
pars = self.parameters
sigma = pars["sigma"].quantity
t_ref = Time(pars["t_ref"].quantity, format="mjd")
norm = np.sqrt(np.pi / 2) * sigma
u_min = (t_min - t_ref) / (np.sqrt(2) * sigma)
u_max = (t_max - t_ref) / (np.sqrt(2) * sigma)
integral = norm * (scipy.special.erf(u_max) - scipy.special.erf(u_min))
return integral / self.time_sum(t_min, t_max)
class PowerlawDecayTemporalModel(TemporalModel):
"""Temporal model with an exponential decay.
..math::
F(t) = exp(t - t_ref)/t0
Parameters
----------
t0 : `~astropy.units.Quantity`
Decay time scale
t_ref: `~astropy.units.Quantity`
The reference time in mjd
"""
tag = "PowerlawDecayTemporalModel"
t0 = Parameter("t0", "1 d", frozen=False)
_t_ref_default = Time("2000-01-01")
t_ref = Parameter("t_ref", _t_ref_default.mjd, unit="day", frozen=True)
@staticmethod
def evaluate(time, t0, t_ref):
"""Evaluate at given times"""
return np.exp(-(time - t_ref) / t0)
def integral(self, t_min, t_max):
"""Evaluate the integrated flux within the given time intervals
Parameters
----------
t_min: `~astropy.time.Time`
Start times of observation
t_max: `~astropy.time.Time`
Stop times of observation
Returns
-------
norm : float
Integrated flux norm on the given time intervals
"""
pars = self.parameters
t0 = pars["t0"].quantity
t_ref = Time(pars["t_ref"].quantity, format="mjd")
value = self.evaluate(t_max, t0, t_ref) - self.evaluate(t_min, t0, t_ref)
return -t0 * value / self.time_sum(t_min, t_max)
class LightCurveTemplateTemporalModel(TemporalModel):
"""Temporal light curve model.
The lightcurve is given as a table with columns ``time`` and ``norm``.
The ``norm`` is supposed to be a unit-less multiplicative factor in the model,
to be multiplied with a spectral model.
The model does linear interpolation for times between the given ``(time, norm)`` values.
The implementation currently uses `scipy.interpolate.InterpolatedUnivariateSpline`,
using degree ``k=1`` to get linear interpolation.
This class also contains an ``integral`` method, making the computation of
mean fluxes for a given time interval a one-liner.
Parameters
----------
table : `~astropy.table.Table`
A table with 'TIME' vs 'NORM'
Examples
--------
Read an example light curve object:
>>> from gammapy.modeling.models import LightCurveTemplateTemporalModel
>>> path = '$GAMMAPY_DATA/tests/models/light_curve/lightcrv_PKSB1222+216.fits'
>>> light_curve = LightCurveTemplateTemporalModel.read(path)
Show basic information about the lightcurve:
>>> print(light_curve)
LightCurve model summary:
Start time: 59000.5 MJD
End time: 61862.5 MJD
Norm min: 0.01551196351647377
Norm max: 1.0
Compute ``norm`` at a given time:
>>> light_curve.evaluate(46300)
0.49059393580053845
Compute mean ``norm`` in a given time interval:
>>> light_curve.mean_norm_in_time_interval(46300, 46301)
"""
tag = "LightCurveTemplateTemporalModel"
def __init__(self, table, filename=None):
self.table = table
if filename is not None:
filename = str(make_path(filename))
self.filename = filename
super().__init__()
def __str__(self):
norm = self.table["NORM"]
return (
f"{self.__class__.__name__} model summary:\n"
f"Start time: {self._time[0].mjd} MJD\n"
f"End time: {self._time[-1].mjd} MJD\n"
f"Norm min: {norm.min()}\n"
f"Norm max: {norm.max()}\n"
)
@classmethod
def read(cls, path):
"""Read lightcurve model table from FITS file.
TODO: This doesn't read the XML part of the model yet.
"""
filename = str(make_path(path))
return cls(Table.read(filename), filename=filename)
def write(self, path=None, overwrite=False):
if path is None:
path = self.filename
if path is None:
raise ValueError(f"filename is required for {self.tag}")
else:
self.filename = str(make_path(path))
self.table.write(self.filename, overwrite=overwrite)
@lazyproperty
def _interpolator(self, ext=0):
x = self._time.value
y = self.table["NORM"].data
return scipy.interpolate.InterpolatedUnivariateSpline(x, y, k=1, ext=ext)
@lazyproperty
def _time_ref(self):
return time_ref_from_dict(self.table.meta)
@lazyproperty
def _time(self):
return self._time_ref + self.table["TIME"].data * getattr(
u, self.table.meta["TIMEUNIT"]
)
def evaluate(self, time, ext=0):
"""Evaluate for a given time.
Parameters
----------
time : array_like
Time since the ``reference`` time.
ext : int or str, optional, default: 0
Parameter passed to ~scipy.interpolate.InterpolatedUnivariateSpline
Controls the extrapolation mode for GTIs outside the range
0 or "extrapolate", return the extrapolated value.
1 or "zeros", return 0
2 or "raise", raise a ValueError
3 or "const", return the boundary value.
Returns
-------
norm : array_like
Norm at the given times.
"""
return self._interpolator(time, ext=ext)
def integral(self, t_min, t_max):
"""Evaluate the integrated flux within the given time intervals
Parameters
----------
t_min: `~astropy.time.Time`
Start times of observation
t_max: `~astropy.time.Time`
Stop times of observation
Returns
-------
norm: The model integrated flux
"""
n1 = self._interpolator.antiderivative()(t_max.mjd)
n2 = self._interpolator.antiderivative()(t_min.mjd)
return u.Quantity(n1 - n2, "day") / self.time_sum(t_min, t_max)
def sample_time(self, n_events, t_min, t_max, t_delta="1 s", random_state=0):
"""Sample arrival times of events.
Parameters
----------
n_events : int
Number of events to sample.
t_min : `~astropy.time.Time`
Start time of the sampling.
t_max : `~astropy.time.Time`
Stop time of the sampling.
t_delta : `~astropy.units.Quantity`
Time step used for sampling of the temporal model.
random_state : {int, 'random-seed', 'global-rng', `~numpy.random.RandomState`}
Defines random number generator initialisation.
Passed to `~gammapy.utils.random.get_random_state`.
Returns
-------
time : `~astropy.units.Quantity`
Array with times of the sampled events.
"""
time_unit = getattr(u, self.table.meta["TIMEUNIT"])
t_min = Time(t_min)
t_max = Time(t_max)
t_delta = u.Quantity(t_delta)
random_state = get_random_state(random_state)
ontime = u.Quantity((t_max - t_min).sec, "s")
t_stop = ontime.to_value(time_unit)
# TODO: the separate time unit handling is unfortunate, but the quantity support for np.arange and np.interp
# is still incomplete, refactor once we change to recent numpy and astropy versions
t_step = t_delta.to_value(time_unit)
t = np.arange(0, t_stop, t_step)
pdf = self.evaluate(t)
sampler = InverseCDFSampler(pdf=pdf, random_state=random_state)
time_pix = sampler.sample(n_events)[0]
time = np.interp(time_pix, np.arange(len(t)), t) * time_unit
return t_min + time
@classmethod
def from_dict(cls, data):
return cls.read(data["filename"])
def to_dict(self, overwrite=False):
"""Create dict for YAML serilisation"""
return {"type": self.tag, "filename": self.filename}