Exemple #1
0
def main(argv=None):
    parser = argparse.ArgumentParser(description="PINT tool for command-line barycentering calculations.")

    parser.add_argument("time",help="MJD (UTC, by default)")
    parser.add_argument("--timescale",default="utc",
        help="Time scale for MJD argument ('utc', 'tt', 'tdb'), default=utc")
    parser.add_argument("--format",
        help="Format for time argument ('mjd' or any astropy.Time format (e.g. 'isot'), see <http://docs.astropy.org/en/stable/time/#time-format>)",
        default="mjd")
    parser.add_argument("--freq",type=float,default=np.inf,
        help="Frequency to use, MHz")
    parser.add_argument("--obs",default="Geocenter",
        help="Observatory code (default = Geocenter)")
    parser.add_argument("--parfile",help="par file to read model from",default=None)
    parser.add_argument("--ra",
        help="RA to use (e.g. '12h22m33.2s' if not read from par file)")
    parser.add_argument("--dec",
        help="Decl. to use (e.g. '19d21m44.2s' if not read from par file)")
    parser.add_argument("--dm",
        help="DM to use (if not read from par file)",type=float,default=0.0)
    parser.add_argument("--ephem",default="DE421",help="Ephemeris to use")
    parser.add_argument("--use_gps",default=False,action='store_true',help="Apply GPS to UTC clock corrections")
    parser.add_argument("--use_bipm",default=False,action='store_true',help="Use TT(BIPM) instead of TT(TAI)")


    args = parser.parse_args(argv)

    if args.format in ("mjd","jd", "unix"):
        # These formats require conversion from string to longdouble first
        fmt = args.format
        # Never allow format == 'mjd' because it fails when scale is 'utc'
        # Change 'mjd' to 'pulsar_mjd' to deal with this.
        if fmt == "mjd":
            fmt = "pulsar_mjd"
        t = Time(np.longdouble(args.time),scale=args.timescale,format=fmt,
            precision=9)
        print(t)
    else:
        t = Time(args.time,scale=args.timescale,format=args.format, precision=9)
    log.debug(t.iso)

    t = toa.TOA(t,freq=args.freq,obs=args.obs)
    # Build TOAs and compute TDBs and positions from ephemeris
    ts = toa.get_TOAs_list([t],ephem=args.ephem, include_bipm=args.use_bipm,
        include_gps=args.use_gps, planets=False)

    if args.parfile is not None:
        m=pint.models.get_model(args.parfile)
    else:
        # Construct model by hand
        m=pint.models.StandardTimingModel
        # Should check if 12:13:14.2 syntax is used and support that as well!
        m.RAJ.quantity = Angle(args.ra)
        m.DECJ.quantity = Angle(args.dec)
        m.DM.quantity = args.dm*u.parsec/u.cm**3

    tdbtimes = m.get_barycentric_toas(ts)

    print("{0:.16f}".format(tdbtimes[0].value))
    return
def make_toa_list(time, obs, freq, **other_meta):
    """Convert times to a list of PINT TOAs.

    The input timestamps will be flattened to a 1-dimensional array.

    Parameters
    ----------
    time : `~astropy.time.Time`
        Input timestamps.
    obs : str
        The observatory code or names
    freq : float or `~astropy.units.Quantity`
        The observing frequency, the default unit is MHz

    Returns
    -------
    List of `~pint.toa.TOA`.
    """
    toa_list = []
    for t_stamp in time.ravel():
        # This format converting should be done by PINT in the future.
        if t_stamp.scale == 'utc':
            t_stamp = t_stamp.replicate(format='pulsar_mjd')
        toa_entry = toa.TOA(t_stamp, obs=obs, freq=freq, **other_meta)
        toa_list.append(toa_entry)
    return toa_list
Exemple #3
0
def test_generate_polycos(tmpdir, par_file, obs, obsfreq, nspan, ncoeff):
    output_polyco = tmpdir / "B1855_polyco_round_trip_from_par.dat"
    mjd_start, mjd_end = 55000.0, 55001.0

    model = get_model(str(par_file))

    p = Polycos()
    p.generate_polycos(model, mjd_start, mjd_end, obs, nspan, ncoeff, obsfreq)
    p.write_polyco_file(output_polyco)
    q = Polycos()
    q.read_polyco_file(output_polyco)

    mjds = np.linspace(mjd_start, mjd_end, 51)

    t = toa.get_TOAs_list(
        [toa.TOA(mjd, obs=obs, freq=obsfreq) for mjd in mjds],
        ephem=model.EPHEM.value)
    ph1 = p.eval_abs_phase(mjds)
    ph2 = q.eval_abs_phase(mjds)
    ph3 = model.phase(t, abs_phase=True)

    assert np.allclose(ph1.int.value[0], ph3.int.value[0])
    assert np.allclose(ph1.frac.value[0], ph3.frac.value[0])
    assert np.allclose(ph2.int.value[0], ph3.int.value[0])
    assert np.allclose(ph2.frac.value[0], ph3.frac.value[0])
Exemple #4
0
    def __call__(self, time):
        """Create list of TOAs for one or more times.

        Parameters
        ----------
        time : `~astropy.time.Time`
            Input time stamps.

        Returns
        -------
        toas : `~pint.toa.TOAs`
            Combining all TOAs.
        """
        # local import since we cannot count on PINT being present,
        # and doing it globally messes up sphinx.
        from pint import toa

        if time.scale == 'utc':
            time = time.replicate(format='pulsar_mjd')

        freq, _ = np.broadcast_arrays(self.frequency, time.jd1, subok=True)
        time = time._apply(np.broadcast_to, freq.shape)
        toa_list = []
        for t, f in zip(time.ravel(), freq.ravel()):
            # This format converting should be done by PINT in the future.
            toa_entry = toa.TOA(t, obs=self.observatory, freq=f)
            toa_list.append(toa_entry)

        toas = toa.get_TOAs_list(toa_list, **self.control_params)
        toas.shape = time.shape
        return toas
Exemple #5
0
    def get_TZR_toa(self, toas):
        """Get the TOAs class for the TZRMJD.

        We are treating the TZRMJD as a special TOA.
        Note that any observatory clock corrections will be applied
        to this TOA, as with any other TOA. This does not affect the
        value of the TZRMJD parmeter, however.
        """

        if self.tz_cache is not None:
            # This should also verify that the clkc_info is the same so it is valid.
            return self.tz_cache
        else:
            # NOTE: Using TZRMJD.quantity.jd[1,2] so that the time scale can be properly
            # set to the TZRSITE default timescale (e.g. UTC for TopoObs and TDB for SSB)
            TZR_toa = toa.TOA(
                (self.TZRMJD.quantity.jd1 - 2400000.5,
                 self.TZRMJD.quantity.jd2),
                obs=self.TZRSITE.value,
                freq=self.TZRFRQ.quantity,
            )
            clkc_info = toas.clock_corr_info
            tz = toa.get_TOAs_list(
                [TZR_toa],
                include_bipm=clkc_info["include_bipm"],
                include_gps=clkc_info["include_gps"],
                ephem=toas.ephem,
                planets=toas.planets,
            )
            self.tz_cache = tz
            return tz
Exemple #6
0
def get_barycentric_correction(orbfile, parfile, dt=5, ephem='DE421'):
    no = SatelliteObs(name="NuSTAR", FPorbname=orbfile, tt2tdb_mode="pint")
    with fits.open(orbfile) as hdul:
        mjdref = high_precision_keyword_read(hdul[1].header, 'MJDREF')

    mjds = np.arange(no.X.x[1], no.X.x[-2], dt / 86400)
    mets = (mjds - mjdref) * 86400

    obs, scale = 'nustar', "tt"
    toalist = [None] * len(mjds)

    for i in range(len(mjds)):
        # Create TOA list
        toalist[i] = toa.TOA(mjds[i], obs=obs, scale=scale)

    if parfile is not None and os.path.exists(parfile):
        modelin = pint.models.get_model(parfile)
    else:
        modelin = get_dummy_parfile_for_position(orbfile)

    ts = toa.get_TOAs_list(
        toalist,
        ephem=ephem,
        include_bipm=False,
        include_gps=False,
        planets=False,
        tdb_method='default',
    )
    bats = modelin.get_barycentric_toas(ts)
    return interp1d(mets, (bats.value - mjds) * 86400,
                    assume_sorted=True,
                    bounds_error=False,
                    fill_value='extrapolate',
                    kind='quadratic')
