PhaseCurveTemplateTemporalModel

class gammapy.modeling.models.PhaseCurveTemplateTemporalModel(table, time_0, phase_0, f0, f1=0, f2=0)[source]

Bases: gammapy.modeling.models.TemporalModel

Temporal phase curve model.

Phase for a given time is computed as:

\[\phi(t) = \phi_0 + f_0(t-t_0) + (1/2)f_1(t-t_0)^2 + (1/6)f_2(t-t_0)^3\]

Strictly periodic sources such as gamma-ray binaries have f1=0 and f2=0. Sources like some pulsars where the period spins up or down have f1!=0 and / or f2 !=0. For a binary, f0 should be calculated as 1/T, where T is the period of the binary in unit of seconds.

The “phase curve”, i.e. multiplicative flux factor for a given phase is given by a Table of nodes (phase, norm), using linear interpolation and circular behaviour, where norm(phase=0) == norm(phase=1).

Parameters:
table : Table

A table of ‘PHASE’ vs ‘NORM’ should be given

time_0 : float

The MJD value where phase is considered as 0.

phase_0 : float

Phase at the reference MJD

f0, f1, f2 : float

Derivatives of the function phi with time of order 1, 2, 3 in units of s^-1, s^-2 & s^-3, respectively.

Examples

Create an example phase curve object:

from astropy.table import Table
from gammapy.utils.scripts import make_path
from gammapy.modeling.models import PhaseCurveTemplateTemporalModel
filename = make_path('$GAMMAPY_DATA/tests/phasecurve_LSI_DC.fits')
table = Table.read(str(filename))
phase_curve = PhaseCurveTemplateTemporalModel(table, time_0=43366.275, phase_0=0.0, f0=4.367575e-7, f1=0.0, f2=0.0)

Use it to compute a phase and evaluate the phase curve model for a given time:

>>> phase_curve.phase(time=46300.0)
0.7066006737999402
>>> phase_curve.evaluate_norm_at_time(46300)
0.49059393580053845

Attributes Summary

f0
f1
f2
parameters Parameters (Parameters)
phase_0
table
tag
time_0

Methods Summary

copy(self) A deep copy.
create(tag, \*args, \*\*kwargs) Create a model instance.
evaluate_norm_at_phase(self, phase)
evaluate_norm_at_time(self, time) Evaluate for a given time.
from_dict(data)
phase(self, time) Evaluate phase for a given time.
sample_time(self, n_events, t_min, t_max[, …]) Sample arrival times of events.
to_dict(self)

Attributes Documentation

f0
f1
f2
parameters

Parameters (Parameters)

phase_0
table
tag = 'PhaseCurveTemplateTemporalModel'
time_0

Methods Documentation

copy(self)

A deep copy.

static create(tag, *args, **kwargs)

Create a model instance.

Examples

>>> from gammapy.modeling import Model
>>> spectral_model = Model.create("PowerLaw2SpectralModel", amplitude="1e-10 cm-2 s-1", index=3)
>>> type(spectral_model)
gammapy.modeling.models.spectral.PowerLaw2SpectralModel
evaluate_norm_at_phase(self, phase)[source]
evaluate_norm_at_time(self, time)[source]

Evaluate for a given time.

Parameters:
time : array_like

Time since the reference time.

Returns:
norm : array_like
classmethod from_dict(data)
phase(self, time)[source]

Evaluate phase for a given time.

Parameters:
time : array_like
Returns:
phase : array_like
sample_time(self, n_events, t_min, t_max, t_delta='1 s', random_state=0)[source]

Sample arrival times of events.

Parameters:
n_events : int

Number of events to sample.

t_min : Time

Start time of the sampling.

t_max : Time

Stop time of the sampling.

t_delta : Quantity

Time step used for sampling of the temporal model.

random_state : {int, ‘random-seed’, ‘global-rng’, RandomState}

Defines random number generator initialisation. Passed to get_random_state.

Returns:
time : Quantity

Array with times of the sampled events.

to_dict(self)