Utilities
Several useful utilities used throughout various modules
Function tools
- tatpulsar.utils.functions.cal_2dchisquare(data, f, F1, pepoch, nbins, F2=0, F3=0, F4=0)[source]
Calculate the chisquare distribution for 2-D frequency search on the pepoch time. For example, search in a two-dimensianal parameter space (\(M \times N\), as \(M\)-length frequency array, and \(N\)-length frequency derivative array).
\[\chi^2 = f_{0} \cdot (t-T_{\mathrm{ref}}) + \frac{1}{2} \cdot f_{1} \cdot (t-T_{\mathrm{ref}})^2 + \frac{1}{6} \cdot f_{2} \cdot (t-T_{\mathrm{ref}})^3 + \cdots,\]where \(T_{\mathrm{ref}\) is the reference time, \(f_{0}\), \(f_{1}\), \(f_{2}\), …, are the parameters of pulsar.
- Parameters
data : array-like
The time array of photons to calculate the chisquare
f : array-like
A set of frequencies to calculate the chisquare for event array
F1 : float, optional, default 0
The frequency derivative
pepoch : float
The reference time of pulsar timing parameters
nbins : int
The number of bins to fold profile
F2 : float, optional, default 0
The second frequency derivative
F3 : float, optional, default 0
The third frequency derivative
F4 : float, optional, default 0
The forth frequency derivative
- Returns
chi_square : array-like
An \(M \times N\) array, as \(M\) is the length of frequency f, \(N\) is the length of frequency derivative F1
- tatpulsar.utils.functions.cal_chisquare(data, f, pepoch, nbins, F1=0, F2=0, F3=0, F4=0, parallel=False)[source]
Calculate the Pearson-Chisquare value for given spinning parameters at given epoch time.
\[\chi^2 = f_{0} \cdot (t-T_{\mathrm{ref}}) + \frac{1}{2} \cdot f_{1} \cdot (t-T_{\mathrm{ref}})^2 + \frac{1}{6} \cdot f_{2} \cdot (t-T_{\mathrm{ref}})^3 + \cdots,\]where \(T_{\mathrm{ref}\) is the reference time, \(f_{0}\), \(f_{1}\), \(f_{2}\), …, are the parameters of pulsar.
- Parameters
data : array-like
The time array of photons to calculate the chisquare
f : array-like
A set of frequencies to calculate the chisquare for event array
pepoch : float
The reference time of pulsar timing parameters
nbins : int
The number of bins to fold profile
F1 : float, optional, default 0
The frequency derivative
F2 : float, optional, default 0
The second frequency derivative
F3 : float, optional, default 0
The third frequency derivative
F4 : float, optional, default 0
The forth frequency derivative
parallel : boolean, optional, default
False
whether to use multi-core CPU to calculate the chisquare
- Returns
chi_square : array-like
The calculated \(\chi^2\) array
- tatpulsar.utils.functions.cal_event_gti(data, tgap=1)[source]
calculate the gti edges of given event data. if the time gap between two adjacent event is larger than tgap, it split the event into two intervals. Otherwise, we take the event as continous observation.
- data: array-like
the event array
- tgap: float
the critical time gap to split GTI
- Returns
gtis: ndarray
the list of GTI array, example [[gti0_0, gti0_1], [gti1_0, gti1_1], …]
- tatpulsar.utils.functions.ccf(f1, f2)[source]
Calculate the cross-correlation function for given data. f1 is the original signal f2 is probe signal(shift and test)
- Returns
y : array-like
the ccf function distribution
delay : float
the index of the delay between the input data f2 and the f1
- tatpulsar.utils.functions.get_parameters(kwargs)[source]
get the parameters for searching The format of input could be a name of parfile, a dictionary, or standard python function arguments.
- tatpulsar.utils.functions.met2mjd(data, telescope='fermi')[source]
Convert Mission Elapse Time (MET) to Modified Julian Date (MJD).
\[T_{\mathrm{MJD}} = T_{\mathrm{MET}}/86400 + \mathrm{MJDREF},\]where MJDREF is the reference time for each mission.
- Parameters
data : float
The MET time
telescope : str, default ‘fermi’
The name of the mission, support mission are {‘fermi’, ‘hxmt’, ‘nicer’, ‘gecam’, ‘nustar’, ‘ixpe’}
- Returns
mjd : float
The MJD time
- tatpulsar.utils.functions.mjd2met(data, telescope='fermi')[source]
Convert Modified Julian Date (MJD) to Mission Elapse Time (MET)
\[T_{\mathrm{MJD}} = T_{\mathrm{MET}}/86400 + \mathrm{MJDREF},\]where MJDREF is the reference time for each mission.
- Parameters
data : float
The MJD time
telescope : str, default ‘fermi’
The name of the mission, support mission are {‘fermi’, ‘hxmt’, ‘nicer’, ‘gecam’, ‘nustar’, ‘ixpe’}
- Returns
met : float
The MET time
- tatpulsar.utils.functions.print_loop_percentage(iterator_i, total, printstr='')[source]
print the percentage in a loop
Functions for GTI manipulation
- tatpulsar.utils.gti.create_gti_fits(gti_template, outfile, tstart, tstop)[source]
write tstart and tstop array to a FITS file, based on the GTI FITS template
Parameters:
- gti_templateFITS file
a GTI fits file as template
- outfilefilename
the name of the output FITS file
- tstartarray-like
array of GTIs start time
- tstoparray-like
array of GTIs stop time
- tatpulsar.utils.gti.create_gti_txt(outfile, tstart, tstop)[source]
create a text file storing the TSTART and TSTOP value, it’s normally used for hxmtscreen
- Parameters
outfile : filename
the name of the output FITS file
- tstartarray-like
array of GTIs start time
- tstoparray-like
array of GTIs stop time
- tatpulsar.utils.gti.gti_intersection(gti1, gti2)[source]
find the intersection for given gti1 and gti2.
- Parameters
gti1: list
the list of input GTI, the format of gti list is a 2D list: [[start, stop], [start, stop], … ,[start, stop]]
gti2: list
the list of input GTI, the format of gti list is a 2D list: [[start, stop], [start, stop], … ,[start, stop]]
- Returns
gti: list
the list of output GTI, the format of gti list is a 2D list: [[start, stop], [start, stop], … ,[start, stop]]
- tatpulsar.utils.gti.gti_union(gti1, gti2)[source]
find the intersection for given gti1 and gti2.
- Parameters
gti1: list
the list of input GTI, the format of gti list is a 2D list: [[start, stop], [start, stop], … ,[start, stop]]
gti2: list
the list of input GTI, the format of gti list is a 2D list: [[start, stop], [start, stop], … ,[start, stop]]
- Returns
gti: list
the list of output GTI, the format of gti list is a 2D list: [[start, stop], [start, stop], … ,[start, stop]]
Class to read the TEMPO2 parameter file (.par)
- class tatpulsar.utils.readpar.readpar(filepath)[source]
Class to parse the TEMPO2 parameter file (.par) The parameters in the ‘.par’ file might be capitalized but those parameters stored in this object are case INSENSITIVE (see examples below).
Examples
>>> eph = readpar('test.par') >>> print("F0 = ", eph.F0.value, eph.f0.value) >>> F0 = 29.636679699921209437 >>> print("F0 error = ", eph.F0.error) >>> F0 error = 1.7247236495e-09 >>> print("PEPOCH = ", eph.PEPOCH.value) >>> PEPOCH = 58066.18087539147382 >>> print("PEPOCH error = ", eph.PEPOCH.error) >>> PEPOCH error = None