Exemple #7
0
def load_Polar_TOAs(eventname,weightcolumn="",minunix=0.0, maxunix=2500000000.0):
    '''
    TOAlist = load_Polar_TOAs(eventname)
      Read photon event times out of a POLAR event root file and return
      a list of PINT TOA objects.

      weightcolumn specifies ROOT Branch name to read the photon weights
      from. 

    '''
    import ROOT
    # Load photon times from POLAR event ROOT file
    f = ROOT.TFile.Open(eventname,"read")
    t = f.Get("t_trigger")

    # Read time Branch from ROOT file
    filter = '(tunix > '+str(minunix)+') && (tunix < '+str(maxunix)+')'
    if weightcolumn == "":
        array = tree2array(t,branches=['tunix'],selection=filter,start=0,stop=10000,step=1)
        ut=array['tunix']
    else:
        array = tree2array(t,branches=['tunix',weightcolumn],selection=filter,start=0,stop=10000,step=1)
        ut=array['tunix']
        weights=array[weightcolumn]
    if len(ut) == 0:
        log.error('No unix time read from file!')
        raise

    e=50*u.keV
    mjd1=Time(ut[0],format="unix")
    mjd2=Time(ut[-1],format="unix")
    log.info('Building spacecraft local TOAs, with MJDs in range {0} to {1}'.format(
            mjd1.tt.mjd,mjd2.tt.mjd))
    polarobs = get_observatory('Polar')

    try:
        if weightcolumn == "":
            toalist=[toa.TOA(Time(m,format="unix").tt.mjd, obs='Polar', scale='tt',energy=e)
                            for m in ut]
        else:
            toalist=[toa.TOA(Time(m,format="unix").tt.mjd, obs='Polar', scale='tt',energy=e,weight=w)
                            for m,w in zip(ut,weights)]
    except KeyError:
        log.error('Error processing POLAR TOAs. You may have forgotten to specify an ENG file file with --eng')
        raise

    return toalist
Exemple #8
0
def load_event_TOAs(eventname, mission, weights=None):
    '''
    Read photon event times out of a FITS file as PINT TOA objects.

    Correctly handles raw event files, or ones processed with axBary to have
    barycentered  TOAs. Different conditions may apply to different missions.

    Parameters
    ----------
    eventname : str
        File name of the FITS event list
    mission : str
        Name of the mission (e.g. RXTE, XMM)
    weights : array or None
        The array has to be of the same size as the event list. Overwrites
        possible weight lists from mission-specific FITS files

    Returns
    -------
    toalist : list of TOA objects
    '''
    # Load photon times from event file
    hdulist = pyfits.open(eventname)

    extension = mission_config[mission]["fits_extension"]

    if hdulist[1].name not in extension.split(','):
        raise RuntimeError(
            'First table in FITS file' +
            'must be {}. Found {}'.format(extension, hdulist[1].name))

    timesys, timeref = _get_timesys_and_timeref(hdulist[1])

    if not mission_config[mission]['allow_local'] \
            and timesys != 'TDB':
        log.error('Raw spacecraft TOAs not yet supported for ' + mission)

    obs, scale = _default_obs_and_scale(mission, timesys, timeref)

    # Read time column from FITS file
    mjds = read_fits_event_mjds_tuples(hdulist[1])

    new_kwargs = _get_columns_from_fits(
        hdulist[1], mission_config[mission]["fits_columns"])

    hdulist.close()

    if weights is not None:
        new_kwargs["weights"] = weights

    toalist = [None] * len(mjds)
    kw = {}
    for i in range(len(mjds)):
        # Create TOA list
        for key in new_kwargs.keys():
            kw[key] = new_kwargs[key][i]
        toalist[i] = toa.TOA(mjds[i], obs=obs, scale=scale, **kw)

    return toalist
Exemple #9
0
    def _gen_polyco(self, parfile, MJD_start, segLength = 60.0, ncoeff = 15, \
                       maxha=12.0, method="TEMPO", numNodes=20, usePINT = True):
        """
        This will be a convenience function to generate polycos and subsequent parameters to replace the values in a 
        PSRFITS file header. The default way to do this will be to use PINT (usePINT = True), with other methods
        currently unsupported. The input values are:
        parfile [string] : path to par file used to generate the polycos. The observing frequency, and observatory will 
                            come from the par file
        MJD_start [float] : Start MJD of the polyco. Should start no later than the beginning of the observation
        segLength [float] : Length in minutes of the range covered by the polycos generated. Default is 60 minutes
        ncoeff [int] : number of polyco coeffeicients to generate. Default is 15, the same as in the PSRFITS file
        maxha [float] : max hour angle needed by PINT. Default is 12.0
        method [string] : Method PINT uses to generate the polyco. Currently only TEMPO is supported.
        numNodes [int] : Number of nodes PINT will use to fit the polycos. Must be larger than ncoeff
        usePINT [bool] : Method used to generate polycos. Currently only PINT is supported.
        """
        if usePINT:
            # Define dictionary to put parameters into
            polyco_dict = {'NSPAN': segLength, 'NCOEF': ncoeff}
            # load parfile to PINT model object
            m = models.get_model(parfile)
            # Determine MJD_end based on segLength
            MJD_end = MJD_start + np.double(
                make_quant(segLength, 'min').to('day').value)  # MJD
            # Determine obseratory and observing frequency
            obsFreq = m.TZRFRQ.value  # in MHz
            polyco_dict['REF_FREQ'] = obsFreq
            obs = m.TZRSITE.value  # string
            polyco_dict['NSITE'] = obs.encode(
                'utf-8')  # observatory code needs to be in binary
            # Get pulsar frequency
            polyco_dict['REF_F0'] = m.F0.value
            # get the polycos
            pcs = polycos.Polycos()
            pcs.generate_polycos(m, MJD_start, MJD_end, obs, segLength, ncoeff, obsFreq, maxha=12.0, method="TEMPO", \
                                     numNodes=20)
            coeffs = pcs.polycoTable['entry'][-1].coeffs
            polyco_dict['COEFF'] = coeffs
            # Now we need to determine the reference MJD, and phase
            REF_MJD = np.double(pcs.polycoTable['tmid'][-1])
            polyco_dict['REF_MJD'] = REF_MJD
            # Now find the phase difference
            tmid_toa = toa.get_TOAs_list(
                [toa.TOA(REF_MJD, obs=obs, freq=obsFreq)])
            ref_phase = m.phase(tmid_toa)
            # Need to force positive value
            if ref_phase.frac.value[0] < 0.0:
                ref_frac_phase = 1.0 - abs(ref_phase.frac.value[0])
            else:
                ref_frac_phase = ref_phase.frac.value[0]
            polyco_dict['REF_PHS'] = ref_frac_phase

            return polyco_dict

        else:
            print("Only PINT is currently supported for generating polycos")
            raise NotImplementedError()
Exemple #10
0
 def test_roundtrip_bary_toa_Tempo2format(self):
     # Create a barycentric TOA
     t1time = Time(58534.0, 0.0928602471130208, format="mjd", scale="tdb")
     t1 = toa.TOA(t1time, obs="Barycenter", freq=0.0)
     ts = toa.get_TOAs_list([t1], ephem="DE421")
     ts.write_TOA_file("testbary.tim", format="Tempo2")
     ts2 = toa.get_TOAs("testbary.tim")
     print(ts.table, ts2.table)
     assert np.abs(ts.table["mjd"][0] - ts2.table["mjd"][0]) < 1.0e-15 * u.d
     assert np.abs(ts.table["tdb"][0] - ts2.table["tdb"][0]) < 1.0e-15 * u.d
