def load_nustar_orbit(orb_filename): """Load data from a NuSTAR orbit file Parameters ---------- orb_filename : str Name of file to load Returns ------- astropy.table.Table containing Time, x, y, z, v_x, v_y, v_z data """ # Load photon times from FT1 file if "_orb" in orb_filename: log.warning("The NuSTAR orbit file you are providing is known to give" "a solution precise only to the ~0.5ms level. Use the " "pipeline-produced attitude-orbit file ('*.attorb.gz') for" "better precision.") hdulist = pyfits.open(orb_filename) orb_hdr = hdulist[1].header orb_dat = hdulist[1].data log.info("Opened orb FITS file {0}".format(orb_filename)) # TIMESYS should be 'TT' # TIMEREF should be 'LOCAL', since no delays are applied timesys = orb_hdr["TIMESYS"] log.info("orb TIMESYS {0}".format(timesys)) try: timeref = orb_hdr["TIMEREF"] except KeyError: timeref = "LOCAL" log.info("orb TIMEREF {0}".format(timeref)) # The X, Y, Z position are for the START time mjds_TT = read_fits_event_mjds(hdulist[1]) mjds_TT = mjds_TT * u.d # SC_POS is in meters in X,Y,Z Earth-centered Inertial (ECI) coordinates SC_POS = orb_dat.field("POSITION") X = SC_POS[:, 0] * u.km Y = SC_POS[:, 1] * u.km Z = SC_POS[:, 2] * u.km SC_VEL = orb_dat.field("VELOCITY") Vx = SC_VEL[:, 0] * u.km / u.s Vy = SC_VEL[:, 1] * u.km / u.s Vz = SC_VEL[:, 2] * u.km / u.s log.info("Building orb table covering MJDs {0} to {1}".format( mjds_TT.min(), mjds_TT.max())) orb_table = Table( [mjds_TT, X, Y, Z, Vx, Vy, Vz], names=("MJD_TT", "X", "Y", "Z", "Vx", "Vy", "Vz"), meta={"name": "orb"}, ) return orb_table
def load_FPorbit(orbit_filename): """Load data from an (RXTE or NICER) FPorbit file Reads a FPorbit FITS file Parameters ---------- orbit_filename : str Name of file to load Returns ------- astropy Table containing Time, x, y, z, v_x, v_y, v_z data """ # Load orbit FITS file hdulist = pyfits.open(orbit_filename) # log.info('orb file HDU name is {0}'.format(hdulist[1].name)) if hdulist[1].name not in ("ORBIT", "XTE_PE"): log.error( "NICER orb file first extension is {0}. It should be ORBIT".format( hdulist[1].name ) ) FPorbit_hdr = hdulist[1].header FPorbit_dat = hdulist[1].data log.info("Opened FPorbit FITS file {0}".format(orbit_filename)) # TIMESYS should be 'TT' # TIMEREF should be 'LOCAL', since no delays are applied timesys = FPorbit_hdr["TIMESYS"] log.debug("FPorbit TIMESYS {0}".format(timesys)) timeref = FPorbit_hdr["TIMEREF"] log.debug("FPorbit TIMEREF {0}".format(timeref)) mjds_TT = read_fits_event_mjds(hdulist[1]) mjds_TT = mjds_TT * u.d log.debug("FPorbit spacing is {0}".format((mjds_TT[1] - mjds_TT[0]).to(u.s))) X = FPorbit_dat.field("X") * u.m Y = FPorbit_dat.field("Y") * u.m Z = FPorbit_dat.field("Z") * u.m Vx = FPorbit_dat.field("Vx") * u.m / u.s Vy = FPorbit_dat.field("Vy") * u.m / u.s Vz = FPorbit_dat.field("Vz") * u.m / u.s log.info( "Building FPorbit table covering MJDs {0} to {1}".format( mjds_TT.min(), mjds_TT.max() ) ) FPorbit_table = Table( [mjds_TT, X, Y, Z, Vx, Vy, Vz], names=("MJD_TT", "X", "Y", "Z", "Vx", "Vy", "Vz"), meta={"name": "FPorbit"}, ) # Make sure table is sorted by time log.debug("Sorting FPorbit table") FPorbit_table.sort("MJD_TT") # Now delete any bad entries where the positions are 0.0 idx = np.where( np.logical_and(FPorbit_table["X"] != 0.0, FPorbit_table["Y"] != 0.0) )[0] if len(idx) != len(FPorbit_table): log.warning( "Dropping {0} zero entries from FPorbit table".format( len(FPorbit_table) - len(idx) ) ) FPorbit_table = FPorbit_table[idx] return FPorbit_table
def load_Fermi_FT2(ft2_filename): """Load data from a Fermi FT2 file The contents of the FT2 file are described here: https://fermi.gsfc.nasa.gov/ssc/data/analysis/documentation/Cicerone/Cicerone_Data/LAT_Data_Columns.html#SpacecraftFile The coordinates are X, Y, Z in the ECI (Earth-Centered Inertial) frame. I (@paulray) **believe** this is the same as astropy's GCRS <http://docs.astropy.org/en/stable/api/astropy.coordinates.GCRS.html>, but this should be confirmed. Parameters ---------- ft2_filename : str Name of file to load Returns ------- astropy Table containing Time, x, y, z, v_x, v_y, v_z data """ # Load photon times from FT1 file hdulist = pyfits.open(ft2_filename) FT2_hdr = hdulist[1].header FT2_dat = hdulist[1].data log.info("Opened FT2 FITS file {0}".format(ft2_filename)) # TIMESYS should be 'TT' # TIMEREF should be 'LOCAL', since no delays are applied timesys = FT2_hdr["TIMESYS"] log.info("FT2 TIMESYS {0}".format(timesys)) timeref = FT2_hdr["TIMEREF"] log.info("FT2 TIMEREF {0}".format(timeref)) # The X, Y, Z position are for the START time mjds_TT = read_fits_event_mjds(hdulist[1], timecolumn="START") mjds_TT = mjds_TT * u.d # SC_POS is in meters in X,Y,Z Earth-centered Inertial (ECI) coordinates SC_POS = FT2_dat.field("SC_POSITION") X = SC_POS[:, 0] * u.m Y = SC_POS[:, 1] * u.m Z = SC_POS[:, 2] * u.m try: # If available, get the velocities from the FT2 file SC_VEL = FT2_dat.field("SC_VELOCITY") Vx = SC_VEL[:, 0] * u.m / u.s Vy = SC_VEL[:, 1] * u.m / u.s Vz = SC_VEL[:, 2] * u.m / u.s except: # Otherwise, compute velocities by differentiation because FT2 does not have velocities # This is not the best way. Should fit an orbit and determine velocity from that. dt = mjds_TT[1] - mjds_TT[0] log.info("FT2 spacing is " + str(dt.to(u.s))) # Use "spacing" argument for gradient to handle nonuniform entries tt = mjds_TT.to(u.s).value Vx = np.gradient(X.value, tt) * u.m / u.s Vy = np.gradient(Y.value, tt) * u.m / u.s Vz = np.gradient(Z.value, tt) * u.m / u.s log.info( "Building FT2 table covering MJDs {0} to {1}".format( mjds_TT.min(), mjds_TT.max() ) ) FT2_table = Table( [mjds_TT, X, Y, Z, Vx, Vy, Vz], names=("MJD_TT", "X", "Y", "Z", "Vx", "Vy", "Vz"), meta={"name": "FT2"}, ) return FT2_table
def main(argv=None): import argparse parser = argparse.ArgumentParser( description= "Use PINT to compute event phases and make plots of photon event files." ) parser.add_argument( "eventfile", help= "Photon event FITS file name (e.g. from NICER, RXTE, XMM, Chandra).", ) parser.add_argument("parfile", help="par file to construct model from") parser.add_argument("--orbfile", help="Name of orbit file", default=None) parser.add_argument("--maxMJD", help="Maximum MJD to include in analysis", default=None) parser.add_argument("--minMJD", help="Minimum MJD to include in analysis", default=None) parser.add_argument("--plotfile", help="Output figure file name (default=None)", default=None) parser.add_argument( "--addphase", help="Write FITS file with added phase column", default=False, action="store_true", ) parser.add_argument( "--addorbphase", help="Write FITS file with added orbital phase column", default=False, action="store_true", ) parser.add_argument( "--absphase", help="Write FITS file with integral portion of pulse phase (ABS_PHASE)", default=False, action="store_true", ) parser.add_argument( "--barytime", help= "Write FITS file with a column containing the barycentric time as double precision MJD.", default=False, action="store_true", ) parser.add_argument( "--outfile", help="Output FITS file name (default=same as eventfile)", default=None, ) parser.add_argument("--ephem", help="Planetary ephemeris to use (default=DE421)", default="DE421") parser.add_argument( "--tdbmethod", help="Method for computing TT to TDB (default=astropy)", default="default", ) parser.add_argument("--plot", help="Show phaseogram plot.", action="store_true", default=False) 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)", ) # parser.add_argument("--fix",help="Apply 1.0 second offset for NICER", action='store_true', default=False) args = parser.parse_args(argv) # If outfile is specified, that implies addphase if args.outfile is not None: args.addphase = True # If plotfile is specified, that implies plot if args.plotfile is not None: args.plot = True # set MJD ranges maxmjd = np.inf if (args.maxMJD is None) else float(args.maxMJD) minmjd = 0.0 if (args.minMJD is None) else float(args.minMJD) # Read event file header to figure out what instrument is is from hdr = pyfits.getheader(args.eventfile, ext=1) log.info("Event file TELESCOPE = {0}, INSTRUMENT = {1}".format( hdr["TELESCOP"], hdr["INSTRUME"])) if hdr["TELESCOP"] == "NICER": # Instantiate NICER observatory once so it gets added to the observatory registry if args.orbfile is not None: log.info("Setting up NICER observatory") get_satellite_observatory("NICER", args.orbfile) # Read event file and return list of TOA objects try: tl = load_NICER_TOAs(args.eventfile, minmjd=minmjd, maxmjd=maxmjd) except KeyError: log.error( "Observatory not recognized. This probably means you need to provide an orbit file or barycenter the event file." ) sys.exit(1) elif hdr["TELESCOP"] == "XTE": # Instantiate RXTE observatory once so it gets added to the observatory registry if args.orbfile is not None: # Determine what observatory type is. log.info("Setting up RXTE observatory") get_satellite_observatory("RXTE", args.orbfile) # Read event file and return list of TOA objects tl = load_RXTE_TOAs(args.eventfile, minmjd=minmjd, maxmjd=maxmjd) elif hdr["TELESCOP"].startswith("XMM"): # Not loading orbit file here, since that is not yet supported. tl = load_XMM_TOAs(args.eventfile, minmjd=minmjd, maxmjd=maxmjd) elif hdr["TELESCOP"].lower().startswith("nustar"): if args.orbfile is not None: log.info("Setting up NuSTAR observatory") get_satellite_observatory("NuSTAR", args.orbfile) tl = load_NuSTAR_TOAs(args.eventfile, minmjd=minmjd, maxmjd=maxmjd) else: log.error( "FITS file not recognized, TELESCOPE = {0}, INSTRUMENT = {1}". format(hdr["TELESCOP"], hdr["INSTRUME"])) sys.exit(1) # Now convert to TOAs object and compute TDBs and posvels if len(tl) == 0: log.error("No TOAs, exiting!") sys.exit(0) # Read in model modelin = pint.models.get_model(args.parfile) use_planets = False if "PLANET_SHAPIRO" in modelin.params: if modelin.PLANET_SHAPIRO.value: use_planets = True if "AbsPhase" not in modelin.components: log.error( "TimingModel does not include AbsPhase component, which is required " "for computing phases. Make sure you have TZR* parameters in your par file!" ) raise ValueError("Model missing AbsPhase component.") if args.addorbphase and (not hasattr(modelin, "binary_model_name")): log.error( "TimingModel does not include a binary model, which is required for " "computing orbital phases. Make sure you have BINARY and associated " "model parameters in your par file!") raise ValueError("Model missing BINARY component.") ts = toa.get_TOAs_list( tl, ephem=args.ephem, include_bipm=args.use_bipm, include_gps=args.use_gps, planets=use_planets, tdb_method=args.tdbmethod, ) ts.filename = args.eventfile # if args.fix: # ts.adjust_TOAs(TimeDelta(np.ones(len(ts.table))*-1.0*u.s,scale='tt')) print(ts.get_summary()) mjds = ts.get_mjds() print(mjds.min(), mjds.max()) # Compute model phase for each TOA iphss, phss = modelin.phase(ts, abs_phase=True) phases = phss.value % 1 h = float(hm(phases)) print("Htest : {0:.2f} ({1:.2f} sigma)".format(h, h2sig(h))) if args.plot: phaseogram_binned(mjds, phases, bins=100, plotfile=args.plotfile) # Compute orbital phases for each photon TOA if args.addorbphase: delay = modelin.delay(ts) orbits = modelin.binary_instance.orbits() # These lines are already in orbits.orbit_phase() in binary_orbits.py. # What is the correct syntax is to call this function here? norbits = np.array(np.floor(orbits), dtype=np.long) orbphases = orbits - norbits # fractional phase if args.addphase or args.addorbphase: # 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) datacol = [] data_to_add = {} # Handle case where minMJD/maxMJD do not exceed length of events mjds_float = read_fits_event_mjds(hdulist[1]) time_mask = np.logical_and((mjds_float > minmjd), (mjds_float < maxmjd)) if args.addphase: if time_mask.sum() != len(phases): raise RuntimeError( "Mismatch between data selection {0} and length of phase array ({1})!" .format(time_mask.sum(), len(phases))) data_to_add["PULSE_PHASE"] = [phases, "D"] if args.absphase: data_to_add["ABS_PHASE"] = [iphss - negmask, "K"] if args.barytime: bats = modelin.get_barycentric_toas(ts) data_to_add["BARY_TIME"] = [bats, "D"] if args.addorbphase: if time_mask.sum() != len(orbphases): raise RuntimeError( "Mismatch between data selection ({0}) and length of orbital phase array ({1})!" .format(time_mask.sum(), len(orbphases))) data_to_add["ORBIT_PHASE"] = [orbphases, "D"] # End if args.addorbphase for key in data_to_add.keys(): if key in hdulist[1].columns.names: log.info("Found existing %s column, overwriting..." % key) # Overwrite values in existing Column hdulist[1].data[key][time_mask] = data_to_add[key][0] else: # Construct and append new column, preserving HDU header and name log.info("Adding new %s column." % key) new_dat = np.full(time_mask.shape, -1, dtype=data_to_add[key][0].dtype) new_dat[time_mask] = data_to_add[key][0] datacol.append( pyfits.ColDefs([ pyfits.Column(name=key, format=data_to_add[key][1], array=new_dat) ])) if len(datacol) > 0: cols = hdulist[1].columns for c in datacol: cols = cols + c bt = pyfits.BinTableHDU.from_columns(cols, header=hdulist[1].header, name=hdulist[1].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")
def main(argv=None): import argparse parser = argparse.ArgumentParser( description= "Use PINT to compute event phases and make plots of photon event files.", formatter_class=argparse.ArgumentDefaultsHelpFormatter, ) parser.add_argument( "eventfile", help= "Photon event FITS file name (e.g. from NICER, RXTE, XMM, Chandra).", ) parser.add_argument("parfile", help="par file to construct model from") parser.add_argument("--orbfile", help="Name of orbit file", default=None) parser.add_argument("--maxMJD", help="Maximum MJD to include in analysis", default=None) parser.add_argument("--minMJD", help="Minimum MJD to include in analysis", default=None) parser.add_argument("--plotfile", help="Output figure file name", default=None) parser.add_argument( "--addphase", help="Write FITS file with added phase column", default=False, action="store_true", ) parser.add_argument( "--addorbphase", help="Write FITS file with added orbital phase column", default=False, action="store_true", ) parser.add_argument( "--absphase", help="Write FITS file with integral portion of pulse phase (ABS_PHASE)", default=False, action="store_true", ) parser.add_argument( "--barytime", help= "Write FITS file with a column containing the barycentric time as double precision MJD.", default=False, action="store_true", ) parser.add_argument( "--outfile", help="Output FITS file name (default=same as eventfile)", default=None, ) parser.add_argument("--ephem", help="Planetary ephemeris to use", default="DE421") parser.add_argument( "--tdbmethod", help="Method for computing TT to TDB (default=astropy)", default="default", ) parser.add_argument("--plot", help="Show phaseogram plot.", action="store_true", default=False) 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)", ) parser.add_argument( "--log-level", type=str, choices=("TRACE", "DEBUG", "INFO", "WARNING", "ERROR"), default="WARNING", help="Logging level", dest="loglevel", ) # parser.add_argument("--fix",help="Apply 1.0 second offset for NICER", action='store_true', default=False) parser.add_argument( "--polycos", default=False, action="store_true", help= "Use polycos to calculate phases; use when working with very large event files", ) args = parser.parse_args(argv) log.remove() log.add( sys.stderr, level=args.loglevel, colorize=True, format=pint.logging.format, filter=pint.logging.LogFilter(), ) # If outfile is specified, that implies addphase if args.outfile is not None: args.addphase = True # If plotfile is specified, that implies plot if args.plotfile is not None: args.plot = True # set MJD ranges maxmjd = np.inf if (args.maxMJD is None) else float(args.maxMJD) minmjd = 0.0 if (args.minMJD is None) else float(args.minMJD) # Read event file header to figure out what instrument is is from hdr = pyfits.getheader(args.eventfile, ext=1) log.info("Event file TELESCOPE = {0}, INSTRUMENT = {1}".format( hdr["TELESCOP"], hdr["INSTRUME"])) telescope = hdr["TELESCOP"].lower() # Instantiate observatory once so it gets added to the observatory registry if args.orbfile is not None: log.info(f"Setting up {telescope.upper()} observatory") try: get_satellite_observatory(telescope, args.orbfile) except Exception: log.error( "The orbit file is not recognized. It is likely that this mission is not supported. " "Please barycenter the event file using the official mission tools before processing with PINT" ) # Read event file and return list of TOA objects, if not using polycos if args.polycos == False: try: tl = load_event_TOAs(args.eventfile, telescope, minmjd=minmjd, maxmjd=maxmjd) except KeyError: log.error( "Observatory not recognized. This probably means you need to provide an orbit file or barycenter the event file." ) sys.exit(1) # Now convert to TOAs object and compute TDBs and posvels if len(tl) == 0: log.error("No TOAs, exiting!") sys.exit(0) # Read in model modelin = pint.models.get_model(args.parfile) use_planets = False if "PLANET_SHAPIRO" in modelin.params: if modelin.PLANET_SHAPIRO.value: use_planets = True if "AbsPhase" not in modelin.components: log.error( "TimingModel does not include AbsPhase component, which is required " "for computing phases. Make sure you have TZR* parameters in your par file!" ) raise ValueError("Model missing AbsPhase component.") if args.addorbphase and (not hasattr(modelin, "binary_model_name")): log.error( "TimingModel does not include a binary model, which is required for " "computing orbital phases. Make sure you have BINARY and associated " "model parameters in your par file!") raise ValueError("Model missing BINARY component.") # Use polycos to calculate pulse phases if args.polycos: log.info("Using polycos to get pulse phases.") if args.addorbphase: raise ValueError("Cannot use orbphase with polycos.") # Polycos parameters segLength = 120 # in minutes ncoeff = 10 obsfreq = 0 # Open event file and get start and end mjds hdulist = pyfits.open(args.eventfile) data = hdulist[1].data mjds = read_fits_event_mjds(hdulist[1]) minmjd = min(mjds) maxmjd = max(mjds) # Check if event file is barycentered if hdulist[1].header["TIMESYS"] != "TDB": raise ValueError( "The event file has not been barycentered. Polycos can only be used with barycentered data." ) # Create polycos table log.debug("Generating polycos") telescope_n = "@" p = polycos.Polycos() ptable = p.generate_polycos(modelin, minmjd, maxmjd, telescope_n, segLength, ncoeff, obsfreq) # Calculate phases log.debug("Evaluating polycos") phases = p.eval_phase(mjds) phases[phases < 0] += 1.0 h = float(hm(phases)) print("Htest : {0:.2f} ({1:.2f} sigma)".format(h, h2sig(h))) else: # Normal mode, not polycos ts = toa.get_TOAs_list( tl, ephem=args.ephem, include_bipm=args.use_bipm, include_gps=args.use_gps, planets=use_planets, tdb_method=args.tdbmethod, ) ts.filename = args.eventfile # if args.fix: # ts.adjust_TOAs(TimeDelta(np.ones(len(ts.table))*-1.0*u.s,scale='tt')) print(ts.get_summary()) mjds = ts.get_mjds() print(mjds.min(), mjds.max()) # Compute model phase for each TOA iphss, phss = modelin.phase(ts, abs_phase=True) phases = phss.value % 1 h = float(hm(phases)) print("Htest : {0:.2f} ({1:.2f} sigma)".format(h, h2sig(h))) if args.plot: phaseogram_binned(mjds, phases, bins=100, plotfile=args.plotfile) # Compute orbital phases for each photon TOA if args.addorbphase: delay = modelin.delay(ts) orbits = modelin.binary_instance.orbits() # These lines are already in orbits.orbit_phase() in binary_orbits.py. # What is the correct syntax is to call this function here? norbits = np.array(np.floor(orbits), dtype=int) orbphases = orbits - norbits # fractional phase if args.addphase or args.addorbphase: # 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) datacol = [] data_to_add = {} # Handle case where minMJD/maxMJD do not exceed length of events mjds_float = read_fits_event_mjds(hdulist[1]) time_mask = np.logical_and((mjds_float > minmjd), (mjds_float < maxmjd)) if args.polycos: time_mask = np.logical_and((mjds_float >= minmjd), (mjds_float <= maxmjd)) if args.addphase: if time_mask.sum() != len(phases): raise RuntimeError( "Mismatch between data selection {0} and length of phase array ({1})!" .format(time_mask.sum(), len(phases))) data_to_add["PULSE_PHASE"] = [phases, "D"] if args.absphase: data_to_add["ABS_PHASE"] = [iphss, "K"] if args.barytime: bats = modelin.get_barycentric_toas(ts) data_to_add["BARY_TIME"] = [bats, "D"] if args.addorbphase: if time_mask.sum() != len(orbphases): raise RuntimeError( "Mismatch between data selection ({0}) and length of orbital phase array ({1})!" .format(time_mask.sum(), len(orbphases))) data_to_add["ORBIT_PHASE"] = [orbphases, "D"] # End if args.addorbphase for key in data_to_add.keys(): if key in hdulist[1].columns.names: log.info("Found existing %s column, overwriting..." % key) # Overwrite values in existing Column hdulist[1].data[key][time_mask] = data_to_add[key][0] else: # Construct and append new column, preserving HDU header and name log.info("Adding new %s column." % key) new_dat = np.full(time_mask.shape, -1, dtype=data_to_add[key][0].dtype) new_dat[time_mask] = data_to_add[key][0] datacol.append( pyfits.ColDefs([ pyfits.Column(name=key, format=data_to_add[key][1], array=new_dat) ])) if len(datacol) > 0: cols = hdulist[1].columns for c in datacol: cols = cols + c bt = pyfits.BinTableHDU.from_columns(cols, header=hdulist[1].header, name=hdulist[1].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")
evhdr = hdulist[1].header evdat = hdulist[1].data # Hack for XMM barycentered data if evhdr['TELESCOP'].startswith('XMM'): log.info('XMM data, setting MET0 to 50814.0') MET0 = Time(50814.0, format='mjd', scale='tdb') # Hack for NuSTAR barycentered data if evhdr['TELESCOP'].startswith('NuSTAR'): log.info('NuSTAR data, setting MET0 to 55197.00076601852') MET0 = Time(55197.00076601852, format='mjd', scale='tdb') #mets = evdat.field('TIME') #mets += evhdr['TIMEZERO'] mjds = read_fits_event_mjds(hdulist[1]) #log.info('MJDs {0}'.format(mjds[:10])) #evtimes = MET0 + mets*u.s #mjds = evtimes.mjd #log.info('Evtimes {0}'.format(evtimes[:10])) try: phases = evdat.field('PULSE_PHASE') except: phases = evdat.field('PHASE') if options.weights is not None: weights = evdat.field(options.weights) # FILTER ON ENERGY try: en = evdat.field('PI') * PI_TO_KEV