Exemple #11
0
def test_roundtrip_topo_toa_TEMPOformat(tmpdir):
    # Create a barycentric TOA
    t1time = Time(58534.0, 0.0928602471130208, format="mjd", scale="utc")
    t1 = toa.TOA(t1time, obs="gbt", freq=0.0)
    ts = toa.get_TOAs_list([t1], ephem="DE421")
    outnm = os.path.join(tmpdir, "testtopot1.tim")
    ts.write_TOA_file(outnm, format="TEMPO")
    ts2 = toa.get_TOAs(outnm)
    assert np.abs(ts.table["mjd"][0] - ts2.table["mjd"][0]) < 1.0e-15 * u.d
    assert np.abs(ts.table["tdb"][0] - ts2.table["tdb"][0]) < 1.0e-15 * u.d
Exemple #12
0
 def test_roundtrip_topo_toa_TEMPOformat(self):
     # Create a barycentric TOA
     t1time = Time(58534.0, 0.0928602471130208, format="mjd", scale="utc")
     t1 = toa.TOA(t1time, obs="gbt", freq=0.0)
     ts = toa.get_TOAs_list([t1], ephem="DE421")
     ts.write_TOA_file("testtopot1.tim", format="TEMPO")
     ts2 = toa.get_TOAs("testtopot1.tim")
     print(ts.table, ts2.table)
     assert ts.table["mjd"][0] - ts2.table["mjd"][0] < 1.0e-15
     assert ts.table["tdb"][0] - ts2.table["tdb"][0] < 1.0e-15
Exemple #13
0
def _load_and_prepare_TOAs(mjds, ephem="DE405"):
    toalist = [None] * len(mjds)
    for i, m in enumerate(mjds):
        toalist[i] = toa.TOA(m, obs='Barycenter', scale='tdb')

    toalist = toa.TOAs(toalist=toalist)
    if 'tdb' not in toalist.table.colnames:
        toalist.compute_TDBs()
    if 'ssb_obs_pos' not in toalist.table.colnames:
        toalist.compute_posvels(ephem, False)
    return toalist
Exemple #14
0
def get_orbital_correction_from_ephemeris_file(mjdstart, mjdstop, parfile,
                                               ntimes=1000, ephem="DE405"):
    """Get a correction for orbital motion from pulsar parameter file.

    Parameters
    ----------
    mjdstart, mjdstop : float
        Start and end of the time interval where we want the orbital solution
    parfile : str
        Any parameter file understood by PINT (Tempo or Tempo2 format)

    Other parameters
    ----------------
    ntimes : int
        Number of time intervals to use for interpolation. Default 1000

    Returns
    -------
    correction_sec : function
        Function that accepts in input an array of times in seconds and a
        floating-point MJDref value, and returns the deorbited times
    correction_mjd : function
        Function that accepts times in MJDs and returns the deorbited times.
    """
    from scipy.interpolate import interp1d
    simon("Assuming events are already referred to the solar system "
          "barycenter (timescale is TDB)")
    if not HAS_PINT:
        raise ImportError("You need the optional dependency PINT to use this "
                          "functionality: github.com/nanograv/pint")

    mjds = np.linspace(mjdstart, mjdstop, ntimes)
    toalist = [None] * len(mjds)
    for i, m in enumerate(mjds):
        toalist[i] = toa.TOA(m, obs='Barycenter', scale='tdb')

    toalist = toa.TOAs(toalist = toalist)
    if 'tdb' not in toalist.table.colnames:
        toalist.compute_TDBs()
    if 'ssb_obs_pos' not in toalist.table.colnames:
        toalist.compute_posvels(ephem, False)

    m = pint.models.get_model(parfile)
    tdbtimes = m.get_barycentric_toas(toalist.table)

    correction_mjd = interp1d(mjds, tdbtimes)

    def correction_sec(times, mjdref):
        deorb_mjds = correction_mjd(times / 86400 + mjdref)
        return np.array((deorb_mjds - mjdref) * 86400)

    return correction_sec, correction_mjd
Exemple #15
0
def test_commenting_toas(tmpdir):
    # Create a barycentric TOA
    t1time = Time(58534.0, 0.0928602471130208, format="mjd", scale="utc")
    t1 = toa.TOA(t1time, obs="gbt", freq=0.0)
    t2 = toa.TOA(t1time, obs="gbt", freq=0.0)
    ts = toa.get_TOAs_list([t1, t2], ephem="DE421")
    assert ts.ntoas == 2
    outnm = os.path.join(tmpdir, "testtopo.tim")
    ts.write_TOA_file(outnm)
    ts2 = toa.get_TOAs(outnm)
    assert ts2.ntoas == 2  # none should be commented
    ts.table[0]["flags"]["cut"] = "do_not_like"  # cut flag
    ts.table[1]["flags"]["ignore"] = str(1)  # ignore flag
    ts.write_TOA_file(outnm)
    ts3 = toa.get_TOAs(outnm)  # defaut is to not comment
    assert ts3.ntoas == 2  # none should be commented by default
    ts.write_TOA_file(outnm, commentflag="cut")
    ts4 = toa.get_TOAs(outnm)
    assert ts4.ntoas == 1  # one should be commented
    ts.write_TOA_file(outnm, commentflag="ignore")
    ts5 = toa.get_TOAs(outnm)
    assert ts5.ntoas == 1  # one should be commented
Exemple #16
0
def _load_and_prepare_TOAs(mjds, errs_us=None, ephem="DE405"):
    errs_us = assign_value_if_none(errs_us, np.zeros_like(mjds))

    toalist = [None] * len(mjds)
    for i, m in enumerate(mjds):
        toalist[i] = \
            toa.TOA(m, error=errs_us[i], obs='Barycenter', scale='tdb')

    toalist = toa.TOAs(toalist=toalist)
    if 'tdb' not in toalist.table.colnames:
        toalist.compute_TDBs()
    if 'ssb_obs_pos' not in toalist.table.colnames:
        toalist.compute_posvels(ephem, False)
    return toalist
Exemple #17
0
 def get_TZR_toa(self, toas):
     """ Get the TOAs class for the TZRMJD. We are treating the TZRMJD as a
         special TOA.
     """
     TZR_toa = toa.TOA(self.TZRMJD.quantity,
                       obs=self.TZRSITE.value,
                       freq=self.TZRFRQ.quantity)
     clkc_info = toas.clock_corr_info
     tz = toa.get_TOAs_list([
         TZR_toa,
     ],
                            include_bipm=clkc_info['include_bipm'],
                            include_gps=clkc_info['include_gps'],
                            ephem=toas.ephem,
                            planets=toas.planets)
     return tz
Exemple #18
0
    def get_TZR_toa(self, toas):
        """Get the TOAs class for the TZRMJD.

        We are treating the TZRMJD as a special TOA.
        Note that any observatory clock corrections will be applied
        to this TOA, as with any other TOA. This does not affect the
        value of the TZRMJD parmeter, however.
        """
        clkc_info = toas.clock_corr_info
        # If we have cached the TZR TOA and all the TZR* and clock info has not changed, then don't rebuild it
        if self.tz_cache is not None:
            if (self.tz_clkc_info["include_bipm"] == clkc_info["include_bipm"]
                    and self.tz_clkc_info["include_gps"]
                    == clkc_info["include_gps"]
                    and self.tz_planets == toas.planets
                    and self.tz_ephem == toas.ephem and self.tz_hash == hash(
                        (self.TZRMJD.value, self.TZRSITE.value,
                         self.TZRFRQ.value))):
                return self.tz_cache
        # Otherwise we have to build the TOA and apply clock corrections
        # NOTE: Using TZRMJD.quantity.jd[1,2] so that the time scale can be properly
        # set to the TZRSITE default timescale (e.g. UTC for TopoObs and TDB for SSB)
        log.debug(
            "Creating and dealing with the single TZR_toa for absolute phase")
        TZR_toa = toa.TOA(
            (self.TZRMJD.quantity.jd1 - 2400000.5, self.TZRMJD.quantity.jd2),
            obs=self.TZRSITE.value,
            freq=self.TZRFRQ.quantity,
        )
        tz = toa.get_TOAs_list(
            [TZR_toa],
            include_bipm=clkc_info["include_bipm"],
            include_gps=clkc_info["include_gps"],
            ephem=toas.ephem,
            planets=toas.planets,
        )
        log.debug("Done with TZR_toa")
        self.tz_cache = tz
        self.tz_hash = hash(
            (self.TZRMJD.value, self.TZRSITE.value, self.TZRFRQ.value))
        self.tz_clkc_info = clkc_info
        self.tz_planets = toas.planets
        self.tz_ephem = toas.ephem
        return tz
Exemple #19
0
 def get_TZR_toa(self, toas):
     """ Get the TOAs class for the TZRMJD. We are treating the TZRMJD as a
         special TOA.
     """
     # NOTE: Using TZRMJD.quantity.jd[1,2] so that the time scale can be properly
     # set to the TZRSITE default timescale (e.g. UTC for TopoObs and TDB for SSB)
     TZR_toa = toa.TOA(
         (self.TZRMJD.quantity.jd1 - 2400000.5, self.TZRMJD.quantity.jd2),
         obs=self.TZRSITE.value,
         freq=self.TZRFRQ.quantity)
     clkc_info = toas.clock_corr_info
     tz = toa.get_TOAs_list([
         TZR_toa,
     ],
                            include_bipm=clkc_info['include_bipm'],
                            include_gps=clkc_info['include_gps'],
                            ephem=toas.ephem,
                            planets=toas.planets)
     return tz
Exemple #20
0
    def change_pepoch(self, new_epoch, toas=None, delay=None):
        """Move PEPOCH to a new time and change the related paramters.

        Parameters
        ----------
        new_epoch: float or `astropy.Time` object
            The new PEPOCH value.
        toas: `toa` object, optional.
            If current PEPOCH is not provided, the first pulsar frame toa will
            be treated as PEPOCH.
        delay: `numpy.array` object
            If current PEPOCH is not provided, it is required for computing the
            first pulsar frame toa.
        """
        if isinstance(new_epoch, Time):
            new_epoch = Time(new_epoch, scale='tdb', precision=9)
        else:
            new_epoch = Time(new_epoch, scale='tdb', format='mjd', precision=9)
        # make new_epoch a toa for delay calculation.
        new_epoch_toa = toa.get_TOAs_list([
            toa.TOA(new_epoch),
        ],
                                          ephem=toas.ephem)

        if self.PEPOCH.value is None:
            if toas is None or delay is None:
                raise ValueError(
                    "`PEPOCH` is not in the model, thus, 'toa' and"
                    " 'delay' shoule be givne.")
            tbl = toas.table
            phsepoch_ld = time_to_longdouble(tbl['tdb'][0] - delay[0])
        else:
            phsepoch_ld = time_to_longdouble(self.PEPOCH.quantity)
        dt = ((time_to_longdouble(new_epoch) - phsepoch_ld) * u.day)
        fterms = [0.0 * u.Unit("")] + self.get_spin_terms()
        # rescale the fterms
        for n in range(len(fterms) - 1):
            f_par = getattr(self, 'F{}'.format(n))
            f_par.value = taylor_horner_deriv(dt.to(u.second),
                                              fterms,
                                              deriv_order=n + 1)
        self.PEPOCH.value = new_epoch
Exemple #21
0
    def test_roundtrip_topo_toa_TEMPOformat(self):
        # Create a barycentric TOA
        t1time = Time(58534.0, 0.0928602471130208, format="mjd", scale="utc")
        t1 = toa.TOA(t1time, obs="gbt", freq=0.0)
        ts = toa.get_TOAs_list([t1], ephem="DE421")
        ts.write_TOA_file("testtopot1.tim", format="TEMPO")
        ts2 = toa.get_TOAs("testtopot1.tim")
        print(ts.table, ts2.table)
        assert np.abs(ts.table["mjd"][0] - ts2.table["mjd"][0]) < 1.0e-15 * u.d
        assert np.abs(ts.table["tdb"][0] - ts2.table["tdb"][0]) < 1.0e-15 * u.d

        # Comment out because TEMPO2 distro doesn't include gmrt2gps.clk so far
        # def test_roundtrip_gmrt_toa_Tempo2format(self):
        #     if os.getenv("TEMPO2") is None:
        #         pytest.skip("TEMPO2 evnironment variable is not set, can't run this test")
        #     # Create a barycentric TOA
        #     t1time = Time(58534.0, 0.0928602471130208, format="mjd", scale="utc")
        #     t1 = toa.TOA(t1time, obs="gmrt", freq=0.0)
        #     ts = toa.get_TOAs_list([t1], ephem="DE421")
        #     ts.write_TOA_file("testgmrt.tim", format="Tempo2")
        #     ts2 = toa.get_TOAs("testgmrt.tim")
        #     print(ts.table, ts2.table)
        assert np.abs(ts.table["mjd"][0] - ts2.table["mjd"][0]) < 1.0e-15 * u.d
        assert np.abs(ts.table["tdb"][0] - ts2.table["tdb"][0]) < 1.0e-15 * u.d
Exemple #22
0
    def generate_polycos(
        self,
        model,
        mjdStart,
        mjdEnd,
        obs,
        segLength,
        ncoeff,
        obsFreq,
        maxha=12.0,
        method="TEMPO",
        numNodes=20,
    ):
        """
        Generate the polyco data.

        Parameters
        ----------
        model : TimingModel
            TimingModel to generate the Polycos with parameters
            setup.

        mjdStart : float / numpy longdouble
            Start time of polycos in mjd

        mjdEnd : float / numpy longdouble
            Ending time of polycos in mjd

        obs : str
            Observatory code

        segLength : int
            Length of polyco segement [minutes]

        ncoeff : int
            Number of coefficents

        obsFreq : float
            Observing frequency [MHz]

        maxha : float optional. Default 12.0
            Maximum hour angle. Only 12.0 is supported for now.

        method : string optional ["TEMPO", "TEMPO2", ...] Default TEMPO
            Method to generate polycos. Only the TEMPO method is supported for now.

        numNodes : int optional. Default 20
            Number of nodes for fitting. It cannot be less then the number of
            coefficents.

        Return
        ---------
        A polyco table.
        """
        mjdStart = data2longdouble(mjdStart)
        mjdEnd = data2longdouble(mjdEnd)
        segLength = int(segLength)
        obsFreq = float(obsFreq)

        # Use the planetary ephemeris specified in the model, if available.
        if model.EPHEM.value is not None:
            ephem = model.EPHEM.value
        else:
            log.info(
                "No ephemeris specified in model, using DE421 to generate polycos"
            )
            ephem = "DE421"

        if maxha != 12.0:
            raise ValueError("Maximum hour angle != 12.0 is not supported.")

        # Make sure the number of nodes is bigger than number of coeffs.
        if numNodes < ncoeff:
            numNodes = ncoeff + 1

        mjdSpan = data2longdouble(segLength / MIN_PER_DAY)
        # Generate "nice" MJDs for consistency with what tempo2 does
        tmids = np.arange(
            int(mjdStart * 24) * 60,
            int(mjdEnd * 24) * 60 + segLength, segLength)
        tmids = data2longdouble(tmids) / MIN_PER_DAY

        # generate the ploynomial coefficents
        if method == "TEMPO":
            entryList = []
            # Using tempo1 method to create polycos
            # If you want to disable the progress bar, add disable=True to the tqdm() call.
            for tmid in tqdm(tmids):
                tStart = tmid - mjdSpan / 2
                tStop = tmid + mjdSpan / 2
                nodes = np.linspace(tStart, tStop, numNodes)

                toaMid = toa.get_TOAs_list(
                    [
                        toa.TOA((np.modf(tmid)[1], np.modf(tmid)[0]),
                                obs=obs,
                                freq=obsFreq)
                    ],
                    ephem=ephem,
                )

                refPhase = model.phase(toaMid, abs_phase=True)

                # Create node toas(Time sample using TOA class)
                toaList = [
                    toa.TOA(
                        (np.modf(toaNode)[1], np.modf(toaNode)[0]),
                        obs=obs,
                        freq=obsFreq,
                    ) for toaNode in nodes
                ]

                toas = toa.get_TOAs_list(toaList, ephem=ephem)

                ph = model.phase(toas, abs_phase=True)
                dt = (nodes - tmid) * MIN_PER_DAY
                rdcPhase = ph - refPhase
                rdcPhase = rdcPhase.int - (dt * model.F0.value *
                                           60.0) + rdcPhase.frac
                dtd = dt.astype(float)  # Truncate to double
                rdcPhased = rdcPhase.astype(float)
                coeffs = np.polyfit(dtd, rdcPhased, ncoeff - 1)[::-1]

                date, hms = Time(tmid, format="mjd", scale="utc").iso.split()
                yy, mm, dd = date.split("-")
                date = dd + "-" + MONTHS[int(mm) - 1] + "-" + yy[-2:]
                hms = float(hms.replace(":", ""))

                entry = PolycoEntry(
                    tmid,
                    segLength,
                    refPhase.int,
                    refPhase.frac,
                    model.F0.value,
                    ncoeff,
                    coeffs,
                )

                entry_dict = OrderedDict()
                entry_dict["psr"] = model.PSR.value
                entry_dict["date"] = date
                entry_dict["utc"] = hms
                entry_dict["tmid"] = tmid
                entry_dict["dm"] = model.DM.value
                entry_dict["doppler"] = 0.0
                entry_dict["logrms"] = 0.0
                entry_dict["mjd_span"] = segLength
                entry_dict["t_start"] = entry.tstart
                entry_dict["t_stop"] = entry.tstop
                entry_dict["obs"] = obs
                entry_dict["obsfreq"] = obsFreq

                if model.is_binary:
                    binphase = model.orbital_phase(toaMid, radians=False)[0]
                    entry_dict["binary_phase"] = binphase
                    b = model.get_components_by_category()["pulsar_system"][0]
                    entry_dict["f_orbit"] = 1 / b.PB.value

                entry_dict["entry"] = entry
                entryList.append(entry_dict)

            pTable = table.Table(entryList, meta={"name": "Polyco Data Table"})

            self.polycoTable = pTable
            if len(self.polycoTable) == 0:
                raise ValueError("Zero polycos found for table")

        else:
            raise NotImplementedError(
                "Only TEMPO method has been implemented.")
Exemple #23
0
from astropy.visualization import quantity_support

quantity_support()

# %%
# Here is how to create a single TOA in Python
# The first argument is an MJD(UTC) as a 2-double tuple to allow extended precision
# and the second argument is the TOA uncertainty
# Wherever possible, it is good to use astropy units on the values,
# but there are sensible defaults if you leave them out (us for uncertainty, MHz for freq)
import pint.toa as toa

a = toa.TOA(
    (54567, 0.876876876876876),
    4.5 * u.us,
    freq=1400.0 * u.MHz,
    obs="GBT",
    backend="GUPPI",
    name="guppi_56789.fits",
)
print(a)

# %%
# An example of reading a TOA file
import pint.toa as toa

t = toa.get_TOAs("NGC6440E.tim", usepickle=False)

# %%
#  You can print a summary of the loaded TOAs
t.print_summary()
Exemple #24
0
def main(argv=None):

    parser = argparse.ArgumentParser(description="Use PINT to compute H-test and plot Phaseogram from a Fermi FT1 event file.")
    parser.add_argument("eventfile",help="Fermi event FITS file name.")
    parser.add_argument("parfile",help="par file to construct model from")
    parser.add_argument("weightcol",help="Column name for event weights (or 'CALC' to compute them)")
    parser.add_argument("--ft2",help="Path to FT2 file.",default=None)
    parser.add_argument("--addphase",help="Write FT1 file with added phase column",
        default=False,action='store_true')
    parser.add_argument("--plot",help="Show phaseogram plot.", action='store_true', default=False)
    parser.add_argument("--plotfile",help="Output figure file name (default=None)", default=None)
    parser.add_argument("--maxMJD",help="Maximum MJD to include in analysis", default=None)
    parser.add_argument("--outfile",help="Output figure file name (default is to overwrite input file)", default=None)
    parser.add_argument("--planets",help="Use planetary Shapiro delay in calculations (default=False)", default=False, action="store_true")
    parser.add_argument("--ephem",help="Planetary ephemeris to use (default=DE421)", default="DE421")
    args = parser.parse_args(argv)

    # If outfile is specified, that implies addphase
    if args.outfile is not None:
        args.addphase = True

    # Read in model
    modelin = pint.models.get_model(args.parfile)
    if 'ELONG' in modelin.params:
        tc = SkyCoord(modelin.ELONG.quantity,modelin.ELAT.quantity,
            frame='barycentrictrueecliptic')
    else:
        tc = SkyCoord(modelin.RAJ.quantity,modelin.DECJ.quantity,frame='icrs')

    # Read TZR parameters from parfile separately
    tzrmjd = None
    tzrsite = '@'
    tzrfrq = np.inf*u.MHz
    for line in open(args.parfile):
        if line.startswith('TZRMJD'):
            tzrmjd = np.longdouble(line.split()[1])
        if line.startswith('TZRSITE'):
            tzrsite = line.split()[1]
        if line.startswith('TZRFRQ'):
            tzrfrq = np.float(line.split()[1])*u.MHz

    if tzrmjd is None:
        tzrmjd = tl[0].mjd

    tztoa = toa.TOA(tzrmjd,obs=tzrsite,freq=tzrfrq)
    tz = toa.get_TOAs_list([tztoa],include_bipm=False,include_gps=True,
        ephem=args.ephem, planets=True)

    if args.ft2 is not None:
        # Instantiate FermiObs once so it gets added to the observatory registry
        FermiObs(name='Fermi',ft2name=args.ft2)

    # Read event file and return list of TOA objects
    tl  = load_Fermi_TOAs(args.eventfile, weightcolumn=args.weightcol,
                          targetcoord=tc)

    # Discard events outside of MJD range
    if args.maxMJD is not None:
        tlnew = []
        print("pre len : ",len(tl))
        maxT = Time(float(args.maxMJD),format='mjd')
        print("maxT : ",maxT)
        for tt in tl:
            if tt.mjd < maxT:
                tlnew.append(tt)
        tl=tlnew
        print("post len : ",len(tlnew))

    # Now convert to TOAs object and compute TDBs and posvels
    # For Fermi, we are not including GPS or TT(BIPM) corrections
    ts = toa.get_TOAs_list(tl,include_gps=False,include_bipm=False,planets=args.planets,ephem=args.ephem)
    ts.filename = args.eventfile

    print(ts.get_summary())
    mjds = ts.get_mjds()
    print(mjds.min(),mjds.max())

    # Compute model phase for each TOA
    # Subtracts off model phase at TZRMJD so that point defines phase 0.0
    iphss,phss = modelin.phase(ts.table) - modelin.phase(tz.table)
    # ensure all postive
    phases = np.where(phss < 0.0 * u.cycle, phss + 1.0 * u.cycle, phss)
    mjds = ts.get_mjds()
    weights = np.array([w['weight'] for w in ts.table['flags']])
    h = float(hmw(phases,weights))
    print("Htest : {0:.2f} ({1:.2f} sigma)".format(h,h2sig(h)))
    if args.plot:
        log.info("Making phaseogram plot with {0} photons".format(len(mjds)))
        phaseogram(mjds,phases,weights,bins=100,plotfile = args.plotfile)

    if args.addphase:
        # Read input FITS file (again).
        # If overwriting, open in 'update' mode
        if args.outfile is None:
            hdulist = pyfits.open(args.eventfile,mode='update')
        else:
            hdulist = pyfits.open(args.eventfile)
        event_hdu = hdulist[1]
        event_hdr=event_hdu.header
        event_dat=event_hdu.data
        if len(event_dat) != len(phases):
            raise RuntimeError('Mismatch between length of FITS table ({0}) and length of phase array ({1})!'.format(len(event_dat),len(phases)))
        if 'PULSE_PHASE' in event_hdu.columns.names:
            log.info('Found existing PULSE_PHASE column, overwriting...')
            # Overwrite values in existing Column
            event_dat['PULSE_PHASE'] = phases
        else:
            # Construct and append new column, preserving HDU header and name
            log.info('Adding new PULSE_PHASE column.')
            phasecol = pyfits.ColDefs([pyfits.Column(name='PULSE_PHASE', format='D',
                array=phases)])
            bt = pyfits.BinTableHDU.from_columns( event_hdu.columns + phasecol,
                header=event_hdr,name=event_hdu.name)
            hdulist[1] = bt
        if args.outfile is None:
            # Overwrite the existing file
            log.info('Overwriting existing FITS file '+args.eventfile)
            hdulist.flush(verbose=True, output_verify='warn')
        else:
            # Write to new output file
            log.info('Writing output FITS file '+args.outfile)
            hdulist.writeto(args.outfile,overwrite=True, checksum=True, output_verify='warn')

    return 0
Exemple #25
0
def main(argv=None):
    import argparse
    parser = argparse.ArgumentParser(
        description="PINT tool for simulating TOAs")
    parser.add_argument("parfile", help="par file to read model from")
    parser.add_argument("timfile", help="Output TOA file name")
    parser.add_argument("--startMJD",
                        help="MJD of first fake TOA (default=56000.0)",
                        type=float,
                        default=56000.0)
    parser.add_argument("--ntoa",
                        help="Number of fake TOAs to generate",
                        type=int,
                        default=100)
    parser.add_argument("--duration",
                        help="Span of TOAs to generate (days)",
                        type=int,
                        default=400)
    parser.add_argument("--obs",
                        help="Observatory code (default: GBT)",
                        default="GBT")
    parser.add_argument("--freq",
                        help="Frequency for TOAs (MHz) (default: 1400)",
                        type=float,
                        default=1400.0)
    parser.add_argument(
        "--error",
        help="Random error to apply to each TOA (us, default=1.0)",
        type=float,
        default=1.0)
    parser.add_argument("--plot",
                        help="Plot residuals",
                        action="store_true",
                        default=False)
    parser.add_argument("--ephem", help="Ephemeris to use", default="DE421")
    parser.add_argument("--planets",
                        help="Use planetary Shapiro delay",
                        action="store_true",
                        default=False)
    args = parser.parse_args(argv)

    log.info("Reading model from {0}".format(args.parfile))
    m = pint.models.get_model(args.parfile)

    duration = args.duration * u.day
    start = Time(args.startMJD, scale='utc', format='mjd', precision=9)
    error = args.error * u.microsecond
    freq = args.freq * u.MHz
    scale = 'utc'

    times = np.linspace(0, duration.to(u.day).value, args.ntoa) * u.day + start

    tl = [
        toa.TOA(t, error=error, obs=args.obs, freq=freq, scale=scale)
        for t in times
    ]

    ts = toa.TOAs(toalist=tl)

    # WARNING! I'm not sure how clock corrections should be handled here!
    # Do we apply them, or not?
    if not any(['clkcorr' in f for f in ts.table['flags']]):
        log.info("Applying clock corrections.")
        ts.apply_clock_corrections()
    if 'tdb' not in ts.table.colnames:
        log.info("Getting IERS params and computing TDBs.")
        ts.compute_TDBs()
    if 'ssb_obs_pos' not in ts.table.colnames:
        log.info("Computing observatory positions and velocities.")
        ts.compute_posvels(args.ephem, args.planets)

    # F_local has units of Hz; discard cycles unit in phase to get a unit
    # that TimeDelta understands
    F_local = m.d_phase_d_toa(ts)
    rs = m.phase(ts.table).frac.value / F_local

    # Adjust the TOA times to put them where their residuals will be 0.0
    ts.adjust_TOAs(TimeDelta(-1.0 * rs))
    rspost = m.phase(ts.table).frac.value / F_local

    # Do a second iteration
    ts.adjust_TOAs(TimeDelta(-1.0 * rspost))

    # Write TOAs to a file
    ts.write_TOA_file(args.timfile, name='fake', format='Tempo2')

    if args.plot:
        # This should be a very boring plot with all residuals flat at 0.0!
        import matplotlib.pyplot as plt
        rspost2 = m.phase(ts.table).frac / F_local
        plt.errorbar(ts.get_mjds(),
                     rspost2.to(u.us).value,
                     yerr=ts.get_errors().to(u.us).value)
        newts = pint.toa.get_TOAs(args.timfile)
        rsnew = m.phase(newts.table).frac / F_local
        plt.errorbar(newts.get_mjds(),
                     rsnew.to(u.us).value,
                     yerr=newts.get_errors().to(u.us).value)
        #plt.plot(ts.get_mjds(),rspost.to(u.us),'x')
        plt.xlabel('MJD')
        plt.ylabel('Residual (us)')
        plt.show()
Exemple #26
0
def load_Fermi_TOAs(ft1name,
                    weightcolumn=None,
                    targetcoord=None,
                    logeref=4.1,
                    logesig=0.5,
                    minweight=0.0,
                    minmjd=0.0,
                    maxmjd=np.inf):
    '''
    TOAlist = load_Fermi_TOAs(ft1name)
      Read photon event times out of a Fermi FT1 file and return
      a list of PINT TOA objects.
      Correctly handles raw FT1 files, or ones processed with gtbary
      to have barycentered or geocentered TOAs.

      weightcolumn specifies the FITS column name to read the photon weights
      from.  The special value 'CALC' causes the weights to be computed empirically
      as in Philippe Bruel's SearchPulsation code.
      logeref and logesig are parameters for the weight computation and are only
      used when weightcolumn='CALC'.

      When weights are loaded, or computed, events are filtered by weight >= minweight
    '''
    import astropy.io.fits as pyfits
    # Load photon times from FT1 file
    hdulist = pyfits.open(ft1name)
    ft1hdr = hdulist[1].header
    ft1dat = hdulist[1].data

    # TIMESYS will be 'TT' for unmodified Fermi LAT events (or geocentered), and
    #                 'TDB' for events barycentered with gtbary
    # TIMEREF will be 'GEOCENTER' for geocentered events,
    #                 'SOLARSYSTEM' for barycentered,
    #             and 'LOCAL' for unmodified events

    timesys = ft1hdr['TIMESYS']
    log.info("TIMESYS {0}".format(timesys))
    timeref = ft1hdr['TIMEREF']
    log.info("TIMEREF {0}".format(timeref))

    # Read time column from FITS file
    mjds = read_fits_event_mjds_tuples(hdulist[1])
    if len(mjds) == 0:
        log.error('No MJDs read from file!')
        raise

    energies = ft1dat.field('ENERGY') * u.MeV
    if weightcolumn is not None:
        if weightcolumn == 'CALC':
            photoncoords = SkyCoord(ft1dat.field('RA') * u.degree,
                                    ft1dat.field('DEC') * u.degree,
                                    frame='icrs')
            weights = calc_lat_weights(ft1dat.field('ENERGY'),
                                       photoncoords.separation(targetcoord),
                                       logeref=logeref,
                                       logesig=logesig)
        else:
            weights = ft1dat.field(weightcolumn)
        if minweight > 0.0:
            idx = np.where(weights > minweight)[0]
            mjds = mjds[idx]
            energies = energies[idx]
            weights = weights[idx]

    # limit the TOAs to ones in selected MJD range
    mjds_float = np.array([r[0] + r[1] for r in mjds])
    idx = np.logical_and((mjds_float > minmjd), (mjds_float < maxmjd))
    mjds = mjds[idx]
    energies = energies[idx]
    if weightcolumn is not None:
        weights = weights[idx]

    if timesys == 'TDB':
        log.info("Building barycentered TOAs")
        if weightcolumn is None:
            toalist = [
                toa.TOA(m,
                        obs='Barycenter',
                        scale='tdb',
                        energy=e,
                        error=1.0 * u.us) for m, e in zip(mjds, energies)
            ]
        else:
            toalist = [
                toa.TOA(m,
                        obs='Barycenter',
                        scale='tdb',
                        energy=e,
                        weight=w,
                        error=1.0 * u.us)
                for m, e, w in zip(mjds, energies, weights)
            ]
    else:
        if timeref == 'LOCAL':
            log.info(
                'Building spacecraft local TOAs, with MJDs in range {0} to {1}'
                .format(mjds[0], mjds[-1]))
            assert timesys == 'TT'
            try:
                fermiobs = get_observatory('Fermi')
            except KeyError:
                log.error(
                    'Fermi observatory not defined. Make sure you have specified an FT2 file!'
                )
                raise

            try:
                if weightcolumn is None:
                    toalist = [
                        toa.TOA(m,
                                obs='Fermi',
                                scale='tt',
                                energy=e,
                                error=1.0 * u.us)
                        for m, e in zip(mjds, energies)
                    ]
                else:
                    toalist = [
                        toa.TOA(m,
                                obs='Fermi',
                                scale='tt',
                                energy=e,
                                weight=w,
                                error=1.0 * u.us)
                        for m, e, w in zip(mjds, energies, weights)
                    ]
            except KeyError:
                log.error(
                    'Error processing Fermi TOAs. You may have forgotten to specify an FT2 file with --ft2'
                )
                raise
        else:
            log.info("Building geocentered TOAs")
            if weightcolumn is None:
                toalist = [
                    toa.TOA(m, obs='Geocenter', scale='tt', energy=e)
                    for m, e in zip(mjds, energies)
                ]
            else:
                toalist = [
                    toa.TOA(m, obs='Geocenter', scale='tt', energy=e, weight=w)
                    for m, e, w in zip(mjds, energies, weights)
                ]

    return toalist
from astropy.visualization import quantity_support

quantity_support()

# %%
# Here is how to create a single TOA in Python
# The first argument is an MJD(UTC) as a 2-double tuple to allow extended precision
# and the second argument is the TOA uncertainty
# Wherever possible, it is good to use astropy units on the values,
# but there are sensible defaults if you leave them out (us for uncertainty, MHz for freq)
import pint.toa as toa

a = toa.TOA(
    (54567, 0.876876876876876),
    4.5 * u.us,
    freq=1400.0 * u.MHz,
    obs="GBT",
    backend="GUPPI",
)
print(a)

# %%
# An example of reading a TOA file
import pint.toa as toa

t = toa.get_TOAs("NGC6440E.tim", usepickle=False)

# %%
#  You can print a summary of the loaded TOAs
t.print_summary()
Exemple #28
0
def load_fits_TOAs(
    eventname,
    mission,
    weights=None,
    extension=None,
    timesys=None,
    timeref=None,
    minmjd=-np.inf,
    maxmjd=np.inf,
):
    """
    Read photon event times out of a FITS file as PINT TOA objects.

    Correctly handles raw event files, or ones processed with axBary to have
    barycentered TOAs. Different conditions may apply to different missions.

    The minmjd/maxmjd parameters can be used to avoid instantiation of TOAs
    we don't want, which can otherwise be very slow.

    Parameters
    ----------
    eventname : str
        File name of the FITS event list
    mission : str
        Name of the mission (e.g. RXTE, XMM)
    weights : array or None
        The array has to be of the same size as the event list. Overwrites
        possible weight lists from mission-specific FITS files
    extension : str
        FITS extension to read
    timesys : str, default None
        Force this time system
    timeref : str, default None
        Forse this time reference
    minmjd : float, default "-infinity"
        minimum MJD timestamp to return
    maxmjd : float, default "infinity"
        maximum MJD timestamp to return

    Returns
    -------
    toalist : list of TOA objects
    """
    # Load photon times from event file
    hdulist = pyfits.open(eventname)
    if mission not in mission_config:
        log.warning("Mission not recognized. Using generic")
        mission = "generic"

    if (extension is not None and isinstance(extension, str)
            and hdulist[1].name not in extension.split(",")):
        raise RuntimeError(
            "First table in FITS file" +
            "must be {}. Found {}".format(extension, hdulist[1].name))
    if isinstance(extension, int) and extension != 1:
        raise ValueError(
            "At the moment, only data in the first FITS extension is supported"
        )

    if timesys is None:
        timesys = _get_timesys(hdulist[1])
    if timeref is None:
        timeref = _get_timeref(hdulist[1])
    check_timesys(timesys)
    check_timeref(timeref)

    if not mission_config[mission]["allow_local"] and timesys != "TDB":
        log.error("Raw spacecraft TOAs not yet supported for " + mission)

    obs, scale = _default_obs_and_scale(mission, timesys, timeref)

    # Read time column from FITS file
    mjds = read_fits_event_mjds_tuples(hdulist[1])

    new_kwargs = _get_columns_from_fits(
        hdulist[1], mission_config[mission]["fits_columns"])

    hdulist.close()

    if weights is not None:
        new_kwargs["weights"] = weights

    # mask out times/columns outside of mjd range
    mjds_float = np.asarray([r[0] + r[1] for r in mjds])
    idx = (minmjd < mjds_float) & (mjds_float < maxmjd)
    mjds = mjds[idx]
    for key in new_kwargs.keys():
        new_kwargs[key] = new_kwargs[key][idx]

    toalist = [None] * len(mjds)
    kw = {}
    for i in range(len(mjds)):
        # Create TOA list
        for key in new_kwargs.keys():
            kw[key] = new_kwargs[key][i]
        toalist[i] = toa.TOA(mjds[i], obs=obs, scale=scale, **kw)

    return toalist
Exemple #29
0
    def generate_polycos(self,
                         model,
                         mjdStart,
                         mjdEnd,
                         obs,
                         segLength,
                         ncoeff,
                         obsFreq,
                         maxha,
                         method="TEMPO",
                         numNodes=20):
        """
        Generate the polyco file data file. 

        Parameters
        ---------
        model : TimingModel
            TimingModel for generate the Polycos with parameters
            setup.

        mjdStart : float / nump longdouble
            Start time of polycos in mjd

        mjdEnd : float / nump longdouble
            Ending time of polycos in mjd

        obs : str
            Observatory code

        segLength : 
            Length of polyco segement [unit: minutes]

        ncoeff : 
            number of coefficents

        obsFreq : 
            observing frequency

        maxha :
            Maximum hour angle

        method : sting optional ['TEMPO','TEMPO2',...] Default TEMPO
            Method to generate polycos. Now it is only support the TEMPO method. 

        numNodes : int optional. Default 20
            Number of nodes for fitting. It can not be less then the number of 
            coefficents.
        Return 
        ---------

        A polyco table. 


        """
        mjdStart = np.longdouble(mjdStart) * u.day
        mjdEnd = np.longdouble(mjdEnd) * u.day
        timeLength = mjdEnd - mjdStart
        segLength = np.longdouble(segLength) * u.min
        obsFreq = float(obsFreq)
        month = [
            'Jan', 'Feb', 'Mar', 'Apr', 'May', 'Jun', 'Jul', 'Aug', 'Sep',
            'Oct', 'Nov', 'Dec'
        ]
        # Alocate memery
        coeffs = np.longdouble(np.zeros(ncoeff))
        entryList = []
        entryIntvl = np.arange(mjdStart.value, mjdEnd.value,
                               segLength.to('day').value)
        if entryIntvl[-1] < mjdEnd.value:
            entryIntvl = np.append(entryIntvl, mjdEnd.value)

        # Make sure the number of nodes is bigger then number of coeffs.
        if numNodes < ncoeff:
            numNodes = ncoeff + 1

    # generate the ploynomial coefficents
        if method == "TEMPO":
            # Using tempo1 method to create polycos
            for i in range(len(entryIntvl) - 1):
                tStart = entryIntvl[i]
                tStop = entryIntvl[i + 1]
                nodes = np.linspace(tStart, tStop, numNodes)
                tmid = ((tStart + tStop) / 2.0) * u.day
                toaMid = toa.get_TOAs_list([
                    toa.TOA((np.modf(tmid.value)[1], np.modf(tmid.value)[0]),
                            obs=obs,
                            freq=obsFreq),
                ])
                refPhase = model.phase(toaMid.table)
                mjdSpan = ((tStop - tStart) * u.day).to('min')
                # Create node toas(Time sample using TOA class)
                toaList = [
                    toa.TOA((np.modf(toaNode)[1], np.modf(toaNode)[0]),
                            obs=obs,
                            freq=obsFreq) for toaNode in nodes
                ]

                toas = toa.get_TOAs_list(toaList)

                ph = model.phase(toas.table)
                dt = (nodes * u.day - tmid).to('min')  # Use constant
                rdcPhase = ph - refPhase
                rdcPhase = rdcPhase.int - dt.value * model.F0.value * 60.0 + rdcPhase.frac
                dtd = dt.value.astype(float)  # Trancate to double
                rdcPhased = rdcPhase.astype(float)
                coeffs = np.polyfit(dtd, rdcPhased, ncoeff - 1)
                coeffs = coeffs[::-1]
                midTime = at.Time(int(tmid.value),
                                  np.modf(tmid.value)[0],
                                  format='mjd',
                                  scale='utc')
                date, hms = midTime.iso.split()
                yy, mm, dd = date.split('-')
                date = dd + '-' + month[int(mm) - 1] + '-' + yy[2:4]
                hms = hms.replace(':', "")
                entry = polycoEntry(tmid.value,
                                    mjdSpan.to('day').value, refPhase.int,
                                    refPhase.frac, model.F0.value, ncoeff,
                                    coeffs, obs)
                entryList.append(
                    (model.PSR.value, date, hms, tmid.value, model.DM.value,
                     0.0, 0.0, 0.0, obsFreq, entry))

            pTable = table.Table(rows=entryList,
                                 names=('psr', 'date', 'utc', 'tmid', 'dm',
                                        'dopper', 'logrms', 'binary_phase',
                                        'obsfreq', 'entry'),
                                 meta={'name': 'Ployco Data Table'})
            self.polycoTable = pTable

        else:
            #  Reading from an old polycofile
            pass
Exemple #30
0
    def generate_polycos(
        self,
        model,
        mjdStart,
        mjdEnd,
        obs,
        segLength,
        ncoeff,
        obsFreq,
        maxha=12.0,
        method="TEMPO",
        numNodes=20,
    ):
        """
        Generate the polyco data.

        Parameters
        ---------
        model : TimingModel
            TimingModel to generate the Polycos with parameters
            setup.

        mjdStart : float / nump longdouble
            Start time of polycos in mjd

        mjdEnd : float / nump longdouble
            Ending time of polycos in mjd

        obs : str
            Observatory code

        segLength : float
            Length of polyco segement [unit: minutes]

        ncoeff : int
            number of coefficents

        obsFreq : float
            observing frequency [unit: MHz]

        maxha : float optional. Default 12.0
            Maximum hour angle

        method : string optional ['TEMPO','TEMPO2',...] Default TEMPO
            Method to generate polycos. Only the TEMPO method is supported for now.

        numNodes : int optional. Default 20
            Number of nodes for fitting. It cannot be less then the number of
            coefficents.

        Return
        ---------
        A polyco table.


        """
        mjdStart = data2longdouble(mjdStart) * u.day
        mjdEnd = data2longdouble(mjdEnd) * u.day
        timeLength = mjdEnd - mjdStart
        segLength = data2longdouble(segLength) * u.min
        obsFreq = float(obsFreq)
        month = [
            "Jan",
            "Feb",
            "Mar",
            "Apr",
            "May",
            "Jun",
            "Jul",
            "Aug",
            "Sep",
            "Oct",
            "Nov",
            "Dec",
        ]
        # Alocate memery
        coeffs = data2longdouble(np.zeros(ncoeff))
        entryList = []
        entryIntvl = np.arange(mjdStart.value, mjdEnd.value,
                               segLength.to("day").value)
        if entryIntvl[-1] < mjdEnd.value:
            entryIntvl = np.append(entryIntvl, mjdEnd.value)

        # Make sure the number of nodes is bigger then number of coeffs.
        if numNodes < ncoeff:
            numNodes = ncoeff + 1

        # generate the ploynomial coefficents
        if method == "TEMPO":
            # Using tempo1 method to create polycos
            for i in range(len(entryIntvl) - 1):
                tStart = entryIntvl[i]
                tStop = entryIntvl[i + 1]
                nodes = np.linspace(tStart, tStop, numNodes)
                tmid = ((tStart + tStop) / 2.0) * u.day
                toaMid = toa.get_TOAs_list([
                    toa.TOA(
                        (np.modf(tmid.value)[1], np.modf(tmid.value)[0]),
                        obs=obs,
                        freq=obsFreq,
                    )
                ])
                refPhase = model.phase(toaMid)
                mjdSpan = ((tStop - tStart) * u.day).to("min")
                # Create node toas(Time sample using TOA class)
                toaList = [
                    toa.TOA(
                        (np.modf(toaNode)[1], np.modf(toaNode)[0]),
                        obs=obs,
                        freq=obsFreq,
                    ) for toaNode in nodes
                ]

                toas = toa.get_TOAs_list(toaList)

                ph = model.phase(toas)
                dt = (nodes * u.day - tmid).to("min")  # Use constant
                rdcPhase = ph - refPhase
                rdcPhase = (rdcPhase.int -
                            (dt.value * model.F0.value * 60.0) * u.cycle +
                            rdcPhase.frac)
                dtd = dt.value.astype(float)  # Truncate to double
                rdcPhased = rdcPhase.astype(float)
                coeffs = np.polyfit(dtd, rdcPhased, ncoeff - 1)
                coeffs = coeffs[::-1]
                midTime = Time(int(tmid.value),
                               np.modf(tmid.value)[0],
                               format="mjd",
                               scale="utc")
                date, hms = midTime.iso.split()
                yy, mm, dd = date.split("-")
                date = dd + "-" + month[int(mm) - 1] + "-" + yy[2:4]
                hms = hms.replace(":", "")
                entry = PolycoEntry(
                    tmid.value,
                    mjdSpan.to("day").value,
                    refPhase.int,
                    refPhase.frac,
                    model.F0.value,
                    ncoeff,
                    coeffs,
                    obs,
                )
                entryList.append((
                    model.PSR.value,
                    date,
                    hms,
                    tmid.value,
                    model.DM.value,
                    0.0,
                    0.0,
                    0.0,
                    mjdSpan.to("day").value,
                    tStart,
                    tStop,
                    obs,
                    obsFreq,
                    entry,
                ))

            pTable = table.Table(
                rows=entryList,
                names=(
                    "psr",
                    "date",
                    "utc",
                    "tmid",
                    "dm",
                    "doppler",
                    "logrms",
                    "binary_phase",
                    "mjd_span",
                    "t_start",
                    "t_stop",
                    "obs",
                    "obsfreq",
                    "entry",
                ),
                meta={"name": "Polyco Data Table"},
            )
            self.polycoTable = pTable
            if len(self.polycoTable) == 0:
                raise ValueError("Zero polycos found for table")
        else:
            #  Reading from an old polycofile
            pass